Fork me on GitHub

source: git/display/DelphesEventDisplay.cc@ 8b04b31

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

Working summary table

Not all collections are in. To be checked why.

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