Fork me on GitHub

source: git/display/DelphesEventDisplay.cc@ 84dd1c8

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

Debugging ongoing - signs of corruption

  • Property mode set to 100644
File size: 18.8 KB
Line 
1#include <cassert>
2#include <iostream>
3#include <utility>
4#include <algorithm>
5#include "TGeoManager.h"
6#include "TGeoVolume.h"
7#include "external/ExRootAnalysis/ExRootConfReader.h"
8#include "external/ExRootAnalysis/ExRootTreeReader.h"
9#include "display/DelphesCaloData.h"
10#include "display/DelphesBranchElement.h"
11#include "display/Delphes3DGeometry.h"
12#include "display/DelphesEventDisplay.h"
13#include "classes/DelphesClasses.h"
14#include "TEveElement.h"
15#include "TEveJetCone.h"
16#include "TEveTrack.h"
17#include "TEveTrackPropagator.h"
18#include "TEveCalo.h"
19#include "TEveManager.h"
20#include "TEveGeoNode.h"
21#include "TEveTrans.h"
22#include "TEveViewer.h"
23#include "TEveBrowser.h"
24#include "TEveArrow.h"
25#include "TMath.h"
26#include "TSystem.h"
27#include "TRootBrowser.h"
28#include "TGButton.h"
29#include "TClonesArray.h"
30
31DelphesEventDisplay::DelphesEventDisplay()
32{
33 event_id = 0;
34 tkRadius_ = 1.29;
35 totRadius_ = 2.0;
36 tkHalfLength_ = 3.0;
37 bz_ = 3.8;
38 chain_ = new TChain("Delphes");
39 gTreeReader = 0;
40 gDelphesDisplay = 0;
41}
42
43DelphesEventDisplay::~DelphesEventDisplay()
44{
45 delete chain_;
46}
47
48DelphesEventDisplay::DelphesEventDisplay(const char *configFile, const char *inputFile, Delphes3DGeometry& det3D)
49{
50 event_id = 0;
51 tkRadius_ = 1.29;
52 totRadius_ = 2.0;
53 tkHalfLength_ = 3.0;
54 bz_ = 3.8;
55 chain_ = new TChain("Delphes");
56 gTreeReader = 0;
57 gDelphesDisplay = 0;
58
59 // initialize the application
60 TEveManager::Create(kTRUE, "IV");
61 TGeoManager* geom = gGeoManager;
62
63 // build the detector
64 tkRadius_ = det3D.getTrackerRadius();
65 totRadius_ = det3D.getDetectorRadius();
66 tkHalfLength_ = det3D.getTrackerHalfLength();
67 bz_ = det3D.getBField();
68
69 //TODO specific to some classical detector... could use better the det3D
70 TGeoVolume* top = det3D.getDetector(false);
71 geom->SetTopVolume(top);
72 TEveElementList *geometry = new TEveElementList("Geometry");
73 TEveGeoTopNode* trk = new TEveGeoTopNode(gGeoManager, top->FindNode("tracker_1"));
74 trk->SetVisLevel(6);
75 geometry->AddElement(trk);
76 TEveGeoTopNode* calo = new TEveGeoTopNode(gGeoManager, top->FindNode("Calorimeter_barrel_1"));
77 calo->SetVisLevel(3);
78 geometry->AddElement(calo);
79 calo = new TEveGeoTopNode(gGeoManager, top->FindNode("Calorimeter_endcap_1"));
80 calo->SetVisLevel(3);
81 calo->UseNodeTrans();
82 geometry->AddElement(calo);
83 calo = new TEveGeoTopNode(gGeoManager, top->FindNode("Calorimeter_endcap_2"));
84 calo->SetVisLevel(3);
85 calo->UseNodeTrans();
86 geometry->AddElement(calo);
87 TEveGeoTopNode* muon = new TEveGeoTopNode(gGeoManager, top->FindNode("muons_barrel_1"));
88 muon->SetVisLevel(4);
89 geometry->AddElement(muon);
90 muon = new TEveGeoTopNode(gGeoManager, top->FindNode("muons_endcap_1"));
91 muon->SetVisLevel(4);
92 muon->UseNodeTrans();
93 geometry->AddElement(muon);
94 muon = new TEveGeoTopNode(gGeoManager, top->FindNode("muons_endcap_2"));
95 muon->SetVisLevel(4);
96 muon->UseNodeTrans();
97 geometry->AddElement(muon);
98 //gGeoManager->DefaultColors();
99
100 // Create chain of root trees
101 chain_->Add(inputFile);
102
103 // Create object of class ExRootTreeReader
104 printf("*** Opening Delphes data file ***\n");
105 gTreeReader = new ExRootTreeReader(chain_);
106
107 // prepare data collections
108 readConfig(configFile, det3D, gElements, gArrays);
109 for(std::vector<DelphesBranchBase*>::iterator element = gElements.begin(); element<gElements.end(); ++element) {
110 DelphesBranchElement<TEveTrackList>* item_v1 = dynamic_cast<DelphesBranchElement<TEveTrackList>*>(*element);
111 DelphesBranchElement<DelphesCaloData>* item_v2 = dynamic_cast<DelphesBranchElement<DelphesCaloData>*>(*element);
112 DelphesBranchElement<TEveElementList>* item_v3 = dynamic_cast<DelphesBranchElement<TEveElementList>*>(*element);
113 if(item_v1) gEve->AddElement(item_v1->GetContainer());
114 if(item_v2) gEve->AddElement(item_v2->GetContainer());
115 if(item_v3) gEve->AddElement(item_v3->GetContainer());
116 }
117
118 // viewers and scenes
119 gDelphesDisplay = new DelphesDisplay;
120 gEve->AddGlobalElement(geometry);
121 gDelphesDisplay->ImportGeomRPhi(geometry);
122 gDelphesDisplay->ImportGeomRhoZ(geometry);
123 // find the first calo data and use that to initialize the calo display
124 for(std::vector<DelphesBranchBase*>::iterator data=gElements.begin();data<gElements.end();++data) {
125 if(TString((*data)->GetType())=="tower") {
126 DelphesCaloData* container = dynamic_cast<DelphesBranchElement<DelphesCaloData>*>((*data))->GetContainer();
127 assert(container);
128 TEveCalo3D *calo3d = new TEveCalo3D(container);
129 calo3d->SetBarrelRadius(tkRadius_);
130 calo3d->SetEndCapPos(tkHalfLength_);
131 gEve->AddGlobalElement(calo3d);
132 gDelphesDisplay->ImportCaloRPhi(calo3d);
133 gDelphesDisplay->ImportCaloRhoZ(calo3d);
134 TEveCaloLego *lego = new TEveCaloLego(container);
135 lego->InitMainTrans();
136 lego->RefMainTrans().SetScale(TMath::TwoPi(), TMath::TwoPi(), TMath::Pi());
137 lego->SetAutoRebin(kFALSE);
138 lego->Set2DMode(TEveCaloLego::kValSizeOutline);
139 gDelphesDisplay->ImportCaloLego(lego);
140 break;
141 }
142 }
143
144 //make_gui();
145 //load_event();
146 gEve->Redraw3D(kTRUE);
147
148}
149// function that parses the config to extract the branches of interest and prepare containers
150void DelphesEventDisplay::readConfig(const char *configFile, Delphes3DGeometry& det3D, std::vector<DelphesBranchBase*>& elements, std::vector<TClonesArray*>& arrays) {
151 ExRootConfReader *confReader = new ExRootConfReader;
152 confReader->ReadFile(configFile);
153 Double_t tk_radius = det3D.getTrackerRadius();
154 Double_t tk_length = det3D.getTrackerHalfLength();
155 Double_t tk_Bz = det3D.getBField();
156 Double_t mu_radius = det3D.getDetectorRadius();
157 Double_t mu_length = det3D.getDetectorHalfLength();
158 TAxis* etaAxis = det3D.getCaloAxes().first;
159 TAxis* phiAxis = det3D.getCaloAxes().second;
160 ExRootConfParam branches = confReader->GetParam("TreeWriter::Branch");
161 Int_t nBranches = branches.GetSize()/3;
162 DelphesBranchElement<TEveTrackList>* tlist;
163 DelphesBranchElement<DelphesCaloData>* clist;
164 DelphesBranchElement<TEveElementList>* elist;
165 for(Int_t b = 0; b<nBranches; ++b) {
166 TString input = branches[b*3].GetString();
167 TString name = branches[b*3+1].GetString();
168 TString className = branches[b*3+2].GetString();
169 if(className=="Track") {
170 if(input.Contains("eflow",TString::kIgnoreCase) || name.Contains("eflow",TString::kIgnoreCase)) continue; //no eflow
171 tlist = new DelphesBranchElement<TEveTrackList>(name,"track",kBlue);
172 elements.push_back(tlist);
173 TEveTrackPropagator *trkProp = tlist->GetContainer()->GetPropagator();
174 trkProp->SetMagField(0., 0., -tk_Bz);
175 trkProp->SetMaxR(tk_radius);
176 trkProp->SetMaxZ(tk_length);
177 } else if(className=="Tower") {
178 if(input.Contains("eflow",TString::kIgnoreCase) || name.Contains("eflow",TString::kIgnoreCase)) continue; //no eflow
179 clist = new DelphesBranchElement<DelphesCaloData>(name,"tower",kBlack);
180 clist->GetContainer()->SetEtaBins(etaAxis);
181 clist->GetContainer()->SetPhiBins(phiAxis);
182 elements.push_back(clist);
183 } else if(className=="Jet") {
184 if(input.Contains("GenJetFinder")) {
185 elist = new DelphesBranchElement<TEveElementList>(name,"jet",kCyan);
186 elist->GetContainer()->SetRnrSelf(false);
187 elist->GetContainer()->SetRnrChildren(false);
188 elements.push_back(elist);
189 } else {
190 elements.push_back(new DelphesBranchElement<TEveElementList>(name,"jet",kYellow));
191 }
192 } else if(className=="Electron") {
193 tlist = new DelphesBranchElement<TEveTrackList>(name,"track",kRed);
194 elements.push_back(tlist);
195 TEveTrackPropagator *trkProp = tlist->GetContainer()->GetPropagator();
196 trkProp->SetMagField(0., 0., -tk_Bz);
197 trkProp->SetMaxR(tk_radius);
198 trkProp->SetMaxZ(tk_length);
199 } else if(className=="Photon") {
200 tlist = new DelphesBranchElement<TEveTrackList>(name,"photon",kYellow);
201 elements.push_back(tlist);
202 TEveTrackPropagator *trkProp = tlist->GetContainer()->GetPropagator();
203 trkProp->SetMagField(0., 0., 0.);
204 trkProp->SetMaxR(tk_radius);
205 trkProp->SetMaxZ(tk_length);
206 } else if(className=="Muon") {
207 tlist = new DelphesBranchElement<TEveTrackList>(name,"track",kGreen);
208 elements.push_back(tlist);
209 TEveTrackPropagator *trkProp = tlist->GetContainer()->GetPropagator();
210 trkProp->SetMagField(0., 0., -tk_Bz);
211 trkProp->SetMaxR(mu_radius);
212 trkProp->SetMaxZ(mu_length);
213 } else if(className=="MissingET") {
214 elements.push_back(new DelphesBranchElement<TEveElementList>(name,"vector",kViolet));
215 } else if(className=="GenParticle") {
216 tlist = new DelphesBranchElement<TEveTrackList>(name,"track",kCyan);
217 elements.push_back(tlist);
218 tlist->GetContainer()->SetRnrSelf(false);
219 tlist->GetContainer()->SetRnrChildren(false);
220 TEveTrackPropagator *trkProp = tlist->GetContainer()->GetPropagator();
221 trkProp->SetMagField(0., 0., -tk_Bz);
222 trkProp->SetMaxR(tk_radius);
223 trkProp->SetMaxZ(tk_length);
224 }
225//TODO one possible simplification could be to add the array to the element class.
226 arrays.push_back(gTreeReader->UseBranch(name));
227 }
228}
229
230
231//______________________________________________________________________________
232void DelphesEventDisplay::load_event()
233{
234 // Load event specified in global event_id.
235 // The contents of previous event are removed.
236
237 //TODO move this to the status bar ???
238 printf("Loading event %d.\n", event_id);
239
240
241 // clear the previous event
242 gEve->GetViewers()->DeleteAnnotations();
243 for(std::vector<DelphesBranchBase*>::iterator data=gElements.begin();data<gElements.end();++data) {
244 (*data)->Reset();
245 }
246
247 // read the new event
248 delphes_read();
249
250//TODO: it blocks somewhere below....
251//TODO: also the event content has one weird entry "TEveCalData" -> corruption????
252//other observation: projections appear in the 3D view ! -> seems to indicate that we have a common problem there.
253//should check the content of the elements vector when filling them
254 // update display
255 TEveElement* top = (TEveElement*)gEve->GetCurrentEvent();
256 gDelphesDisplay->DestroyEventRPhi();
257 gDelphesDisplay->ImportEventRPhi(top);
258 gDelphesDisplay->DestroyEventRhoZ();
259 gDelphesDisplay->ImportEventRhoZ(top);
260 //update_html_summary();
261 gEve->Redraw3D(kFALSE, kTRUE);
262}
263
264void DelphesEventDisplay::delphes_read()
265{
266
267 // safety
268 if(event_id >= gTreeReader->GetEntries() || event_id<0 ) return;
269
270 // Load selected branches with data from specified event
271 gTreeReader->ReadEntry(event_id);
272
273 // loop over selected branches, and apply the proper recipe to fill the collections.
274 // this is basically to loop on gArrays to fill gElements.
275
276//TODO: one option would be to have templated methods in the element classes. We could simply call "element.fill()"
277 std::vector<TClonesArray*>::iterator data = gArrays.begin();
278 std::vector<DelphesBranchBase*>::iterator element = gElements.begin();
279 std::vector<TClonesArray*>::iterator data_tracks = gArrays.begin();
280 std::vector<DelphesBranchBase*>::iterator element_tracks = gElements.begin();
281 Int_t nTracks = 0;
282 for(; data<gArrays.end() && element<gElements.end(); ++data, ++element) {
283 TString type = (*element)->GetType();
284 // keep the most generic track collection for the end
285 if(type=="track" && TString((*element)->GetClassName())=="Track" && nTracks==0) {
286 data_tracks = data;
287 element_tracks = element;
288 nTracks = (*data_tracks)->GetEntries();
289 continue;
290 }
291 // branch on the element type
292 if(type=="tower") delphes_read_towers(*data,*element);
293 else if(type=="track" || type=="photon") delphes_read_tracks(*data,*element);
294 else if(type=="jet") delphes_read_jets(*data,*element);
295 else if(type=="vector") delphes_read_vectors(*data,*element);
296 }
297 // finish whith what we consider to be the main track collection
298 if(nTracks>0) delphes_read_tracks(*data,*element);
299}
300
301void DelphesEventDisplay::delphes_read_towers(TClonesArray* data, DelphesBranchBase* element) {
302 DelphesCaloData* container = dynamic_cast<DelphesBranchElement<DelphesCaloData>*>(element)->GetContainer();
303 assert(container);
304 // Loop over all towers
305 TIter itTower(data);
306 Tower *tower;
307 while((tower = (Tower *) itTower.Next()))
308 {
309 container->AddTower(tower->Edges[0], tower->Edges[1], tower->Edges[2], tower->Edges[3]);
310 container->FillSlice(0, tower->Eem);
311 container->FillSlice(1, tower->Ehad);
312 }
313 container->DataChanged();
314}
315
316void DelphesEventDisplay::delphes_read_tracks(TClonesArray* data, DelphesBranchBase* element) {
317 TEveTrackList* container = dynamic_cast<DelphesBranchElement<TEveTrackList>*>(element)->GetContainer();
318 assert(container);
319 TString className = element->GetClassName();
320 TIter itTrack(data);
321 Int_t counter = 0;
322 TEveTrack *eveTrack;
323 TEveTrackPropagator *trkProp = container->GetPropagator();
324 if(className=="Track") {
325 // Loop over all tracks
326 Track *track;
327 while((track = (Track *) itTrack.Next())) {
328 TParticle pb(track->PID, 1, 0, 0, 0, 0,
329 track->P4().Px(), track->P4().Py(),
330 track->P4().Pz(), track->P4().E(),
331 track->X, track->Y, track->Z, 0.0);
332
333 eveTrack = new TEveTrack(&pb, counter, trkProp);
334 eveTrack->SetName(Form("%s [%d]", pb.GetName(), counter++));
335 eveTrack->SetStdTitle();
336 eveTrack->SetAttLineAttMarker(container);
337 container->AddElement(eveTrack);
338 eveTrack->SetLineColor(element->GetColor());
339 eveTrack->MakeTrack();
340 }
341 } else if(className=="Electron") {
342 // Loop over all electrons
343 Electron *electron;
344 while((electron = (Electron *) itTrack.Next())) {
345 TParticle pb(electron->Charge<0?11:-11, 1, 0, 0, 0, 0,
346 electron->P4().Px(), electron->P4().Py(),
347 electron->P4().Pz(), electron->P4().E(),
348 0., 0., 0., 0.);
349
350 eveTrack = new TEveTrack(&pb, counter, trkProp);
351 eveTrack->SetName(Form("%s [%d]", pb.GetName(), counter++));
352 eveTrack->SetStdTitle();
353 eveTrack->SetAttLineAttMarker(container);
354 container->AddElement(eveTrack);
355 eveTrack->SetLineColor(element->GetColor());
356 eveTrack->MakeTrack();
357 }
358 } else if(className=="Muon") {
359 // Loop over all muons
360 Muon *muon;
361 while((muon = (Muon *) itTrack.Next())) {
362 TParticle pb(muon->Charge<0?13:-13, 1, 0, 0, 0, 0,
363 muon->P4().Px(), muon->P4().Py(),
364 muon->P4().Pz(), muon->P4().E(),
365 0., 0., 0., 0.);
366
367 eveTrack = new TEveTrack(&pb, counter, trkProp);
368 eveTrack->SetName(Form("%s [%d]", pb.GetName(), counter++));
369 eveTrack->SetStdTitle();
370 eveTrack->SetAttLineAttMarker(container);
371 container->AddElement(eveTrack);
372 eveTrack->SetLineColor(element->GetColor());
373 eveTrack->MakeTrack();
374 }
375 } else if(className=="Photon") {
376 // Loop over all photons
377 Photon *photon;
378 while((photon = (Photon *) itTrack.Next())) {
379 TParticle pb(22, 1, 0, 0, 0, 0,
380 photon->P4().Px(), photon->P4().Py(),
381 photon->P4().Pz(), photon->P4().E(),
382 0., 0., 0., 0.);
383
384 eveTrack = new TEveTrack(&pb, counter, trkProp);
385 eveTrack->SetName(Form("%s [%d]", pb.GetName(), counter++));
386 eveTrack->SetStdTitle();
387 eveTrack->SetAttLineAttMarker(container);
388 container->AddElement(eveTrack);
389 eveTrack->SetLineColor(element->GetColor());
390 eveTrack->MakeTrack();
391 }
392 }
393}
394
395void DelphesEventDisplay::delphes_read_jets(TClonesArray* data, DelphesBranchBase* element) {
396 TEveElementList* container = dynamic_cast<DelphesBranchElement<TEveElementList>*>(element)->GetContainer();
397 assert(container);
398 TIter itJet(data);
399 Jet *jet;
400 TEveJetCone *eveJetCone;
401 // Loop over all jets
402 Int_t counter = 0;
403 while((jet = (Jet *) itJet.Next()))
404 {
405 eveJetCone = new TEveJetCone();
406 eveJetCone->SetTitle(Form("jet [%d]: Pt=%f, Eta=%f, \nPhi=%f, M=%f",counter,jet->PT, jet->Eta, jet->Phi, jet->Mass));
407 eveJetCone->SetName(Form("jet [%d]", counter++));
408 eveJetCone->SetMainTransparency(60);
409 eveJetCone->SetLineColor(element->GetColor());
410 eveJetCone->SetFillColor(element->GetColor());
411 eveJetCone->SetCylinder(tkRadius_ - 10, tkHalfLength_ - 10);
412 eveJetCone->SetPickable(kTRUE);
413 eveJetCone->AddEllipticCone(jet->Eta, jet->Phi, jet->DeltaEta, jet->DeltaPhi);
414 container->AddElement(eveJetCone);
415 }
416}
417
418void DelphesEventDisplay::delphes_read_vectors(TClonesArray* data, DelphesBranchBase* element) {
419 TEveElementList* container = dynamic_cast<DelphesBranchElement<TEveElementList>*>(element)->GetContainer();
420 assert(container);
421 TIter itMet(data);
422 MissingET *MET;
423 TEveArrow *eveMet;
424 // Missing Et
425 Double_t maxPt = 50.;
426 // TODO to be changed as we don't have access to maxPt anymore. MET scale could be a general parameter set in GUI
427 while((MET = (MissingET*) itMet.Next())) {
428 eveMet = new TEveArrow((tkRadius_ * MET->MET/maxPt)*cos(MET->Phi), (tkRadius_ * MET->MET/maxPt)*sin(MET->Phi), 0., 0., 0., 0.);
429 eveMet->SetMainColor(element->GetColor());
430 eveMet->SetTubeR(0.04);
431 eveMet->SetConeR(0.08);
432 eveMet->SetConeL(0.10);
433 eveMet->SetPickable(kTRUE);
434 eveMet->SetName("Missing Et");
435 eveMet->SetTitle(Form("Missing Et (%.1f GeV)",MET->MET));
436 container->AddElement(eveMet);
437 }
438}
439
440/******************************************************************************/
441// GUI
442/******************************************************************************/
443
444void DelphesEventDisplay::make_gui()
445{
446 // Create minimal GUI for event navigation.
447 // TODO: better GUI could be made based on the ch15 of the manual (Writing a GUI)
448
449 // add a tab on the left
450 TEveBrowser* browser = gEve->GetBrowser();
451 browser->StartEmbedding(TRootBrowser::kLeft);
452
453 // set the main title
454 TGMainFrame* frmMain = new TGMainFrame(gClient->GetRoot(), 1000, 600);
455 frmMain->SetWindowName("Delphes Event Display");
456 frmMain->SetCleanup(kDeepCleanup);
457
458 // build the navigation menu
459 TGHorizontalFrame* hf = new TGHorizontalFrame(frmMain);
460 {
461 TString icondir;
462 if(gSystem->Getenv("ROOTSYS"))
463 icondir = Form("%s/icons/", gSystem->Getenv("ROOTSYS"));
464 if(!gSystem->OpenDirectory(icondir))
465 icondir = Form("%s/icons/", (const char*)gSystem->GetFromPipe("root-config --etcdir") );
466 TGPictureButton* b = 0;
467 EvNavHandler *fh = new EvNavHandler(this);
468
469 b = new TGPictureButton(hf, gClient->GetPicture(icondir+"GoBack.gif"));
470 hf->AddFrame(b);
471 b->Connect("Clicked()", "EvNavHandler", fh, "Bck()");
472
473 b = new TGPictureButton(hf, gClient->GetPicture(icondir+"GoForward.gif"));
474 hf->AddFrame(b);
475 b->Connect("Clicked()", "EvNavHandler", fh, "Fwd()");
476 }
477 frmMain->AddFrame(hf);
478 frmMain->MapSubwindows();
479 frmMain->Resize();
480 frmMain->MapWindow();
481 browser->StopEmbedding();
482 browser->SetTabTitle("Event Control", 0);
483}
484
Note: See TracBrowser for help on using the repository browser.