Fork me on GitHub

source: git/examples/geometry.C@ b3c42d3

ImprovedOutputFile Timing dual_readout llp
Last change on this file since b3c42d3 was b3c42d3, checked in by Christophe Delaere <christophe.delaere@…>, 10 years ago

Broken version on the way to automatization

  • Property mode set to 100644
File size: 41.4 KB
Line 
1#include <set>
2#include <map>
3#include <utility>
4#include <vector>
5#include <algorithm>
6#include <sstream>
7#include <exception>
8#include "TGeoManager.h"
9#include "TGeoVolume.h"
10#include "TGeoMedium.h"
11#include "TGeoNode.h"
12#include "TGeoCompositeShape.h"
13#include "TGeoMatrix.h"
14#include "TGeoTube.h"
15#include "TGeoCone.h"
16#include "TGeoArb8.h"
17//#include "external/ExRootAnalysis/ExRootConfReader.h"
18//#include "external/ExRootAnalysis/ExRootTreeReader.h"
19//#include "display/DelphesCaloData.h"
20//#include "display/DelphesDisplay.h"
21//#include "display/DelphesBranchElement.h"
22//#include "classes/DelphesClasses.h"
23#include "TF2.h"
24#include "TH1F.h"
25#include "TChain.h"
26#include "TEveElement.h"
27#include "TEveJetCone.h"
28#include "TEveTrack.h"
29#include "TEveTrackPropagator.h"
30#include "TEveCalo.h"
31#include "TMath.h"
32#include "TSystem.h"
33#include "TEveManager.h"
34#include "TEveGeoNode.h"
35#include "TEveTrans.h"
36#include "TEveViewer.h"
37#include "TEveBrowser.h"
38#include "TRootBrowser.h"
39#include "TGLViewer.h"
40#include "TGButton.h"
41#include "TCollection.h"
42#include "TClonesArray.h"
43#include "TGLClip.h"
44#include "TEveArrow.h"
45
46/*
47 * assembly.C: sauvegarde as shape-extract -> implement in the geometry class (read/write)
48 * histobrowser.C: intégration d'histogrammes dans le display (on pourrait avoir Pt, eta, phi pour les principales collections)
49 * also from alice_esd: summary html table
50 *
51 */
52using namespace std;
53
54// Forward declarations.
55class Delphes3DGeometry;
56class ExRootTreeReader;
57class DelphesCaloData;
58class DelphesDisplay;
59class DelphesBranchBase;
60template<typename EveContainer> class DelphesBranchElement;
61void make_gui();
62void load_event();
63void delphes_read();
64void readConfig(const char *configFile, const Delphes3DGeometry& det3D, std::vector<DelphesBranchBase*>& elements, std::vector<TClonesArray*>& arrays);
65
66// Configuration and global variables.
67Int_t event_id = 0; // Current event id.
68Double_t gRadius = 1.29;
69Double_t gTotRadius = 2.0;
70Double_t gHalfLength = 3.0;
71Double_t gBz = 3.8;
72
73TChain gChain("Delphes");
74
75ExRootTreeReader *gTreeReader = 0;
76
77std::vector<DelphesBranchBase*> gElements;
78std::vector<TClonesArray*> gArrays;
79
80DelphesDisplay *gDelphesDisplay = 0;
81
82/******************************************************************************/
83// Construction of the geometry
84/******************************************************************************/
85
86// TODO: asymmetric detector
87
88class Delphes3DGeometry {
89 public:
90 Delphes3DGeometry(TGeoManager *geom = NULL);
91 ~Delphes3DGeometry() {}
92
93 void readFile(const char* filename, const char* ParticlePropagator="ParticlePropagator",
94 const char* TrackingEfficiency="ChargedHadronTrackingEfficiency",
95 const char* MuonEfficiency="MuonEfficiency",
96 const char* Calorimeters="Calorimeter");
97
98 void loadFromFile(const char* filename, const char* name="DelphesGeometry");
99 void save(const char* filename, const char* name="DelphesGeometry");
100
101 void setContingency(Double_t contingency) { contingency_ = contingency; }
102 void setCaloBarrelThickness(Double_t thickness) { calo_barrel_thickness_ = thickness; }
103 void setCaloEndcapThickness(Double_t thickness) { calo_endcap_thickness_ = thickness; }
104 void setMuonSystemThickness(Double_t thickness) { muonSystem_thickn_ = thickness; }
105
106 TGeoVolume* getDetector(bool withTowers = true);
107
108 Double_t getTrackerRadius() const { return tk_radius_; }
109 Double_t getDetectorRadius() const { return muonSystem_radius_; }
110 Double_t getTrackerHalfLength() const { return tk_length_; }
111 Double_t getDetectorHalfLength() const { return muonSystem_length_; }
112 Double_t getBField() const { return tk_Bz_; }
113 std::pair<TAxis*, TAxis*> getCaloAxes() { return std::make_pair(etaAxis_,phiAxis_); }
114
115 private:
116 std::pair<Double_t, Double_t> addTracker(TGeoVolume *top);
117 std::pair<Double_t, Double_t> addCalorimeter(TGeoVolume *top, const char* name, Double_t innerBarrelRadius, Double_t innerBarrelLength, set< pair<Double_t, Int_t> >& caloBinning);
118 std::pair<Double_t, Double_t> addMuonDets(TGeoVolume *top, const char* name, Double_t innerBarrelRadius, Double_t innerBarrelLength);
119 void addCaloTowers(TGeoVolume *top, const char* name, Double_t innerBarrelRadius, Double_t innerBarrelLength, set< pair<Double_t, Int_t> >& caloBinning);
120
121 private:
122
123 TGeoManager *geom_;
124
125 TGeoMedium *vacuum_;
126 TGeoMedium *tkmed_;
127 TGeoMedium *calomed_;
128 TGeoMedium *mudetmed_;
129
130 TAxis* etaAxis_;
131 TAxis* phiAxis_;
132
133 Double_t contingency_;
134 Double_t calo_barrel_thickness_;
135 Double_t calo_endcap_thickness_;
136 Double_t muonSystem_thickn_;
137 Double_t muonSystem_radius_;
138 Double_t muonSystem_length_;
139 Double_t tk_radius_;
140 Double_t tk_length_;
141 Double_t tk_etamax_;
142 Double_t tk_Bz_;
143
144 std::vector<std::string> calorimeters_;
145 std::vector<std::string> muondets_;
146
147 std::map<std::string, Double_t> muonSystem_etamax_;
148 std::map<std::string, set< pair<Double_t, Int_t> > > caloBinning_;
149
150};
151
152Delphes3DGeometry::Delphes3DGeometry(TGeoManager *geom) {
153
154 //--- the geometry manager
155 geom_ = geom==NULL? gGeoManager : geom;
156 //gGeoManager->DefaultColors();
157
158 //--- define some materials
159 TGeoMaterial *matVacuum = new TGeoMaterial("Vacuum", 0,0,0);
160 TGeoMaterial *matAl = new TGeoMaterial("Al", 26.98,13,2.7); // placeholder
161 //TODO: create different materials for different subdetectors???
162 matVacuum->SetTransparency(85); //TODO: tune
163 matAl->SetTransparency(85); //TODO: tune
164
165 //--- define some media
166 TGeoMedium *Vacuum = new TGeoMedium("Vacuum",1, matVacuum);
167 TGeoMedium *Al = new TGeoMedium("Root Material",2, matAl);
168 vacuum_ = Vacuum;
169 tkmed_ = Vacuum; // placeholder
170 calomed_ = Al; // placeholder
171 mudetmed_ = Al; // placeholder
172
173 // custom parameters
174 contingency_ = 10.;
175 calo_barrel_thickness_ = 50.;
176 calo_endcap_thickness_ = 75.;
177 muonSystem_thickn_ = 10.;
178
179 // read these parameters from the Delphes Card (with default values)
180 etaAxis_ = NULL;
181 phiAxis_ = NULL;
182 tk_radius_ = 120.;
183 tk_length_ = 150.;
184 tk_etamax_ = 3.0;
185 tk_Bz_ = 1.;
186 muonSystem_radius_ = 200.;
187}
188
189void Delphes3DGeometry::readFile(const char *configFile,
190 const char* ParticlePropagator, const char* TrackingEfficiency,
191 const char* MuonEfficiency, const char* Calorimeters) {
192
193 ExRootConfReader *confReader = new ExRootConfReader;
194 confReader->ReadFile(configFile);
195
196 tk_radius_ = confReader->GetDouble(Form("%s::Radius",ParticlePropagator), 1.0)*100.; // tk_radius
197 tk_length_ = confReader->GetDouble(Form("%s::HalfLength",ParticlePropagator), 3.0)*100.; // tk_length
198 tk_Bz_ = confReader->GetDouble("ParticlePropagator::Bz", 0.0); // tk_Bz
199
200 {
201 TString tkEffFormula = confReader->GetString(Form("%s::EfficiencyFormula",TrackingEfficiency),"abs(eta)<3.0");
202 tkEffFormula.ReplaceAll("pt","x");
203 tkEffFormula.ReplaceAll("eta","y");
204 tkEffFormula.ReplaceAll("phi","0.");
205 TF2* tkEffFunction = new TF2("tkEff",tkEffFormula,0,1000,-10,10);
206 TH1F etaHisto("eta","eta",100,5.,-5.);
207 Double_t pt,eta;
208 for(int i=0;i<1000;++i) {
209 tkEffFunction->GetRandom2(pt,eta);
210 etaHisto.Fill(eta);
211 }
212 Int_t bin = -1;
213 bin = etaHisto.FindFirstBinAbove(0.5);
214 Double_t etamin = (bin>-1) ? etaHisto.GetBinLowEdge(bin) : -10.;
215 bin = etaHisto.FindLastBinAbove(0.5);
216 Double_t etamax = (bin>-1) ? etaHisto.GetBinLowEdge(bin+1) : -10.;
217 tk_etamax_ = TMath::Max(fabs(etamin),fabs(etamax)); // tk_etamax
218 delete tkEffFunction;
219 }
220
221 {
222 muondets_.push_back("muons");
223 TString muonEffFormula = confReader->GetString(Form("%s::EfficiencyFormula",MuonEfficiency),"abs(eta)<2.0");
224 muonEffFormula.ReplaceAll("pt","x");
225 muonEffFormula.ReplaceAll("eta","y");
226 muonEffFormula.ReplaceAll("phi","0.");
227 TF2* muEffFunction = new TF2("muEff",muonEffFormula,0,1000,-10,10);
228 TH1F etaHisto("eta2","eta2",100,5.,-5.);
229 Double_t pt,eta;
230 for(int i=0;i<1000;++i) {
231 muEffFunction->GetRandom2(pt,eta);
232 etaHisto.Fill(eta);
233 }
234 Int_t bin = -1;
235 bin = etaHisto.FindFirstBinAbove(0.5);
236 Double_t etamin = (bin>-1) ? etaHisto.GetBinLowEdge(bin) : -10.;
237 bin = etaHisto.FindLastBinAbove(0.5);
238 Double_t etamax = (bin>-1) ? etaHisto.GetBinLowEdge(bin+1) : -10.;
239 muonSystem_etamax_["muons"] = TMath::Max(fabs(etamin),fabs(etamax)); // muonSystem_etamax
240 delete muEffFunction;
241 }
242
243 std::string s(Calorimeters);
244 std::replace( s.begin(), s.end(), ',', ' ' );
245 std::istringstream stream( s );
246 std::string word;
247 while (stream >> word) calorimeters_.push_back(word);
248
249 caloBinning_.clear(); // calo binning
250 for(std::vector<std::string>::const_iterator calo=calorimeters_.begin();calo!=calorimeters_.end(); ++calo) {
251 set< pair<Double_t, Int_t> > caloBinning;
252 ExRootConfParam paramEtaBins, paramPhiBins;
253 ExRootConfParam param = confReader->GetParam(Form("%s::EtaPhiBins",calo->c_str()));
254 Int_t size = param.GetSize();
255 for(int i = 0; i < size/2; ++i) {
256 paramEtaBins = param[i*2];
257 paramPhiBins = param[i*2+1];
258 assert(paramEtaBins.GetSize()==1);
259 caloBinning.insert(std::make_pair(paramEtaBins[0].GetDouble(),paramPhiBins.GetSize()-1));
260 }
261 caloBinning_[*calo] = caloBinning;
262 }
263
264 set< pair<Double_t, Int_t> > caloBinning = caloBinning_[*calorimeters_.begin()];
265 Double_t *etaBins = new Double_t[caloBinning.size()]; // note that this is the eta binning of the first calo
266 unsigned int ii = 0;
267 for(set< pair<Double_t, Int_t> >::const_iterator itEtaSet = caloBinning.begin(); itEtaSet != caloBinning.end(); ++itEtaSet) {
268 etaBins[ii++] = itEtaSet->first;
269 }
270 etaAxis_ = new TAxis(caloBinning.size() - 1, etaBins);
271 phiAxis_ = new TAxis(72, -TMath::Pi(), TMath::Pi()); // note that this is fixed while #phibins could vary, also with eta, which doesn't seem possible in ROOT
272
273 muonSystem_radius_ = tk_radius_ + contingency_ + (contingency_+calo_barrel_thickness_)*calorimeters_.size() + muonSystem_thickn_;
274 muonSystem_length_ = tk_length_ + contingency_ + (contingency_+calo_endcap_thickness_)*calorimeters_.size() + muonSystem_thickn_;
275
276 delete confReader;
277
278}
279
280TGeoVolume* Delphes3DGeometry::getDetector(bool withTowers) {
281 // compute the envelope
282 Double_t system_radius = tk_radius_+calo_barrel_thickness_+3*contingency_;
283 Double_t system_length = tk_length_+contingency_+(contingency_+calo_endcap_thickness_)*calorimeters_.size()+contingency_;
284 // the detector volume
285 TGeoVolume *top = geom_->MakeBox("Delphes3DGeometry", vacuum_, system_radius, system_radius, system_length);
286 // build the detector
287 std::pair<Double_t, Double_t> limits = addTracker(top);
288 Double_t radius = limits.first;
289 Double_t length = limits.second;
290 for(std::vector<std::string>::const_iterator calo = calorimeters_.begin(); calo != calorimeters_.end(); ++calo) {
291 limits = addCalorimeter(top,calo->c_str(),radius,length,caloBinning_[*calo]);
292 if (withTowers) {
293 addCaloTowers(top,calo->c_str(),radius,length,caloBinning_[*calo]);
294 }
295 radius = limits.first;
296 length = limits.second;
297 }
298 for(std::vector<std::string>::const_iterator muon = muondets_.begin(); muon != muondets_.end(); ++muon) {
299 limits = addMuonDets(top, muon->c_str(), radius, length);
300 radius = limits.first;
301 length = limits.second;
302 }
303 // return the result
304 return top;
305}
306
307std::pair<Double_t, Double_t> Delphes3DGeometry::addTracker(TGeoVolume *top) {
308 // tracker: a cylinder with two cones substracted
309 new TGeoCone("forwardTkAcceptance",(tk_length_/2.+0.05),0.,tk_radius_,(tk_length_)*2.*exp(-tk_etamax_)/(1-exp(-2.*tk_etamax_)),tk_radius_);
310 TGeoTranslation *tr1 = new TGeoTranslation("tkacc1",0., 0., tk_length_/2.);
311 tr1->RegisterYourself();
312 TGeoRotation *negz = new TGeoRotation("tknegz",0,180,0);
313 negz->RegisterYourself();
314 TGeoCombiTrans *tr2 = new TGeoCombiTrans("tkacc2",0.,0.,-tk_length_/2.,negz);
315 tr2->RegisterYourself();
316 TGeoCompositeShape* tracker_cs = new TGeoCompositeShape("tracker_cs","forwardTkAcceptance:tkacc1+forwardTkAcceptance:tkacc2");
317 TGeoVolume *tracker = new TGeoVolume("tracker",tracker_cs,tkmed_);
318 tracker->SetLineColor(kYellow);
319 top->AddNode(tracker,1);
320 return std::make_pair(tk_radius_,tk_length_);
321}
322
323std::pair<Double_t, Double_t> Delphes3DGeometry::addCalorimeter(TGeoVolume *top, const char* name,
324 Double_t innerBarrelRadius, Double_t innerBarrelLength, set< pair<Double_t, Int_t> >& caloBinning) {
325 // parameters derived from the inputs
326 Double_t calo_endcap_etamax = TMath::Max(fabs(caloBinning.begin()->first),fabs(caloBinning.rbegin()->first));
327 Double_t calo_barrel_innerRadius = innerBarrelRadius+contingency_;
328 Double_t calo_barrel_length = innerBarrelLength + calo_barrel_thickness_;
329 Double_t calo_endcap_etamin = -log(innerBarrelRadius/(2*innerBarrelLength));
330 Double_t calo_endcap_innerRadius1 = innerBarrelLength*2.*exp(-calo_endcap_etamax)/(1-exp(-2.*calo_endcap_etamax));
331 Double_t calo_endcap_innerRadius2 = (innerBarrelLength+calo_endcap_thickness_)*2.*exp(-calo_endcap_etamax)/(1-exp(-2.*calo_endcap_etamax));
332 Double_t calo_endcap_outerRadius1 = innerBarrelRadius;
333 Double_t calo_endcap_outerRadius2 = innerBarrelRadius+calo_barrel_thickness_;
334 Double_t calo_endcap_coneThickness = TMath::Min(calo_barrel_thickness_ * (1-exp(-2.*calo_endcap_etamin)) / (2.*exp(-calo_endcap_etamin)), calo_endcap_thickness_);
335 Double_t calo_endcap_diskThickness = TMath::Max(0.,calo_endcap_thickness_-calo_endcap_coneThickness);
336
337 // calorimeters: tube truncated in eta + cones
338 new TGeoTube(Form("%s_barrel_cylinder",name),calo_barrel_innerRadius,calo_barrel_innerRadius+calo_barrel_thickness_,calo_barrel_length);
339 new TGeoCone(Form("%s_endcap_cone",name),calo_endcap_coneThickness/2.,calo_endcap_innerRadius1,calo_endcap_outerRadius1,calo_endcap_innerRadius2,calo_endcap_outerRadius2);
340 new TGeoTube(Form("%s_endcap_disk",name),calo_endcap_innerRadius2,tk_radius_+calo_barrel_thickness_,calo_endcap_diskThickness/2.);
341 TGeoTranslation *tr1 = new TGeoTranslation(Form("%s_tr1",name),0., 0., (calo_endcap_coneThickness+calo_endcap_diskThickness)/2.);
342 tr1->RegisterYourself();
343 TGeoCompositeShape *calo_endcap_cs = new TGeoCompositeShape(Form("%s_endcap_cs",name),Form("%s_endcap_cone+%s_endcap_disk:%s_tr1",name,name,name));
344 TGeoTranslation *trc1 = new TGeoTranslation(Form("%s_endcap1_position",name),0.,0., innerBarrelLength+calo_endcap_coneThickness/2.);
345 trc1->RegisterYourself();
346 TGeoRotation *negz = new TGeoRotation(Form("%s_negz",name),0,180,0);
347 TGeoCombiTrans *trc2 = new TGeoCombiTrans(Form("%s_endcap2_position",name),0.,0.,-(innerBarrelLength+calo_endcap_coneThickness/2.),negz);
348 trc2->RegisterYourself();
349 TGeoTranslation *trc1c = new TGeoTranslation(Form("%s_endcap1_position_cont",name),0.,0., innerBarrelLength+calo_endcap_coneThickness/2.+contingency_);
350 trc1c->RegisterYourself();
351 TGeoCombiTrans *trc2c = new TGeoCombiTrans(Form("%s_endcap2_position_cont",name),0.,0.,-(innerBarrelLength+calo_endcap_coneThickness/2.)-contingency_,negz);
352 trc2c->RegisterYourself();
353 TGeoVolume *calo_endcap = new TGeoVolume(Form("%s_endcap",name),calo_endcap_cs,calomed_);
354 TGeoCompositeShape *calo_barrel_cs = new TGeoCompositeShape(Form("%s_barrel_cs",name),
355 Form("%s_barrel_cylinder-%s_endcap_cs:%s_endcap1_position-%s_endcap_cs:%s_endcap2_position",name,name,name,name,name));
356 TGeoVolume *calo_barrel = new TGeoVolume(Form("%s_barrel",name),calo_barrel_cs,calomed_);
357 calo_endcap->SetLineColor(kViolet);
358 calo_endcap->SetFillColor(kViolet);
359 calo_barrel->SetLineColor(kRed);
360 top->AddNode(calo_endcap,1,trc1c);
361 top->AddNode(calo_endcap,2,trc2c);
362 top->AddNode(calo_barrel,1);
363 return std::make_pair(calo_barrel_innerRadius+calo_barrel_thickness_,innerBarrelLength+calo_endcap_thickness_+contingency_);
364}
365
366std::pair<Double_t, Double_t> Delphes3DGeometry::addMuonDets(TGeoVolume *top, const char* name, Double_t innerBarrelRadius, Double_t innerBarrelLength) {
367 // muon system: tube + disks
368 Double_t muonSystem_radius = innerBarrelRadius + contingency_;
369 Double_t muonSystem_length = innerBarrelLength + contingency_;
370 Double_t muonSystem_rmin = muonSystem_length*2.*exp(-muonSystem_etamax_[name])/(1-exp(-2.*muonSystem_etamax_[name]));
371 TGeoVolume *muon_barrel = geom_->MakeTube(Form("%s_barrel",name),mudetmed_,muonSystem_radius,muonSystem_radius+muonSystem_thickn_,muonSystem_length);
372 muon_barrel->SetLineColor(kBlue);
373 top->AddNode(muon_barrel,1);
374 TGeoVolume *muon_endcap = geom_->MakeTube(Form("%s_endcap",name),mudetmed_,muonSystem_rmin,muonSystem_radius+muonSystem_thickn_,muonSystem_thickn_/2.);
375 muon_endcap->SetLineColor(kBlue);
376 TGeoTranslation *trm1 = new TGeoTranslation(Form("%sEndcap1_position",name),0.,0.,muonSystem_length);
377 trm1->RegisterYourself();
378 TGeoTranslation *trm2 = new TGeoTranslation(Form("%sEndcap2_position",name),0.,0.,-muonSystem_length);
379 trm1->RegisterYourself();
380 top->AddNode(muon_endcap,1,trm1);
381 top->AddNode(muon_endcap,2,trm2);
382 return std::make_pair(muonSystem_radius,muonSystem_length);
383}
384
385void Delphes3DGeometry::addCaloTowers(TGeoVolume *top, const char* name,
386 Double_t innerBarrelRadius, Double_t innerBarrelLength, set< pair<Double_t, Int_t> >& caloBinning) {
387
388 TGeoVolume* calo_endcap = top->GetNode(Form("%s_endcap_1",name))->GetVolume();
389 TGeoVolume* calo_barrel = top->GetNode(Form("%s_barrel_1",name))->GetVolume();
390 Double_t calo_endcap_etamin = -log(innerBarrelRadius/(2*innerBarrelLength));
391 Double_t calo_endcap_coneThickness = TMath::Min(calo_barrel_thickness_ * (1-exp(-2.*calo_endcap_etamin)) / (2.*exp(-calo_endcap_etamin)), calo_endcap_thickness_);
392
393 // calo towers in the barrel
394 Double_t vertices[16] = {0.,0.,0.,0.,0.,0.,0.,0.}; // summit of the pyramid
395 Double_t R = tk_radius_ + contingency_+(contingency_+calo_barrel_thickness_)*calorimeters_.size(); // radius of the muons system = height of the pyramid
396 Int_t nEtaBins = caloBinning.size();
397 // this rotation is to make the tower point "up"
398 TGeoRotation* initTowerRot = new TGeoRotation(Form("%s_initTowerRot",name),0.,90.,0.);
399 TGeoCombiTrans* initTower = new TGeoCombiTrans(Form("%s_initTower",name),0.,-R/2.,0.,initTowerRot);
400 initTower->RegisterYourself();
401 // eta bins... we build one pyramid per eta slice and then translate it nphi times.
402 // phi bins represented by rotations around z
403 Double_t *y = new Double_t[nEtaBins];
404 Double_t *dx = new Double_t[nEtaBins];
405 Int_t *nphi = new Int_t[nEtaBins];
406 Int_t etaslice = 0;
407 std::map<std::pair<int,int>, TGeoRotation*> phirotations;
408 for(set< pair<Double_t, Int_t> >::const_iterator bin=caloBinning.begin(); bin!=caloBinning.end();++bin) {
409 if(abs(bin->first)>calo_endcap_etamin) continue; // only in the barrel
410 nphi[etaslice] = bin->second;
411 y[etaslice] = 0.5*R*(1-exp(-2*bin->first))/exp(-bin->first);
412 Double_t phiRotationAngle = 360./nphi[etaslice];
413 dx[etaslice] = R*tan(TMath::Pi()*phiRotationAngle/360.);
414 for(int phislice=0;phislice<nphi[etaslice];++phislice) {
415 phirotations[make_pair(etaslice,phislice)] = new TGeoRotation(Form("%s_phi%d_%d",name,etaslice,phislice),phiRotationAngle*phislice,0.,0.);
416 phirotations[make_pair(etaslice,phislice)]->RegisterYourself();
417 }
418 ++etaslice;
419 }
420 nEtaBins = etaslice;
421 for(int i=0;i<nEtaBins-1;++i) { // loop on the eta slices
422 vertices[8] = -dx[i]; vertices[9] = y[i];
423 vertices[10] = -dx[i]; vertices[11] = y[i+1];
424 vertices[12] = dx[i]; vertices[13] = y[i+1];
425 vertices[14] = dx[i]; vertices[15] = y[i];
426 new TGeoArb8(Form("%s_tower%d",name,i),R/2., vertices); // tower in the proper eta slice, at phi=0
427 // intersection between the tower and the calo_barrel
428 TGeoCompositeShape *finaltower_cs = new TGeoCompositeShape(Form("%s_ftower%d_cs",name,i),Form("%s_tower%d:%s_initTower*%s_barrel_cs",name,i,name,name));
429 TGeoVolume *finaltower = new TGeoVolume(Form("%s_ftower%d",name,i),finaltower_cs,calomed_);
430 finaltower->SetLineColor(kRed);
431 for(int j=0;j<nphi[i];++j) { // loop on the phi slices
432 calo_barrel->AddNode(finaltower,j,phirotations[make_pair(i,j)]);
433 }
434 }
435 delete[] y;
436 delete[] dx;
437 delete[] nphi;
438 //the towers in the forward region
439 R = tk_length_+contingency_+(contingency_+calo_endcap_thickness_)*calorimeters_.size(); // Z of the muons system = height of the pyramid
440 nEtaBins = caloBinning.size();
441 // translation to bring the origin of the tower to (0,0,0) (well, not really as the endcap is not yet in place)
442 TGeoTranslation* towerdz = new TGeoTranslation(Form("%s_towerdz",name),0.,0.,R/2.-(innerBarrelLength+calo_endcap_coneThickness/2.));
443 towerdz->RegisterYourself();
444 // eta bins... we build one pyramid per eta slice and then translate it nphi times.
445 Double_t *r = new Double_t[nEtaBins];
446 nphi = new Int_t[nEtaBins];
447 etaslice = 0;
448 phirotations.clear();
449 for(set< pair<Double_t, Int_t> >::const_iterator bin=caloBinning.begin(); bin!=caloBinning.end();++bin) {
450 if(bin->first<calo_endcap_etamin) continue; // only in the + endcap
451 r[etaslice] = R*2*exp(-bin->first)/(1-exp(-2*bin->first));
452 nphi[etaslice] = bin->second;
453 Double_t phiRotationAngle = 360./nphi[etaslice];
454 for(int phislice=0;phislice<nphi[etaslice];++phislice) {
455 phirotations[make_pair(etaslice,phislice)] = new TGeoRotation(Form("%s_forward_phi%d_%d",name,etaslice,phislice),phiRotationAngle*phislice,0.,0.);
456 phirotations[make_pair(etaslice,phislice)]->RegisterYourself();
457 }
458 ++etaslice;
459 }
460 nEtaBins = etaslice;
461 for(int i=0;i<nEtaBins-1;++i) { // loop on the eta slices
462 vertices[8] = -r[i+1]*sin(TMath::Pi()/nphi[i]); vertices[9] = r[i+1]*cos(TMath::Pi()/nphi[i]);
463 vertices[10] = -r[i]*sin(TMath::Pi()/nphi[i]); vertices[11] = r[i]*cos(TMath::Pi()/nphi[i]);
464 vertices[12] = r[i]*sin(TMath::Pi()/nphi[i]); vertices[13] = r[i]*cos(TMath::Pi()/nphi[i]);
465 vertices[14] = r[i+1]*sin(TMath::Pi()/nphi[i]); vertices[15] = r[i+1]*cos(TMath::Pi()/nphi[i]);
466 new TGeoArb8(Form("%sfwdtower%d",name,i),R/2., vertices); // tower in the proper eta slice, at phi=0
467 // intersection between the tower and the calo_endcap
468 TGeoCompositeShape *finalfwdtower_cs = new TGeoCompositeShape(Form("%sffwdtower%d_cs",name,i),Form("%sfwdtower%d:%s_towerdz*%s_endcap_cs",name,i,name,name));
469 TGeoVolume *finalfwdtower = new TGeoVolume(Form("%sffwdtower%d",name,i),finalfwdtower_cs,calomed_);
470 finalfwdtower->SetLineColor(kViolet);
471 for(int j=0;j<nphi[i];++j) { // loop on the phi slices
472 calo_endcap->AddNode(finalfwdtower,j,phirotations[make_pair(i,j)]);
473 }
474 }
475 delete[] r;
476 delete[] nphi;
477}
478
479/******************************************************************************/
480// Initialization and steering functions
481/******************************************************************************/
482
483// function that parses the config to extract the branches of interest and prepare containers
484void readConfig(const char *configFile, Delphes3DGeometry& det3D, std::vector<DelphesBranchBase*>& elements, std::vector<TClonesArray*>& arrays) {
485 ExRootConfReader *confReader = new ExRootConfReader;
486 confReader->ReadFile(configFile);
487 Double_t tk_radius = det3D.getTrackerRadius();
488 Double_t tk_length = det3D.getTrackerHalfLength();
489 Double_t tk_Bz = det3D.getBField();
490 Double_t mu_radius = det3D.getDetectorRadius();
491 Double_t mu_length = det3D.getDetectorHalfLength();
492 TAxis* etaAxis = det3D.getCaloAxes().first;
493 TAxis* phiAxis = det3D.getCaloAxes().second;
494 ExRootConfParam branches = confReader->GetParam("TreeWriter::Branch");
495 Int_t nBranches = branches.GetSize()/3;
496 DelphesBranchElement<TEveTrackList>* tlist;
497 DelphesBranchElement<DelphesCaloData>* clist;
498 DelphesBranchElement<TEveTrackList>* elist;
499 for(Int_t b = 0; b<nBranches; ++b) {
500 TString input = branches[b*3].GetString();
501 TString name = branches[b*3+1].GetString();
502 TString className = branches[b*3+2].GetString();
503 if(className=="Track") {
504 if(input.Contains("eflow",TString::kIgnoreCase) || name.Contains("eflow",TString::kIgnoreCase)) continue; //no eflow
505 elist = new DelphesBranchElement<TEveTrackList>(name,"track",kBlue);
506 elements.push_back(elist);
507 TEveTrackPropagator *trkProp = elist->GetContainer()->GetPropagator();
508 trkProp->SetMagField(0., 0., -tk_Bz);
509 trkProp->SetMaxR(tk_radius);
510 trkProp->SetMaxZ(tk_length);
511 } else if(className=="Tower") {
512 if(input.Contains("eflow",TString::kIgnoreCase) || name.Contains("eflow",TString::kIgnoreCase)) continue; //no eflow
513 clist = new DelphesBranchElement<DelphesCaloData>(name,"tower",kBlack);
514 clist->GetContainer()->SetEtaBins(etaAxis);
515 clist->GetContainer()->SetPhiBins(phiAxis);
516 elements.push_back(clist);
517 } else if(className=="Jet") {
518 if(input.Contains("GenJetFinder")) {
519 elist = new DelphesBranchElement<TEveElementList>(name,"jet",kCyan);
520 elist->GetContainer()->SetRnrSelf(false);
521 elist->GetContainer()->SetRnrChildren(false);
522 elements.push_back(elist);
523 } else {
524 elements.push_back(new DelphesBranchElement<TEveElementList>(name,"jet",kYellow));
525 }
526 } else if(className=="Electron") {
527 tlist = new DelphesBranchElement<TEveTrackList>(name,"track",kRed);
528 elements.push_back(tlist);
529 TEveTrackPropagator *trkProp = tlist->GetContainer()->GetPropagator();
530 trkProp->SetMagField(0., 0., -tk_Bz);
531 trkProp->SetMaxR(tk_radius);
532 trkProp->SetMaxZ(tk_length);
533 } else if(className=="Photon") {
534 tlist = new DelphesBranchElement<TEveTrackList>(name,"photon",kYellow);
535 elements.push_back(tlist);
536 TEveTrackPropagator *trkProp = tlist->GetContainer()->GetPropagator();
537 trkProp->SetMagField(0., 0., 0.);
538 trkProp->SetMaxR(tk_radius);
539 trkProp->SetMaxZ(tk_length);
540 } else if(className=="Muon") {
541 tlist = new DelphesBranchElement<TEveTrackList>(name,"track",kGreen);
542 elements.push_back(tlist);
543 TEveTrackPropagator *trkProp = tlist->GetContainer()->GetPropagator();
544 trkProp->SetMagField(0., 0., -tk_Bz);
545 trkProp->SetMaxR(mu_radius);
546 trkProp->SetMaxZ(mu_length);
547 } else if(className=="MissingET") {
548 elements.push_back(new DelphesBranchElement<TEveElementList>(name,"vector",kViolet));
549 } else if(className=="GenParticle") {
550 tlist = new DelphesBranchElement<TEveTrackList>(name,"track",kCyan);
551 elements.push_back(tlist);
552 tlist->GetContainer()->SetRnrSelf(false);
553 tlist->GetContainer()->SetRnrChildren(false);
554 TEveTrackPropagator *trkProp = tlist->GetContainer()->GetPropagator();
555 trkProp->SetMagField(0., 0., -tk_Bz);
556 trkProp->SetMaxR(tk_radius);
557 trkProp->SetMaxZ(tk_length);
558 }
559//TODO one possible simplification could be to add the array to the element class.
560 arrays.push_back(gTreeReader->UseBranch(name));
561 }
562}
563
564void delphes_event_display(const char *configFile, const char *inputFile, Delphes3DGeometry& det3D)
565{
566
567 // initialize the application
568 TEveManager::Create(kTRUE, "IV");
569 TGeoManager* geom = gGeoManager;
570
571 // build the detector
572 gRadius = det3D.getTrackerRadius();
573 gTotRadius = det3D.getDetectorRadius();
574 gHalfLength = det3D.getTrackerHalfLength();
575 gBz = det3D.getBField();
576
577 //TODO specific to some classical detector... could use better the det3D
578 TGeoVolume* top = det3D.getDetector(false);
579 geom->SetTopVolume(top);
580 TEveElementList *geometry = new TEveElementList("Geometry");
581 TEveGeoTopNode* trk = new TEveGeoTopNode(gGeoManager, top->FindNode("tracker_1"));
582 trk->SetVisLevel(6);
583 geometry->AddElement(trk);
584 TEveGeoTopNode* calo = new TEveGeoTopNode(gGeoManager, top->FindNode("Calorimeter_barrel_1"));
585 calo->SetVisLevel(3);
586 geometry->AddElement(calo);
587 calo = new TEveGeoTopNode(gGeoManager, top->FindNode("Calorimeter_endcap_1"));
588 calo->SetVisLevel(3);
589 calo->UseNodeTrans();
590 geometry->AddElement(calo);
591 calo = new TEveGeoTopNode(gGeoManager, top->FindNode("Calorimeter_endcap_2"));
592 calo->SetVisLevel(3);
593 calo->UseNodeTrans();
594 geometry->AddElement(calo);
595 TEveGeoTopNode* muon = new TEveGeoTopNode(gGeoManager, top->FindNode("muons_barrel_1"));
596 muon->SetVisLevel(4);
597 geometry->AddElement(muon);
598 muon = new TEveGeoTopNode(gGeoManager, top->FindNode("muons_endcap_1"));
599 muon->SetVisLevel(4);
600 muon->UseNodeTrans();
601 geometry->AddElement(muon);
602 muon = new TEveGeoTopNode(gGeoManager, top->FindNode("muons_endcap_2"));
603 muon->SetVisLevel(4);
604 muon->UseNodeTrans();
605 geometry->AddElement(muon);
606 //gGeoManager->DefaultColors();
607
608 // Create chain of root trees
609 gChain.Add(inputFile);
610
611 // Create object of class ExRootTreeReader
612 printf("*** Opening Delphes data file ***\n");
613 gTreeReader = new ExRootTreeReader(&gChain);
614
615 // prepare data collections
616 readConfig(configFile, det3D, gElements, gArrays);
617
618 // viewers and scenes
619 gDelphesDisplay = new DelphesDisplay;
620 gEve->AddGlobalElement(geometry);
621 gDelphesDisplay->ImportGeomRPhi(geometry);
622 gDelphesDisplay->ImportGeomRhoZ(geometry);
623 // find the first calo data and use that to initialize the calo display
624 for(std::vector<DelphesBranchBase*>::iterator data=gElements.begin();data<gElements.end();++Data) {
625 if(TString((*data)->GetType())=="tower") {
626 TEveCalo3D *calo3d = new TEveCalo3D((*data)->GetContainer());
627 calo3d->SetBarrelRadius(gRadius);
628 calo3d->SetEndCapPos(gHalfLength);
629 gEve->AddGlobalElement(calo3d);
630 gDelphesDisplay->ImportCaloRPhi(calo3d);
631 gDelphesDisplay->ImportCaloRhoZ(calo3d);
632 TEveCaloLego *lego = new TEveCaloLego((*data)->GetContainer());
633 lego->InitMainTrans();
634 lego->RefMainTrans().SetScale(TMath::TwoPi(), TMath::TwoPi(), TMath::Pi());
635 lego->SetAutoRebin(kFALSE);
636 lego->Set2DMode(TEveCaloLego::kValSizeOutline);
637 gDelphesDisplay->ImportCaloLego(lego);
638 break;
639 }
640 }
641 gEve->Redraw3D(kTRUE);
642}
643
644//______________________________________________________________________________
645void load_event()
646{
647 // Load event specified in global event_id.
648 // The contents of previous event are removed.
649
650 //TODO move this to the status bar ???
651 printf("Loading event %d.\n", event_id);
652
653 // clear the previous event
654 gEve->GetViewers()->DeleteAnnotations();
655 for(std::vector<DelphesBranchBase*>::iterator data=gElements.begin();data<gElements.end();++Data) {
656 (*data)->Reset()
657 }
658
659 // read the new event
660 delphes_read();
661
662 // update display
663 TEveElement* top = (TEveElement*)gEve->GetCurrentEvent();
664 gDelphesDisplay->DestroyEventRPhi();
665 gDelphesDisplay->ImportEventRPhi(top);
666 gDelphesDisplay->DestroyEventRhoZ();
667 gDelphesDisplay->ImportEventRhoZ(top);
668 //update_html_summary();
669 gEve->Redraw3D(kFALSE, kTRUE);
670}
671
672void delphes_read()
673{
674
675 // safety
676 if(event_id >= gTreeReader->GetEntries() || event_id<0 ) return;
677
678 // Load selected branches with data from specified event
679 gTreeReader->ReadEntry(event_id);
680
681 // loop over selected branches, and apply the proper recipe to fill the collections.
682 // this is basically to loop on gArrays to fill gElements.
683
684//TODO: one option would be to have templated methods in the element classes. We could simply call "element.fill()"
685 std::vector<TClonesArray*>::iterator data = gArrays.begin();
686 std::vector<DelphesBranchBase*>::iterator element = gElements.begin();
687 std::vector<TClonesArray*>::iterator data_tracks = gArrays.begin();
688 std::vector<DelphesBranchBase*>::iterator element_tracks = gElements.begin();
689 Int_t nTracks = 0;
690 for(; data<gArrays.end() && element<gElements.end(); ++data, ++element) {
691 TString type = (*element)->GetType();
692 // keep the most generic track collection for the end
693 if(type=="track" && (*element)->GetClassName()=="Track" && nTracks=0) {
694 data_tracks = data;
695 element_tracks = element;
696 nTracks = (*data_tracks)->GetEntries();
697 continue;
698 }
699 // branch on the element type
700 // TODO : I understand that we will have to cast the elements.
701 if(type=="tower") delphes_read_towers(*data,*element);
702 else if(type=="track" || type=="photon") delphes_read_tracks(*data,*element);
703 else if(type=="jet") delphes_read_jets(*data,*element);
704 else if(type=="vector") delphes_read_vectors(*data,*element);
705 }
706 // finish whith what we consider to be the main track collection
707 if(nTracks>0) delphes_read_tracks(*data,(*element)->GetContainer());
708}
709
710void delphes_read_towers(TClonesArray* data, DelphesBranchElement<TEveElementList>* element) {
711 DelphesCaloData* container = element->GetContainer();
712 // Loop over all towers
713 TIter itTower(data);
714 Tower *tower;
715 while((tower = (Tower *) itTower.Next()))
716 {
717 container->AddTower(tower->Edges[0], tower->Edges[1], tower->Edges[2], tower->Edges[3]);
718 container->FillSlice(0, tower->Eem);
719 container->FillSlice(1, tower->Ehad);
720 }
721 container->DataChanged();
722}
723
724void delphes_read_tracks(TClonesArray* data, DelphesBranchElement<TEveTrackList>* element) {
725 TEveTrackList* container = element->GetContainer();
726 TString className = element->GetClassName();
727 TIter itTrack(data);
728 Int_t counter = 0;
729 TEveTrack *eveTrack;
730 TEveTrackPropagator *trkProp = container->GetPropagator();
731 if(className=="Track") {
732 // Loop over all tracks
733 Track *track;
734 while((track = (Track *) itTrack.Next())) {
735 TParticle pb(track->PID, 1, 0, 0, 0, 0,
736 track->P4().Px(), track->P4().Py(),
737 track->P4().Pz(), track->P4().E(),
738 track->X, track->Y, track->Z, 0.0);
739
740 eveTrack = new TEveTrack(&pb, counter, trkProp);
741 eveTrack->SetName(Form("%s [%d]", pb.GetName(), counter++));
742 eveTrack->SetStdTitle();
743 eveTrack->SetAttLineAttMarker(gTrackList);
744 container->AddElement(eveTrack);
745 eveTrack->SetLineColor(element->GetColor());
746 eveTrack->MakeTrack();
747 }
748 } else if(className=="Electron") {
749 // Loop over all electrons
750 Electron *electron;
751 while((electron = (Electron *) itElectron.Next())) {
752 TParticle pb(electron->Charge<0?11:-11, 1, 0, 0, 0, 0,
753 electron->P4().Px(), electron->P4().Py(),
754 electron->P4().Pz(), electron->P4().E(),
755 0., 0., 0., 0.);
756
757 eveTrack = new TEveTrack(&pb, counter, trkProp);
758 eveTrack->SetName(Form("%s [%d]", pb.GetName(), counter++));
759 eveTrack->SetStdTitle();
760 eveTrack->SetAttLineAttMarker(gElectronList);
761 container->AddElement(eveTrack);
762 eveTrack->SetLineColor(element->GetColor());
763 eveTrack->MakeTrack();
764 }
765 } else if(className=="Muon") {
766 // Loop over all muons
767 Muon *muon;
768 while((muon = (Muon *) itMuon.Next())) {
769 TParticle pb(muon->Charge<0?13:-13, 1, 0, 0, 0, 0,
770 muon->P4().Px(), muon->P4().Py(),
771 muon->P4().Pz(), muon->P4().E(),
772 0., 0., 0., 0.);
773
774 eveTrack = new TEveTrack(&pb, counter, trkProp);
775 eveTrack->SetName(Form("%s [%d]", pb.GetName(), counter++));
776 eveTrack->SetStdTitle();
777 eveTrack->SetAttLineAttMarker(gMuonList);
778 container->AddElement(eveTrack);
779 eveTrack->SetLineColor(element->GetColor());
780 eveTrack->MakeTrack();
781 }
782 } else if(className=="Photon") {
783 // Loop over all photons
784 Photon *photon;
785 while((photon = (Photon *) itPhoton.Next())) {
786 TParticle pb(22, 1, 0, 0, 0, 0,
787 photon->P4().Px(), photon->P4().Py(),
788 photon->P4().Pz(), photon->P4().E(),
789 0., 0., 0., 0.);
790
791 eveTrack = new TEveTrack(&pb, counter, trkProp);
792 eveTrack->SetName(Form("%s [%d]", pb.GetName(), counter++));
793 eveTrack->SetStdTitle();
794 eveTrack->SetAttLineAttMarker(gPhotonList);
795 container->AddElement(eveTrack);
796 eveTrack->SetLineColor(element->GetColor());
797 eveTrack->MakeTrack();
798 }
799 }
800}
801
802void delphes_read_jets(TClonesArray* data, DelphesBranchElement<TEveElementList>* element) {
803 TEveElementList* container = element->GetContainer();
804 TIter itJet(data);
805 Jet *jet;
806 TEveJetCone *eveJetCone;
807 // Loop over all jets
808 Int_t counter = 0;
809 while((jet = (Jet *) itJet.Next()))
810 {
811 eveJetCone = new TEveJetCone();
812 eveJetCone->SetTitle(Form("jet [%d]: Pt=%f, Eta=%f, \nPhi=%f, M=%f",counter,jet->PT, jet->Eta, jet->Phi, jet->Mass));
813 eveJetCone->SetName(Form("jet [%d]", counter++));
814 eveJetCone->SetMainTransparency(60);
815 eveJetCone->SetLineColor(element->GetColor());
816 eveJetCone->SetFillColor(element->GetColor());
817 eveJetCone->SetCylinder(gRadius - 10, gHalfLength - 10);
818 eveJetCone->SetPickable(kTRUE);
819 eveJetCone->AddEllipticCone(jet->Eta, jet->Phi, jet->DeltaEta, jet->DeltaPhi);
820 container->AddElement(eveJetCone);
821 }
822}
823
824void delphes_read_vectors(TClonesArray* data, DelphesBranchElement<TEveElementList>* element) {
825 TEveElementList* container = element->GetContainer();
826 TIter itMet(data);
827 MissingET *MET;
828 TEveArrow *eveMet;
829 // Missing Et
830 Double_t maxPt = 50.;
831 // TODO to be changed as we don't have access to maxPt anymore. MET scale could be a general parameter set in GUI
832 while((MET = (MissingET*) itMet.Next())) {
833 eveMet = new TEveArrow((gRadius * MET->MET/maxPt)*cos(MET->Phi), (gRadius * MET->MET/maxPt)*sin(MET->Phi), 0., 0., 0., 0.);
834 eveMet->SetMainColor(element->GetColor());
835 eveMet->SetTubeR(0.04);
836 eveMet->SetConeR(0.08);
837 eveMet->SetConeL(0.10);
838 eveMet->SetPickable(kTRUE);
839 eveMet->SetName("Missing Et");
840 eveMet->SetTitle(Form("Missing Et (%.1f GeV)",MET->MET));
841 container->AddElement(eveMet);
842 }
843}
844
845/******************************************************************************/
846// GUI
847/******************************************************************************/
848
849//______________________________________________________________________________
850//
851// EvNavHandler class is needed to connect GUI signals.
852
853class EvNavHandler
854{
855public:
856 void Fwd()
857 {
858 if (event_id < gTreeReader->GetEntries() - 1) {
859 ++event_id;
860 load_event();
861 } else {
862 printf("Already at last event.\n");
863 }
864 }
865 void Bck()
866 {
867 if (event_id > 0) {
868 --event_id;
869 load_event();
870 } else {
871 printf("Already at first event.\n");
872 }
873 }
874};
875
876//______________________________________________________________________________
877void make_gui()
878{
879 // Create minimal GUI for event navigation.
880 // TODO: better GUI could be made based on the ch15 of the manual (Writing a GUI)
881
882 // add a tab on the left
883 TEveBrowser* browser = gEve->GetBrowser();
884 browser->StartEmbedding(TRootBrowser::kLeft);
885
886 // set the main title
887 TGMainFrame* frmMain = new TGMainFrame(gClient->GetRoot(), 1000, 600);
888 frmMain->SetWindowName("Delphes Event Display");
889 frmMain->SetCleanup(kDeepCleanup);
890
891 // build the navigation menu
892 TGHorizontalFrame* hf = new TGHorizontalFrame(frmMain);
893 {
894 TString icondir;
895 if(gSystem->Getenv("ROOTSYS"))
896 icondir = Form("%s/icons/", gSystem->Getenv("ROOTSYS"));
897 if(!gSystem->OpenDirectory(icondir))
898 icondir = Form("%s/icons/", (const char*)gSystem->GetFromPipe("root-config --etcdir") );
899 TGPictureButton* b = 0;
900 EvNavHandler *fh = new EvNavHandler;
901
902 b = new TGPictureButton(hf, gClient->GetPicture(icondir+"GoBack.gif"));
903 hf->AddFrame(b);
904 b->Connect("Clicked()", "EvNavHandler", fh, "Bck()");
905
906 b = new TGPictureButton(hf, gClient->GetPicture(icondir+"GoForward.gif"));
907 hf->AddFrame(b);
908 b->Connect("Clicked()", "EvNavHandler", fh, "Fwd()");
909 }
910 frmMain->AddFrame(hf);
911 frmMain->MapSubwindows();
912 frmMain->Resize();
913 frmMain->MapWindow();
914 browser->StopEmbedding();
915 browser->SetTabTitle("Event Control", 0);
916}
917
918/******************************************************************************/
919// MAIN
920/******************************************************************************/
921
922void geometry(const char* filename = "delphes_card_CMS.tcl", const char* ParticlePropagator="ParticlePropagator",
923 const char* TrackingEfficiency="ChargedHadronTrackingEfficiency",
924 const char* MuonEfficiency="MuonEfficiency",
925 const char* Calorimeters="Calorimeter")
926{
927
928 // load the libraries
929 gSystem->Load("libGeom");
930 gSystem->Load("../libDelphes");
931 gSystem->Load("../libDelphesDisplay");
932
933 // create the detector representation
934 Delphes3DGeometry det3D(new TGeoManager("delphes", "Delphes geometry"));
935 det3D.readFile(filename, ParticlePropagator, TrackingEfficiency, MuonEfficiency, Calorimeters);
936
937 // create the application items
938 delphes_event_display("delphes_card_CMS.tcl", "../delphes_output.root", det3D); //TODO root file as input cfg
939 make_gui();
940 load_event();
941 gEve->Redraw3D(kTRUE); // Reset camera after the first event has been shown.
942
943 // EClipType not exported to CINT (see TGLUtil.h):
944 // 0 - no clip, 1 - clip plane, 2 - clip box
945 TGLViewer *v = gEve->GetDefaultGLViewer();
946 //Double_t plane[4] = { 0., 1., 0., 0. };
947 //v->GetClipSet()->SetClipState(1,plane);
948 //v->GetClipSet()->SetClipType(1);
949 //v->ColorSet().Background().SetColor(kMagenta+4);
950 //v->SetGuideState(TGLUtil::kAxesEdge, kTRUE, kFALSE, 0);
951 v->RefreshPadEditor(v);
952 v->CurrentCamera().RotateRad(-1.2, 0.5);
953 v->DoDraw();
954
955}
956
Note: See TracBrowser for help on using the repository browser.