Fork me on GitHub

source: git/display/DelphesEventDisplay.cc@ 2ca23b5

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

Started to work on the summary plots

Investigated the way to add a tab with summary plots. A new class is in
place (dummy for now) to manage these plots.
Other small changes: use the status bar, clean the maxPt for MET (easier
to change later on), filter GenParticles on status==1.

  • Property mode set to 100644
File size: 15.5 KB
Line 
1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
4 *
5 * This program is free software: you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation, either version 3 of the License, or
8 * (at your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program. If not, see <http://www.gnu.org/licenses/>.
17 */
18
19#include <cassert>
20#include <iostream>
21#include <utility>
22#include <algorithm>
23#include "TGeoManager.h"
24#include "TGeoVolume.h"
25#include "external/ExRootAnalysis/ExRootConfReader.h"
26#include "external/ExRootAnalysis/ExRootTreeReader.h"
27#include "display/DelphesCaloData.h"
28#include "display/DelphesBranchElement.h"
29#include "display/Delphes3DGeometry.h"
30#include "display/DelphesEventDisplay.h"
31#include "classes/DelphesClasses.h"
32#include "TEveElement.h"
33#include "TEveJetCone.h"
34#include "TEveTrack.h"
35#include "TEveTrackPropagator.h"
36#include "TEveCalo.h"
37#include "TEveManager.h"
38#include "TEveGeoNode.h"
39#include "TEveTrans.h"
40#include "TEveViewer.h"
41#include "TEveBrowser.h"
42#include "TEveArrow.h"
43#include "TMath.h"
44#include "TSystem.h"
45#include "TRootBrowser.h"
46#include "TGButton.h"
47#include "TRootEmbeddedCanvas.h"
48#include "TClonesArray.h"
49#include "TEveEventManager.h"
50#include "TCanvas.h"
51#include "TH1F.h"
52
53DelphesEventDisplay::DelphesEventDisplay()
54{
55 event_id_ = 0;
56 tkRadius_ = 1.29;
57 totRadius_ = 2.0;
58 tkHalfLength_ = 3.0;
59 muHalfLength_ = 6.0;
60 bz_ = 3.8;
61 chain_ = new TChain("Delphes");
62 treeReader_ = 0;
63 delphesDisplay_ = 0;
64 etaAxis_ = 0;
65 phiAxis_ = 0;
66}
67
68DelphesEventDisplay::~DelphesEventDisplay()
69{
70 delete chain_;
71}
72
73DelphesEventDisplay::DelphesEventDisplay(const char *configFile, const char *inputFile, Delphes3DGeometry& det3D)
74{
75 event_id_ = 0;
76 tkRadius_ = 1.29;
77 totRadius_ = 2.0;
78 tkHalfLength_ = 3.0;
79 bz_ = 3.8;
80 chain_ = new TChain("Delphes");
81 treeReader_ = 0;
82 delphesDisplay_ = 0;
83
84 // initialize the application
85 TEveManager::Create(kTRUE, "IV");
86 fStatusBar_ = gEve->GetBrowser()->GetStatusBar();
87 TGeoManager* geom = gGeoManager;
88
89 // build the detector
90 tkRadius_ = det3D.getTrackerRadius();
91 totRadius_ = det3D.getDetectorRadius();
92 tkHalfLength_ = det3D.getTrackerHalfLength();
93 muHalfLength_ = det3D.getDetectorHalfLength();
94 bz_ = det3D.getBField();
95 etaAxis_ = det3D.getCaloAxes().first;
96 phiAxis_ = det3D.getCaloAxes().second;
97 TGeoVolume* top = det3D.getDetector(false);
98 geom->SetTopVolume(top);
99 TEveElementList *geometry = new TEveElementList("Geometry");
100 TObjArray* nodes = top->GetNodes();
101 TIter itNodes(nodes);
102 TGeoNode* nodeobj;
103 TEveGeoTopNode* node;
104 while((nodeobj = (TGeoNode*)itNodes.Next())) {
105 node = new TEveGeoTopNode(gGeoManager,nodeobj);
106 node->UseNodeTrans();
107 geometry->AddElement(node);
108 }
109
110 // Create chain of root trees
111 chain_->Add(inputFile);
112
113 // Create object of class ExRootTreeReader
114 fStatusBar_->SetText("Opening Delphes data file", 1);
115 treeReader_ = new ExRootTreeReader(chain_);
116
117 // prepare data collections
118 readConfig(configFile, elements_);
119 for(std::vector<DelphesBranchBase*>::iterator element = elements_.begin(); element<elements_.end(); ++element) {
120 DelphesBranchElement<TEveTrackList>* item_v1 = dynamic_cast<DelphesBranchElement<TEveTrackList>*>(*element);
121 DelphesBranchElement<TEveElementList>* item_v2 = dynamic_cast<DelphesBranchElement<TEveElementList>*>(*element);
122 if(item_v1) gEve->AddElement(item_v1->GetContainer());
123 if(item_v2) gEve->AddElement(item_v2->GetContainer());
124 }
125
126 // viewers and scenes
127 delphesDisplay_ = new DelphesDisplay;
128 gEve->AddGlobalElement(geometry);
129 delphesDisplay_->ImportGeomRPhi(geometry);
130 delphesDisplay_->ImportGeomRhoZ(geometry);
131 // find the first calo data and use that to initialize the calo display
132 for(std::vector<DelphesBranchBase*>::iterator data=elements_.begin();data<elements_.end();++data) {
133 if(TString((*data)->GetType())=="Tower") { // we could also use GetClassName()=="DelphesCaloData"
134 DelphesCaloData* container = dynamic_cast<DelphesBranchElement<DelphesCaloData>*>((*data))->GetContainer();
135 assert(container);
136 TEveCalo3D *calo3d = new TEveCalo3D(container);
137 calo3d->SetBarrelRadius(tkRadius_);
138 calo3d->SetEndCapPos(tkHalfLength_);
139 gEve->AddGlobalElement(calo3d);
140 delphesDisplay_->ImportCaloRPhi(calo3d);
141 delphesDisplay_->ImportCaloRhoZ(calo3d);
142 TEveCaloLego *lego = new TEveCaloLego(container);
143 lego->InitMainTrans();
144 lego->RefMainTrans().SetScale(TMath::TwoPi(), TMath::TwoPi(), TMath::Pi());
145 lego->SetAutoRebin(kFALSE);
146 lego->Set2DMode(TEveCaloLego::kValSizeOutline);
147 delphesDisplay_->ImportCaloLego(lego);
148 break;
149 }
150 }
151
152 // the GUI: control panel, summary tab
153 make_gui();
154
155 //ready...
156 fStatusBar_->SetText("Ready.", 1);
157 load_event();
158 gEve->Redraw3D(kTRUE);
159
160}
161
162// function that parses the config to extract the branches of interest and prepare containers
163void DelphesEventDisplay::readConfig(const char *configFile, std::vector<DelphesBranchBase*>& elements) {
164 ExRootConfReader *confReader = new ExRootConfReader;
165 confReader->ReadFile(configFile);
166 ExRootConfParam branches = confReader->GetParam("TreeWriter::Branch");
167 Int_t nBranches = branches.GetSize()/3;
168 DelphesBranchElement<TEveTrackList>* tlist;
169 DelphesBranchElement<DelphesCaloData>* clist;
170 DelphesBranchElement<TEveElementList>* elist;
171 // first loop with all but tracks
172 for(Int_t b = 0; b<nBranches; ++b) {
173 TString input = branches[b*3].GetString();
174 TString name = branches[b*3+1].GetString();
175 TString className = branches[b*3+2].GetString();
176 if(className=="Tower") {
177 if(input.Contains("eflow",TString::kIgnoreCase) || name.Contains("eflow",TString::kIgnoreCase)) continue; //no eflow
178 clist = new DelphesBranchElement<DelphesCaloData>(name,treeReader_->UseBranch(name),kBlack);
179 clist->GetContainer()->SetEtaBins(etaAxis_);
180 clist->GetContainer()->SetPhiBins(phiAxis_);
181 elements.push_back(clist);
182 } else if(className=="Jet") {
183 if(input.Contains("GenJetFinder")) {
184 elist = new DelphesBranchElement<TEveElementList>(name,treeReader_->UseBranch(name),kCyan);
185 elist->GetContainer()->SetRnrSelf(false);
186 elist->GetContainer()->SetRnrChildren(false);
187 elist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
188 elements.push_back(elist);
189 } else {
190 elist = new DelphesBranchElement<TEveElementList>(name,treeReader_->UseBranch(name),kYellow);
191 elist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
192 elements.push_back(elist);
193 }
194 } else if(className=="Electron") {
195 tlist = new DelphesBranchElement<TEveTrackList>(name,treeReader_->UseBranch(name),kRed);
196 tlist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
197 elements.push_back(tlist);
198 } else if(className=="Photon") {
199 tlist = new DelphesBranchElement<TEveTrackList>(name,treeReader_->UseBranch(name),kYellow);
200 tlist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
201 elements.push_back(tlist);
202 } else if(className=="Muon") {
203 tlist = new DelphesBranchElement<TEveTrackList>(name,treeReader_->UseBranch(name),kGreen);
204 tlist->SetTrackingVolume(totRadius_, muHalfLength_, bz_);
205 elements.push_back(tlist);
206 } else if(className=="MissingET") {
207 elist = new DelphesBranchElement<TEveElementList>(name,treeReader_->UseBranch(name),kViolet);
208 elist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
209 elements.push_back(elist);
210 } else if(className=="GenParticle") {
211 tlist = new DelphesBranchElement<TEveTrackList>(name,treeReader_->UseBranch(name),kCyan);
212 tlist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
213 tlist->GetContainer()->SetRnrSelf(false);
214 tlist->GetContainer()->SetRnrChildren(false);
215 elements.push_back(tlist);
216 } else {
217 continue;
218 }
219 }
220 // second loop for tracks
221 for(Int_t b = 0; b<nBranches; ++b) {
222 TString input = branches[b*3].GetString();
223 TString name = branches[b*3+1].GetString();
224 TString className = branches[b*3+2].GetString();
225 if(className=="Track") {
226 if(input.Contains("eflow",TString::kIgnoreCase) || name.Contains("eflow",TString::kIgnoreCase)) continue; //no eflow
227 tlist = new DelphesBranchElement<TEveTrackList>(name,treeReader_->UseBranch(name),kBlue);
228 tlist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
229 elements.push_back(tlist);
230 }
231 }
232}
233
234void DelphesEventDisplay::load_event()
235{
236 // Load event specified in global event_id_.
237 // The contents of previous event are removed.
238
239 // safety
240 if(event_id_ >= treeReader_->GetEntries() || event_id_<0 ) return;
241
242 // message
243 fStatusBar_->SetText(Form("Loading event %d.", event_id_), 1);
244
245 // clear the previous event
246 gEve->GetViewers()->DeleteAnnotations();
247 for(std::vector<DelphesBranchBase*>::iterator data=elements_.begin();data<elements_.end();++data) {
248 (*data)->Reset();
249 }
250
251 // Load selected branches with data from specified event
252 treeReader_->ReadEntry(event_id_);
253 for(std::vector<DelphesBranchBase*>::iterator data=elements_.begin();data<elements_.end();++data) {
254 (*data)->ReadBranch();
255 }
256
257 // update display
258 TEveElement* top = (TEveElement*)gEve->GetCurrentEvent();
259 delphesDisplay_->DestroyEventRPhi();
260 delphesDisplay_->ImportEventRPhi(top);
261 delphesDisplay_->DestroyEventRhoZ();
262 delphesDisplay_->ImportEventRhoZ(top);
263 update_html_summary();
264 //TODO: update plot tab (show current event on top)
265
266 gEve->Redraw3D(kFALSE, kTRUE);
267 fStatusBar_->SetText(Form("Loaded event %d.", event_id_), 1);
268}
269
270void DelphesEventDisplay::update_html_summary()
271{
272 // Update summary of current event.
273
274 TEveElement::List_i i;
275 TEveElement::List_i j;
276 Int_t k;
277 TEveElement *el;
278 DelphesHtmlObjTable *table;
279 TEveEventManager *mgr = gEve ? gEve->GetCurrentEvent() : 0;
280 if (mgr) {
281 htmlSummary_->Clear("D");
282 for (i=mgr->BeginChildren(); i!=mgr->EndChildren(); ++i) {
283 el = ((TEveElement*)(*i));
284 if (el->IsA() == TEvePointSet::Class()) {
285 TEvePointSet *ps = (TEvePointSet *)el;
286 TString ename = ps->GetElementName();
287 TString etitle = ps->GetElementTitle();
288 if (ename.First('\'') != kNPOS)
289 ename.Remove(ename.First('\''));
290 etitle.Remove(0, 2);
291 Int_t nel = atoi(etitle.Data());
292 table = htmlSummary_->AddTable(ename, 0, nel);
293 }
294 else if (el->IsA() == TEveTrackList::Class()) {
295 TEveTrackList *tracks = (TEveTrackList *)el;
296 TString ename = tracks->GetElementName();
297 if (ename.First('\'') != kNPOS)
298 ename.Remove(ename.First('\''));
299 table = htmlSummary_->AddTable(ename.Data(), 5,
300 tracks->NumChildren(), kTRUE, "first");
301 table->SetLabel(0, "Momentum");
302 table->SetLabel(1, "P_t");
303 table->SetLabel(2, "Phi");
304 table->SetLabel(3, "Theta");
305 table->SetLabel(4, "Eta");
306 k=0;
307 for (j=tracks->BeginChildren(); j!=tracks->EndChildren(); ++j) {
308 Float_t p = ((TEveTrack*)(*j))->GetMomentum().Mag();
309 table->SetValue(0, k, p);
310 Float_t pt = ((TEveTrack*)(*j))->GetMomentum().Perp();
311 table->SetValue(1, k, pt);
312 Float_t phi = ((TEveTrack*)(*j))->GetMomentum().Phi();
313 table->SetValue(2, k, phi);
314 Float_t theta = ((TEveTrack*)(*j))->GetMomentum().Theta();
315 table->SetValue(3, k, theta);
316 Float_t eta = theta>0.0005 && theta<3.1413 ? ((TEveTrack*)(*j))->GetMomentum().Eta() : 1e10;
317 table->SetValue(4, k, eta);
318 ++k;
319 }
320 }
321 }
322 htmlSummary_->Build();
323 gHtml_->Clear();
324 gHtml_->ParseText((char*)htmlSummary_->Html().Data());
325 gHtml_->Layout();
326 }
327
328}
329
330/******************************************************************************/
331// GUI
332/******************************************************************************/
333
334void DelphesEventDisplay::make_gui()
335{
336 // Create minimal GUI for event navigation.
337 // TODO: better GUI could be made based on the ch15 of the manual (Writing a GUI)
338
339 // add a tab on the left
340 TEveBrowser* browser = gEve->GetBrowser();
341 browser->StartEmbedding(TRootBrowser::kLeft);
342
343 // set the main title
344 TGMainFrame* frmMain = new TGMainFrame(gClient->GetRoot(), 1000, 600);
345 frmMain->SetWindowName("Delphes Event Display");
346 frmMain->SetCleanup(kDeepCleanup);
347
348 // build the navigation menu
349 TGHorizontalFrame* hf = new TGHorizontalFrame(frmMain);
350 {
351 TString icondir;
352 if(gSystem->Getenv("ROOTSYS"))
353 icondir = Form("%s/icons/", gSystem->Getenv("ROOTSYS"));
354 if(!gSystem->OpenDirectory(icondir))
355 icondir = Form("%s/icons/", (const char*)gSystem->GetFromPipe("root-config --etcdir") );
356 TGPictureButton* b = 0;
357
358 b = new TGPictureButton(hf, gClient->GetPicture(icondir+"GoBack.gif"));
359 hf->AddFrame(b);
360 b->Connect("Clicked()", "DelphesEventDisplay", this, "Bck()");
361
362 b = new TGPictureButton(hf, gClient->GetPicture(icondir+"GoForward.gif"));
363 hf->AddFrame(b);
364 b->Connect("Clicked()", "DelphesEventDisplay", this, "Fwd()");
365 }
366 frmMain->AddFrame(hf);
367 frmMain->MapSubwindows();
368 frmMain->Resize();
369 frmMain->MapWindow();
370 browser->StopEmbedding();
371 browser->SetTabTitle("Event Control", 0);
372
373 // the summary tab
374 htmlSummary_ = new DelphesHtmlSummary("Delphes Event Display Summary Table");
375 TEveWindowSlot* slot = TEveWindow::CreateWindowInTab(gEve->GetBrowser()->GetTabRight());
376 gHtml_ = new TGHtml(0, 100, 100);
377 TEveWindowFrame *wf = slot->MakeFrame(gHtml_);
378 gHtml_->MapSubwindows();
379 wf->SetElementName("Summary tables");
380
381 // plot tab
382 slot = TEveWindow::CreateWindowInTab(gEve->GetBrowser()->GetTabRight());
383 TEveWindowTab* tab = slot->MakeTab();
384 tab->SetElementName("Summary plots");
385 tab->SetShowTitleBar(kFALSE);
386 plotSummary_ = new DelphesPlotSummary(tab);
387 plotSummary_->Init(elements_);
388 plotSummary_->FillSample(treeReader_); //TODO later, control it via a GUI button.
389 plotSummary_->FillEvent(); //TODO later move to event loop
390
391 //for test
392 TH1F* h;
393 TRootEmbeddedCanvas* trec;
394 TCanvas* gCanvas_;
395
396 slot = tab->NewSlot();
397 trec = new TRootEmbeddedCanvas();
398 gCanvas_ = trec->GetCanvas();
399 wf = slot->MakeFrame(trec);
400 wf->SetElementName("Tracks");
401 h = new TH1F("tracks","tracks",100,0,100);
402 gCanvas_->cd();
403 h->Draw();
404
405 slot = tab->NewSlot();
406 trec = new TRootEmbeddedCanvas();
407 gCanvas_ = trec->GetCanvas();
408 wf = slot->MakeFrame(trec);
409 wf->SetElementName("Electrons");
410 h = new TH1F("electrons","electrons",100,0,100);
411 gCanvas_->cd();
412 h->Draw();
413
414 // TODO: here we have, for each collection, Pt,Eta,Phi for all, leading, subleading
415 // for each event, we will then add a marker with the current value and/or a histo for current event.
416 // this means to create one tab with subtabs (one per collection).
417
418}
419
Note: See TracBrowser for help on using the repository browser.