Fork me on GitHub

source: git/examples/geometry.C@ 4d999a57

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

Still broken, but better.

debugging may require to transfer chuncks of the code to compiled code.
I have in mind the geometry and the code that prepares the GUI.
Right now I got a segfault + funny cint errors.

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