Fork me on GitHub

source: git/display/DelphesEventDisplay.cc@ a3b2495

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

Display of automatically selected elements working

The original features are back, but with a much more flexible code that
selects automatically the list of branches from the tcl file. The main
issue left is the GUI (handler not found), and then some more
cleaning. One idea is to further exploit templating to avoid ifs when
leading events.

  • Property mode set to 100644
File size: 19.1 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<TEveElementList>* item_v2 = dynamic_cast<DelphesBranchElement<TEveElementList>*>(*element);
112 if(item_v1) gEve->AddElement(item_v1->GetContainer());
113 if(item_v2) gEve->AddElement(item_v2->GetContainer());
114 }
115
116 // viewers and scenes
117 gDelphesDisplay = new DelphesDisplay;
118 gEve->AddGlobalElement(geometry);
119 gDelphesDisplay->ImportGeomRPhi(geometry);
120 gDelphesDisplay->ImportGeomRhoZ(geometry);
121 // find the first calo data and use that to initialize the calo display
122 for(std::vector<DelphesBranchBase*>::iterator data=gElements.begin();data<gElements.end();++data) {
123 if(TString((*data)->GetType())=="tower") { // we could also use GetClassName()=="DelphesCaloData"
124 DelphesCaloData* container = dynamic_cast<DelphesBranchElement<DelphesCaloData>*>((*data))->GetContainer();
125 assert(container);
126 TEveCalo3D *calo3d = new TEveCalo3D(container);
127 calo3d->SetBarrelRadius(tkRadius_);
128 calo3d->SetEndCapPos(tkHalfLength_);
129 gEve->AddGlobalElement(calo3d);
130 gDelphesDisplay->ImportCaloRPhi(calo3d);
131 gDelphesDisplay->ImportCaloRhoZ(calo3d);
132 TEveCaloLego *lego = new TEveCaloLego(container);
133 lego->InitMainTrans();
134 lego->RefMainTrans().SetScale(TMath::TwoPi(), TMath::TwoPi(), TMath::Pi());
135 lego->SetAutoRebin(kFALSE);
136 lego->Set2DMode(TEveCaloLego::kValSizeOutline);
137 gDelphesDisplay->ImportCaloLego(lego);
138 break;
139 }
140 }
141
142 //make_gui(); //TODO put back!
143 //load_event(); //TODO put back!
144 gEve->Redraw3D(kTRUE);
145
146}
147// function that parses the config to extract the branches of interest and prepare containers
148void DelphesEventDisplay::readConfig(const char *configFile, Delphes3DGeometry& det3D, std::vector<DelphesBranchBase*>& elements, std::vector<TClonesArray*>& arrays) {
149 ExRootConfReader *confReader = new ExRootConfReader;
150 confReader->ReadFile(configFile);
151 Double_t tk_radius = det3D.getTrackerRadius();
152 Double_t tk_length = det3D.getTrackerHalfLength();
153 Double_t tk_Bz = det3D.getBField();
154 Double_t mu_radius = det3D.getDetectorRadius();
155 Double_t mu_length = det3D.getDetectorHalfLength();
156 TAxis* etaAxis = det3D.getCaloAxes().first;
157 TAxis* phiAxis = det3D.getCaloAxes().second;
158 ExRootConfParam branches = confReader->GetParam("TreeWriter::Branch");
159 Int_t nBranches = branches.GetSize()/3;
160 DelphesBranchElement<TEveTrackList>* tlist;
161 DelphesBranchElement<DelphesCaloData>* clist;
162 DelphesBranchElement<TEveElementList>* elist;
163 // first loop with all but tracks
164 for(Int_t b = 0; b<nBranches; ++b) {
165 TString input = branches[b*3].GetString();
166 TString name = branches[b*3+1].GetString();
167 TString className = branches[b*3+2].GetString();
168 if(className=="Tower") {
169 if(input.Contains("eflow",TString::kIgnoreCase) || name.Contains("eflow",TString::kIgnoreCase)) continue; //no eflow
170 clist = new DelphesBranchElement<DelphesCaloData>(name,"tower",kBlack);
171 clist->GetContainer()->SetEtaBins(etaAxis);
172 clist->GetContainer()->SetPhiBins(phiAxis);
173 elements.push_back(clist);
174 } else if(className=="Jet") {
175 if(input.Contains("GenJetFinder")) {
176 elist = new DelphesBranchElement<TEveElementList>(name,"jet",kCyan);
177 elist->GetContainer()->SetRnrSelf(false);
178 elist->GetContainer()->SetRnrChildren(false);
179 elements.push_back(elist);
180 } else {
181 elements.push_back(new DelphesBranchElement<TEveElementList>(name,"jet",kYellow));
182 }
183 } else if(className=="Electron") {
184 tlist = new DelphesBranchElement<TEveTrackList>(name,"electron",kRed);
185 elements.push_back(tlist);
186 TEveTrackPropagator *trkProp = tlist->GetContainer()->GetPropagator();
187 trkProp->SetMagField(0., 0., -tk_Bz);
188 trkProp->SetMaxR(tk_radius);
189 trkProp->SetMaxZ(tk_length);
190 } else if(className=="Photon") {
191 tlist = new DelphesBranchElement<TEveTrackList>(name,"photon",kYellow);
192 elements.push_back(tlist);
193 TEveTrackPropagator *trkProp = tlist->GetContainer()->GetPropagator();
194 trkProp->SetMagField(0., 0., 0.);
195 trkProp->SetMaxR(tk_radius);
196 trkProp->SetMaxZ(tk_length);
197 } else if(className=="Muon") {
198 tlist = new DelphesBranchElement<TEveTrackList>(name,"muon",kGreen);
199 elements.push_back(tlist);
200 TEveTrackPropagator *trkProp = tlist->GetContainer()->GetPropagator();
201 trkProp->SetMagField(0., 0., -tk_Bz);
202 trkProp->SetMaxR(mu_radius);
203 trkProp->SetMaxZ(mu_length);
204 } else if(className=="MissingET") {
205 elements.push_back(new DelphesBranchElement<TEveElementList>(name,"vector",kViolet));
206 } else if(className=="GenParticle") {
207 tlist = new DelphesBranchElement<TEveTrackList>(name,"genparticle",kCyan);
208 elements.push_back(tlist);
209 tlist->GetContainer()->SetRnrSelf(false);
210 tlist->GetContainer()->SetRnrChildren(false);
211 TEveTrackPropagator *trkProp = tlist->GetContainer()->GetPropagator();
212 trkProp->SetMagField(0., 0., -tk_Bz);
213 trkProp->SetMaxR(tk_radius);
214 trkProp->SetMaxZ(tk_length);
215 } else {
216 continue;
217 }
218//TODO one possible simplification could be to add the array to the element class.
219 arrays.push_back(gTreeReader->UseBranch(name));
220 }
221 // second loop for tracks
222 for(Int_t b = 0; b<nBranches; ++b) {
223 TString input = branches[b*3].GetString();
224 TString name = branches[b*3+1].GetString();
225 TString className = branches[b*3+2].GetString();
226 if(className=="Track") {
227 if(input.Contains("eflow",TString::kIgnoreCase) || name.Contains("eflow",TString::kIgnoreCase)) continue; //no eflow
228 tlist = new DelphesBranchElement<TEveTrackList>(name,"track",kBlue);
229 elements.push_back(tlist);
230 TEveTrackPropagator *trkProp = tlist->GetContainer()->GetPropagator();
231 trkProp->SetMagField(0., 0., -tk_Bz);
232 trkProp->SetMaxR(tk_radius);
233 trkProp->SetMaxZ(tk_length);
234 arrays.push_back(gTreeReader->UseBranch(name));
235 }
236 }
237}
238
239
240//______________________________________________________________________________
241void DelphesEventDisplay::load_event()
242{
243 // Load event specified in global event_id.
244 // The contents of previous event are removed.
245
246 //TODO move this to the status bar ???
247 printf("Loading event %d.\n", event_id);
248
249
250 // clear the previous event
251 gEve->GetViewers()->DeleteAnnotations();
252 for(std::vector<DelphesBranchBase*>::iterator data=gElements.begin();data<gElements.end();++data) {
253 (*data)->Reset();
254 }
255
256 // read the new event
257 delphes_read();
258
259 // update display
260 TEveElement* top = (TEveElement*)gEve->GetCurrentEvent();
261 gDelphesDisplay->DestroyEventRPhi();
262 gDelphesDisplay->ImportEventRPhi(top);
263 gDelphesDisplay->DestroyEventRhoZ();
264 gDelphesDisplay->ImportEventRhoZ(top);
265 //update_html_summary();
266 gEve->Redraw3D(kFALSE, kTRUE);
267}
268
269void DelphesEventDisplay::delphes_read()
270{
271
272 // safety
273 if(event_id >= gTreeReader->GetEntries() || event_id<0 ) return;
274
275 // Load selected branches with data from specified event
276 gTreeReader->ReadEntry(event_id);
277
278 // loop over selected branches, and apply the proper recipe to fill the collections.
279 // this is basically to loop on gArrays to fill gElements.
280
281//TODO: one option would be to have templated methods in the element classes. We could simply call "element.fill()"
282 std::vector<TClonesArray*>::iterator data = gArrays.begin();
283 std::vector<DelphesBranchBase*>::iterator element = gElements.begin();
284 for(; data<gArrays.end() && element<gElements.end(); ++data, ++element) {
285 TString type = (*element)->GetType();
286 // branch on the element type
287 if(type=="tower") delphes_read_towers(*data,*element);
288 else if(type=="track" || type=="photon" || type=="electron" || type=="muon" || type=="genparticle") delphes_read_tracks(*data,*element);
289 else if(type=="jet") delphes_read_jets(*data,*element);
290 else if(type=="vector") delphes_read_vectors(*data,*element);
291 }
292}
293
294void DelphesEventDisplay::delphes_read_towers(TClonesArray* data, DelphesBranchBase* element) {
295 DelphesCaloData* container = dynamic_cast<DelphesBranchElement<DelphesCaloData>*>(element)->GetContainer();
296 assert(container);
297 // Loop over all towers
298 TIter itTower(data);
299 Tower *tower;
300 while((tower = (Tower *) itTower.Next()))
301 {
302 container->AddTower(tower->Edges[0], tower->Edges[1], tower->Edges[2], tower->Edges[3]);
303 container->FillSlice(0, tower->Eem);
304 container->FillSlice(1, tower->Ehad);
305 }
306 container->DataChanged();
307}
308
309void DelphesEventDisplay::delphes_read_tracks(TClonesArray* data, DelphesBranchBase* element) {
310 TEveTrackList* container = dynamic_cast<DelphesBranchElement<TEveTrackList>*>(element)->GetContainer();
311 assert(container);
312 TString type = element->GetType();
313 TIter itTrack(data);
314 Int_t counter = 0;
315 TEveTrack *eveTrack;
316 TEveTrackPropagator *trkProp = container->GetPropagator();
317 if(type=="track") {
318 // Loop over all tracks
319 Track *track;
320 while((track = (Track *) itTrack.Next())) {
321 TParticle pb(track->PID, 1, 0, 0, 0, 0,
322 track->P4().Px(), track->P4().Py(),
323 track->P4().Pz(), track->P4().E(),
324 track->X, track->Y, track->Z, 0.0);
325 eveTrack = new TEveTrack(&pb, counter, trkProp);
326 eveTrack->SetName(Form("%s [%d]", pb.GetName(), counter++));
327 eveTrack->SetStdTitle();
328 eveTrack->SetAttLineAttMarker(container);
329 container->AddElement(eveTrack);
330 eveTrack->SetLineColor(element->GetColor());
331 eveTrack->MakeTrack();
332 }
333 } else if(type=="electron") {
334 // Loop over all electrons
335 Electron *electron;
336 while((electron = (Electron *) itTrack.Next())) {
337 TParticle pb(electron->Charge<0?11:-11, 1, 0, 0, 0, 0,
338 electron->P4().Px(), electron->P4().Py(),
339 electron->P4().Pz(), electron->P4().E(),
340 0., 0., 0., 0.);
341
342 eveTrack = new TEveTrack(&pb, counter, trkProp);
343 eveTrack->SetName(Form("%s [%d]", pb.GetName(), counter++));
344 eveTrack->SetStdTitle();
345 eveTrack->SetAttLineAttMarker(container);
346 container->AddElement(eveTrack);
347 eveTrack->SetLineColor(element->GetColor());
348 eveTrack->MakeTrack();
349 }
350 } else if(type=="muon") {
351 // Loop over all muons
352 Muon *muon;
353 while((muon = (Muon *) itTrack.Next())) {
354 TParticle pb(muon->Charge<0?13:-13, 1, 0, 0, 0, 0,
355 muon->P4().Px(), muon->P4().Py(),
356 muon->P4().Pz(), muon->P4().E(),
357 0., 0., 0., 0.);
358
359 eveTrack = new TEveTrack(&pb, counter, trkProp);
360 eveTrack->SetName(Form("%s [%d]", pb.GetName(), counter++));
361 eveTrack->SetStdTitle();
362 eveTrack->SetAttLineAttMarker(container);
363 container->AddElement(eveTrack);
364 eveTrack->SetLineColor(element->GetColor());
365 eveTrack->MakeTrack();
366 }
367 } else if(type=="photon") {
368 // Loop over all photons
369 Photon *photon;
370 while((photon = (Photon *) itTrack.Next())) {
371 TParticle pb(22, 1, 0, 0, 0, 0,
372 photon->P4().Px(), photon->P4().Py(),
373 photon->P4().Pz(), photon->P4().E(),
374 0., 0., 0., 0.);
375 eveTrack = new TEveTrack(&pb, counter, trkProp);
376 eveTrack->SetName(Form("%s [%d]", pb.GetName(), counter++));
377 eveTrack->SetStdTitle();
378 eveTrack->SetAttLineAttMarker(container);
379 eveTrack->SetLineStyle(7);
380 container->AddElement(eveTrack);
381 eveTrack->SetLineColor(element->GetColor());
382 eveTrack->MakeTrack();
383 }
384 } else if(type=="genparticle") {
385 // Loop over all particles
386 GenParticle *particle;
387 while((particle = (GenParticle *) itTrack.Next())) {
388 TParticle pb(particle->PID, particle->Status, particle->M1, particle->M2, particle->D1, particle->D2,
389 particle->P4().Px(), particle->P4().Py(),
390 particle->P4().Pz(), particle->P4().E(),
391 particle->X, particle->Y, particle->Z, particle->T);
392 eveTrack = new TEveTrack(&pb, counter, trkProp);
393 eveTrack->SetName(Form("%s [%d]", pb.GetName(), counter++));
394 eveTrack->SetStdTitle();
395 eveTrack->SetAttLineAttMarker(container);
396 container->AddElement(eveTrack);
397 eveTrack->SetLineColor(element->GetColor());
398 if(particle->Charge==0) eveTrack->SetLineStyle(7);
399 eveTrack->MakeTrack();
400 }
401 }
402}
403
404void DelphesEventDisplay::delphes_read_jets(TClonesArray* data, DelphesBranchBase* element) {
405 TEveElementList* container = dynamic_cast<DelphesBranchElement<TEveElementList>*>(element)->GetContainer();
406 assert(container);
407 TIter itJet(data);
408 Jet *jet;
409 TEveJetCone *eveJetCone;
410 // Loop over all jets
411 Int_t counter = 0;
412 while((jet = (Jet *) itJet.Next()))
413 {
414 eveJetCone = new TEveJetCone();
415 eveJetCone->SetTitle(Form("jet [%d]: Pt=%f, Eta=%f, \nPhi=%f, M=%f",counter,jet->PT, jet->Eta, jet->Phi, jet->Mass));
416 eveJetCone->SetName(Form("jet [%d]", counter++));
417 eveJetCone->SetMainTransparency(60);
418 eveJetCone->SetLineColor(element->GetColor());
419 eveJetCone->SetFillColor(element->GetColor());
420 eveJetCone->SetCylinder(tkRadius_ - 10, tkHalfLength_ - 10);
421 eveJetCone->SetPickable(kTRUE);
422 eveJetCone->AddEllipticCone(jet->Eta, jet->Phi, jet->DeltaEta, jet->DeltaPhi);
423 container->AddElement(eveJetCone);
424 }
425}
426
427void DelphesEventDisplay::delphes_read_vectors(TClonesArray* data, DelphesBranchBase* element) {
428 TEveElementList* container = dynamic_cast<DelphesBranchElement<TEveElementList>*>(element)->GetContainer();
429 assert(container);
430 TIter itMet(data);
431 MissingET *MET;
432 TEveArrow *eveMet;
433 // Missing Et
434 Double_t maxPt = 50.;
435 // TODO to be changed as we don't have access to maxPt anymore. MET scale could be a general parameter set in GUI
436 while((MET = (MissingET*) itMet.Next())) {
437 eveMet = new TEveArrow((tkRadius_ * MET->MET/maxPt)*cos(MET->Phi), (tkRadius_ * MET->MET/maxPt)*sin(MET->Phi), 0., 0., 0., 0.);
438 eveMet->SetMainColor(element->GetColor());
439 eveMet->SetTubeR(0.04);
440 eveMet->SetConeR(0.08);
441 eveMet->SetConeL(0.10);
442 eveMet->SetPickable(kTRUE);
443 eveMet->SetName("Missing Et");
444 eveMet->SetTitle(Form("Missing Et (%.1f GeV)",MET->MET));
445 container->AddElement(eveMet);
446 }
447}
448
449/******************************************************************************/
450// GUI
451/******************************************************************************/
452
453void DelphesEventDisplay::make_gui()
454{
455 // Create minimal GUI for event navigation.
456 // TODO: better GUI could be made based on the ch15 of the manual (Writing a GUI)
457
458 // add a tab on the left
459 TEveBrowser* browser = gEve->GetBrowser();
460 browser->StartEmbedding(TRootBrowser::kLeft);
461
462 // set the main title
463 TGMainFrame* frmMain = new TGMainFrame(gClient->GetRoot(), 1000, 600);
464 frmMain->SetWindowName("Delphes Event Display");
465 frmMain->SetCleanup(kDeepCleanup);
466
467 // build the navigation menu
468 TGHorizontalFrame* hf = new TGHorizontalFrame(frmMain);
469 {
470 TString icondir;
471 if(gSystem->Getenv("ROOTSYS"))
472 icondir = Form("%s/icons/", gSystem->Getenv("ROOTSYS"));
473 if(!gSystem->OpenDirectory(icondir))
474 icondir = Form("%s/icons/", (const char*)gSystem->GetFromPipe("root-config --etcdir") );
475 TGPictureButton* b = 0;
476 EvNavHandler *fh = new EvNavHandler(this);
477
478 b = new TGPictureButton(hf, gClient->GetPicture(icondir+"GoBack.gif"));
479 hf->AddFrame(b);
480 b->Connect("Clicked()", "EvNavHandler", fh, "Bck()");
481
482 b = new TGPictureButton(hf, gClient->GetPicture(icondir+"GoForward.gif"));
483 hf->AddFrame(b);
484 b->Connect("Clicked()", "EvNavHandler", fh, "Fwd()");
485 }
486 frmMain->AddFrame(hf);
487 frmMain->MapSubwindows();
488 frmMain->Resize();
489 frmMain->MapWindow();
490 browser->StopEmbedding();
491 browser->SetTabTitle("Event Control", 0);
492}
493
Note: See TracBrowser for help on using the repository browser.