Fork me on GitHub

source: git/display/DelphesEventDisplay.cc@ fc6300d

ImprovedOutputFile Timing dual_readout llp
Last change on this file since fc6300d was 1fa50c2, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

fix GPLv3 header

  • Property mode set to 100644
File size: 16.8 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 "TGTextEntry.h"
48#include "TGProgressBar.h"
49#include "TGNumberEntry.h"
50#include "TRootEmbeddedCanvas.h"
51#include "TClonesArray.h"
52#include "TEveEventManager.h"
53#include "TCanvas.h"
54#include "TH1F.h"
55
56DelphesEventDisplay::DelphesEventDisplay()
57{
58 event_id_ = 0;
59 tkRadius_ = 1.29;
60 totRadius_ = 2.0;
61 tkHalfLength_ = 3.0;
62 muHalfLength_ = 6.0;
63 bz_ = 3.8;
64 chain_ = new TChain("Delphes");
65 treeReader_ = 0;
66 delphesDisplay_ = 0;
67 etaAxis_ = 0;
68 phiAxis_ = 0;
69}
70
71DelphesEventDisplay::~DelphesEventDisplay()
72{
73 delete chain_;
74}
75
76void DelphesEventDisplay::EventChanged(Int_t e)
77{
78 if( e!= event_id_ ) {
79 event_id_ = e;
80 Emit("EventChanged(Int_t)",e);
81 load_event();
82 }
83}
84
85
86DelphesEventDisplay::DelphesEventDisplay(const char *configFile, const char *inputFile, Delphes3DGeometry& det3D)
87{
88 event_id_ = 0;
89 tkRadius_ = 1.29;
90 totRadius_ = 2.0;
91 tkHalfLength_ = 3.0;
92 bz_ = 3.8;
93 chain_ = new TChain("Delphes");
94 treeReader_ = 0;
95 delphesDisplay_ = 0;
96
97 // initialize the application
98 TEveManager::Create(kTRUE, "IV");
99 fStatusBar_ = gEve->GetBrowser()->GetStatusBar();
100 TGeoManager* geom = gGeoManager;
101
102 // build the detector
103 tkRadius_ = det3D.getTrackerRadius();
104 totRadius_ = det3D.getDetectorRadius();
105 tkHalfLength_ = det3D.getTrackerHalfLength();
106 muHalfLength_ = det3D.getDetectorHalfLength();
107 bz_ = det3D.getBField();
108 etaAxis_ = det3D.getCaloAxes().first;
109 phiAxis_ = det3D.getCaloAxes().second;
110 TGeoVolume* top = det3D.getDetector(false);
111 geom->SetTopVolume(top);
112 TEveElementList *geometry = new TEveElementList("Geometry");
113 TObjArray* nodes = top->GetNodes();
114 TIter itNodes(nodes);
115 TGeoNode* nodeobj;
116 TEveGeoTopNode* node;
117 while((nodeobj = (TGeoNode*)itNodes.Next())) {
118 node = new TEveGeoTopNode(gGeoManager,nodeobj);
119 node->UseNodeTrans();
120 geometry->AddElement(node);
121 }
122
123 // Create chain of root trees
124 chain_->Add(inputFile);
125
126 // Create object of class ExRootTreeReader
127 fStatusBar_->SetText("Opening Delphes data file", 1);
128 gSystem->ProcessEvents();
129 treeReader_ = new ExRootTreeReader(chain_);
130
131 // prepare data collections
132 readConfig(configFile, elements_);
133 for(std::vector<DelphesBranchBase*>::iterator element = elements_.begin(); element<elements_.end(); ++element) {
134 DelphesBranchElement<TEveTrackList>* item_v1 = dynamic_cast<DelphesBranchElement<TEveTrackList>*>(*element);
135 DelphesBranchElement<TEveElementList>* item_v2 = dynamic_cast<DelphesBranchElement<TEveElementList>*>(*element);
136 if(item_v1) gEve->AddElement(item_v1->GetContainer());
137 if(item_v2) gEve->AddElement(item_v2->GetContainer());
138 }
139
140 // viewers and scenes
141 delphesDisplay_ = new DelphesDisplay;
142 gEve->AddGlobalElement(geometry);
143 delphesDisplay_->ImportGeomRPhi(geometry);
144 delphesDisplay_->ImportGeomRhoZ(geometry);
145 // find the first calo data and use that to initialize the calo display
146 for(std::vector<DelphesBranchBase*>::iterator data=elements_.begin();data<elements_.end();++data) {
147 if(TString((*data)->GetType())=="Tower") { // we could also use GetClassName()=="DelphesCaloData"
148 DelphesCaloData* container = dynamic_cast<DelphesBranchElement<DelphesCaloData>*>((*data))->GetContainer();
149 assert(container);
150 TEveCalo3D *calo3d = new TEveCalo3D(container);
151 calo3d->SetBarrelRadius(tkRadius_);
152 calo3d->SetEndCapPos(tkHalfLength_);
153 gEve->AddGlobalElement(calo3d);
154 delphesDisplay_->ImportCaloRPhi(calo3d);
155 delphesDisplay_->ImportCaloRhoZ(calo3d);
156 TEveCaloLego *lego = new TEveCaloLego(container);
157 lego->InitMainTrans();
158 // lego->RefMainTrans().SetScale(TMath::TwoPi(), TMath::TwoPi(), TMath::Pi());
159 lego->RefMainTrans().SetScale(100, 100, TMath::Pi());
160 lego->SetAutoRebin(kFALSE);
161 lego->Set2DMode(TEveCaloLego::kValSizeOutline);
162 delphesDisplay_->ImportCaloLego(lego);
163 break;
164 }
165 }
166
167 // the GUI: control panel, summary tab
168 make_gui();
169
170 //ready...
171 fStatusBar_->SetText("Ready.", 1);
172 gSystem->ProcessEvents();
173 load_event();
174 gEve->Redraw3D(kTRUE);
175
176}
177
178// function that parses the config to extract the branches of interest and prepare containers
179void DelphesEventDisplay::readConfig(const char *configFile, std::vector<DelphesBranchBase*>& elements) {
180 ExRootConfReader *confReader = new ExRootConfReader;
181 confReader->ReadFile(configFile);
182 ExRootConfParam branches = confReader->GetParam("TreeWriter::Branch");
183 Int_t nBranches = branches.GetSize()/3;
184 DelphesBranchElement<TEveTrackList>* tlist;
185 DelphesBranchElement<DelphesCaloData>* clist;
186 DelphesBranchElement<TEveElementList>* elist;
187 // first loop with all but tracks
188 for(Int_t b = 0; b<nBranches; ++b) {
189 TString input = branches[b*3].GetString();
190 TString name = branches[b*3+1].GetString();
191 TString className = branches[b*3+2].GetString();
192 if(className=="Tower") {
193 if(input.Contains("eflow",TString::kIgnoreCase) || name.Contains("eflow",TString::kIgnoreCase)) continue; //no eflow
194 clist = new DelphesBranchElement<DelphesCaloData>(name,treeReader_->UseBranch(name),kBlack);
195 clist->GetContainer()->SetEtaBins(etaAxis_);
196 clist->GetContainer()->SetPhiBins(phiAxis_);
197 elements.push_back(clist);
198 } else if(className=="Jet") {
199 if(input.Contains("GenJetFinder")) {
200 elist = new DelphesBranchElement<TEveElementList>(name,treeReader_->UseBranch(name),kCyan);
201 elist->GetContainer()->SetRnrSelf(false);
202 elist->GetContainer()->SetRnrChildren(false);
203 elist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
204 elements.push_back(elist);
205 } else {
206 elist = new DelphesBranchElement<TEveElementList>(name,treeReader_->UseBranch(name),kYellow);
207 elist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
208 elements.push_back(elist);
209 }
210 } else if(className=="Electron") {
211 tlist = new DelphesBranchElement<TEveTrackList>(name,treeReader_->UseBranch(name),kRed);
212 tlist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
213 elements.push_back(tlist);
214 } else if(className=="Photon") {
215 tlist = new DelphesBranchElement<TEveTrackList>(name,treeReader_->UseBranch(name),kYellow);
216 tlist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
217 elements.push_back(tlist);
218 } else if(className=="Muon") {
219 tlist = new DelphesBranchElement<TEveTrackList>(name,treeReader_->UseBranch(name),kGreen);
220 tlist->SetTrackingVolume(totRadius_, muHalfLength_, bz_);
221 elements.push_back(tlist);
222 } else if(className=="MissingET") {
223 elist = new DelphesBranchElement<TEveElementList>(name,treeReader_->UseBranch(name),kViolet);
224 elist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
225 elements.push_back(elist);
226 } else if(className=="GenParticle") {
227 tlist = new DelphesBranchElement<TEveTrackList>(name,treeReader_->UseBranch(name),kCyan);
228 tlist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
229 tlist->GetContainer()->SetRnrSelf(false);
230 tlist->GetContainer()->SetRnrChildren(false);
231 elements.push_back(tlist);
232 } else {
233 continue;
234 }
235 }
236 // second loop for tracks
237 for(Int_t b = 0; b<nBranches; ++b) {
238 TString input = branches[b*3].GetString();
239 TString name = branches[b*3+1].GetString();
240 TString className = branches[b*3+2].GetString();
241 if(className=="Track") {
242 if(input.Contains("eflow",TString::kIgnoreCase) || name.Contains("eflow",TString::kIgnoreCase)) continue; //no eflow
243 tlist = new DelphesBranchElement<TEveTrackList>(name,treeReader_->UseBranch(name),kBlue);
244 tlist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
245 elements.push_back(tlist);
246 }
247 }
248}
249
250void DelphesEventDisplay::load_event()
251{
252 // Load event specified in global event_id_.
253 // The contents of previous event are removed.
254
255 // safety
256 if(event_id_ >= treeReader_->GetEntries() || event_id_<0 ) return;
257
258 // message
259 fStatusBar_->SetText(Form("Loading event %d.", event_id_), 1);
260 gSystem->ProcessEvents();
261
262 // clear the previous event
263 gEve->GetViewers()->DeleteAnnotations();
264 for(std::vector<DelphesBranchBase*>::iterator data=elements_.begin();data<elements_.end();++data) {
265 (*data)->Reset();
266 }
267
268 // Load selected branches with data from specified event
269 treeReader_->ReadEntry(event_id_);
270 for(std::vector<DelphesBranchBase*>::iterator data=elements_.begin();data<elements_.end();++data) {
271 (*data)->ReadBranch();
272 }
273
274 // update display
275 TEveElement* top = (TEveElement*)gEve->GetCurrentEvent();
276 delphesDisplay_->DestroyEventRPhi();
277 delphesDisplay_->ImportEventRPhi(top);
278 delphesDisplay_->DestroyEventRhoZ();
279 delphesDisplay_->ImportEventRhoZ(top);
280 update_html_summary();
281 plotSummary_->FillEvent();
282 plotSummary_->Draw();
283
284 gEve->Redraw3D(kFALSE, kTRUE);
285 fStatusBar_->SetText(Form("Loaded event %d.", event_id_), 1);
286 gSystem->ProcessEvents();
287}
288
289void DelphesEventDisplay::update_html_summary()
290{
291 // Update summary of current event.
292
293 TEveElement::List_i i;
294 TEveElement::List_i j;
295 Int_t k;
296 TEveElement *el;
297 DelphesHtmlObjTable *table;
298 TEveEventManager *mgr = gEve ? gEve->GetCurrentEvent() : 0;
299 if (mgr) {
300 htmlSummary_->Clear("D");
301 for (i=mgr->BeginChildren(); i!=mgr->EndChildren(); ++i) {
302 el = ((TEveElement*)(*i));
303 if (el->IsA() == TEvePointSet::Class()) {
304 TEvePointSet *ps = (TEvePointSet *)el;
305 TString ename = ps->GetElementName();
306 TString etitle = ps->GetElementTitle();
307 if (ename.First('\'') != kNPOS)
308 ename.Remove(ename.First('\''));
309 etitle.Remove(0, 2);
310 Int_t nel = atoi(etitle.Data());
311 table = htmlSummary_->AddTable(ename, 0, nel);
312 }
313 else if (el->IsA() == TEveTrackList::Class()) {
314 TEveTrackList *tracks = (TEveTrackList *)el;
315 TString ename = tracks->GetElementName();
316 if (ename.First('\'') != kNPOS)
317 ename.Remove(ename.First('\''));
318 table = htmlSummary_->AddTable(ename.Data(), 5,
319 tracks->NumChildren(), kTRUE, "first");
320 table->SetLabel(0, "Momentum");
321 table->SetLabel(1, "P_t");
322 table->SetLabel(2, "Phi");
323 table->SetLabel(3, "Theta");
324 table->SetLabel(4, "Eta");
325 k=0;
326 for (j=tracks->BeginChildren(); j!=tracks->EndChildren(); ++j) {
327 Float_t p = ((TEveTrack*)(*j))->GetMomentum().Mag();
328 table->SetValue(0, k, p);
329 Float_t pt = ((TEveTrack*)(*j))->GetMomentum().Perp();
330 table->SetValue(1, k, pt);
331 Float_t phi = ((TEveTrack*)(*j))->GetMomentum().Phi();
332 table->SetValue(2, k, phi);
333 Float_t theta = ((TEveTrack*)(*j))->GetMomentum().Theta();
334 table->SetValue(3, k, theta);
335 Float_t eta = theta>0.0005 && theta<3.1413 ? ((TEveTrack*)(*j))->GetMomentum().Eta() : 1e10;
336 table->SetValue(4, k, eta);
337 ++k;
338 }
339 }
340 }
341 htmlSummary_->Build();
342 gHtml_->Clear();
343 gHtml_->ParseText((char*)htmlSummary_->Html().Data());
344 gHtml_->Layout();
345 }
346
347}
348
349/******************************************************************************/
350// GUI
351/******************************************************************************/
352
353void DelphesEventDisplay::make_gui()
354{
355 // Create minimal GUI for event navigation.
356
357 // add a tab on the left
358 TEveBrowser* browser = gEve->GetBrowser();
359 browser->SetWindowName("Delphes Event Display");
360 browser->StartEmbedding(TRootBrowser::kLeft);
361
362 // set the main title
363 TGMainFrame* frmMain = new TGMainFrame(gClient->GetRoot(), 1000, 600);
364 frmMain->SetWindowName("Delphes Event Display");
365 frmMain->SetCleanup(kDeepCleanup);
366
367 // build the navigation menu
368 TString icondir;
369 if(gSystem->Getenv("ROOTSYS"))
370 icondir = Form("%s/icons/", gSystem->Getenv("ROOTSYS"));
371 if(!gSystem->OpenDirectory(icondir))
372 icondir = Form("%s/icons/", (const char*)gSystem->GetFromPipe("root-config --etcdir") );
373 TGGroupFrame* vf = new TGGroupFrame(frmMain,"Event navigation",kVerticalFrame | kFitWidth );
374 {
375 TGHorizontalFrame* hf = new TGHorizontalFrame(frmMain);
376 {
377 TGPictureButton* b = 0;
378
379 b = new TGPictureButton(hf, gClient->GetPicture(icondir+"GoBack.gif"));
380 hf->AddFrame(b, new TGLayoutHints(kLHintsLeft | kLHintsCenterY , 10, 2, 10, 10));
381 b->Connect("Clicked()", "DelphesEventDisplay", this, "Bck()");
382
383 TGNumberEntry* numberEntry = new TGNumberEntry(hf,0,9,-1,TGNumberFormat::kNESInteger, TGNumberFormat::kNEANonNegative, TGNumberFormat::kNELLimitMinMax, 0, treeReader_->GetEntries());
384 hf->AddFrame(numberEntry, new TGLayoutHints(kLHintsCenterX | kLHintsCenterY , 2, 0, 10, 10));
385 this->Connect("EventChanged(Int_t)","TGNumberEntry",numberEntry,"SetIntNumber(Long_t)");
386 numberEntry->GetNumberEntry()->Connect("TextChanged(char*)", "DelphesEventDisplay", this, "PreSetEv(char*)");
387 numberEntry->GetNumberEntry()->Connect("ReturnPressed()", "DelphesEventDisplay", this, "GoTo()");
388
389 b = new TGPictureButton(hf, gClient->GetPicture(icondir+"GoForward.gif"));
390 hf->AddFrame(b, new TGLayoutHints(kLHintsRight | kLHintsCenterY , 2, 10, 10, 10));
391 b->Connect("Clicked()", "DelphesEventDisplay", this, "Fwd()");
392
393 }
394 vf->AddFrame(hf, new TGLayoutHints(kLHintsExpandX , 2, 2, 2, 2));
395
396 TGHProgressBar* progress = new TGHProgressBar(frmMain, TGProgressBar::kFancy, 100);
397 progress->SetMax( treeReader_->GetEntries());
398 progress->ShowPosition(kTRUE, kFALSE, "Event %.0f");
399 progress->SetBarColor("green");
400 vf->AddFrame(progress, new TGLayoutHints(kLHintsExpandX, 10, 10, 5, 5));
401 this->Connect("EventChanged(Int_t)","TGHProgressBar",progress,"SetPosition(Float_t)");
402 }
403 frmMain->AddFrame(vf, new TGLayoutHints(kLHintsExpandX , 5, 5, 5, 5));
404 vf = new TGGroupFrame(frmMain,"Batch operations",kVerticalFrame | kFitWidth );
405 {
406 TGTextButton *b = new TGTextButton(vf, "Initialize Summary Plots");
407 vf->AddFrame(b, new TGLayoutHints(kLHintsCenterX | kLHintsCenterY | kLHintsExpandX, 10, 10, 10, 10));
408 b->Connect("Clicked()", "DelphesEventDisplay", this, "InitSummaryPlots()");
409 }
410 frmMain->AddFrame(vf, new TGLayoutHints(kLHintsExpandX , 5, 5, 5, 5));
411
412 frmMain->MapSubwindows();
413 frmMain->Resize();
414 frmMain->MapWindow();
415 browser->StopEmbedding();
416 browser->SetTabTitle("Event Control", 0);
417
418 // the summary tab
419 htmlSummary_ = new DelphesHtmlSummary("Delphes Event Display Summary Table");
420 TEveWindowSlot* slot = TEveWindow::CreateWindowInTab(gEve->GetBrowser()->GetTabRight());
421 gHtml_ = new TGHtml(0, 100, 100);
422 TEveWindowFrame *wf = slot->MakeFrame(gHtml_);
423 gHtml_->MapSubwindows();
424 wf->SetElementName("Summary tables");
425
426 // plot tab
427 slot = TEveWindow::CreateWindowInTab(gEve->GetBrowser()->GetTabRight());
428 TEveWindowTab* tab = slot->MakeTab();
429 tab->SetElementName("Summary plots");
430 tab->SetShowTitleBar(kFALSE);
431 plotSummary_ = new DelphesPlotSummary(tab);
432 plotSummary_->Init(elements_);
433 plotSummary_->Connect("Progress(Int_t)", "DelphesEventDisplay", this, "DisplayProgress(Int_t)");
434
435}
436
437
Note: See TracBrowser for help on using the repository browser.