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 |
|
---|
56 | DelphesEventDisplay::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 |
|
---|
71 | DelphesEventDisplay::~DelphesEventDisplay()
|
---|
72 | {
|
---|
73 | delete chain_;
|
---|
74 | }
|
---|
75 |
|
---|
76 | void 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 |
|
---|
86 | DelphesEventDisplay::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->SetAutoRebin(kFALSE);
|
---|
160 | lego->Set2DMode(TEveCaloLego::kValSizeOutline);
|
---|
161 | delphesDisplay_->ImportCaloLego(lego);
|
---|
162 | break;
|
---|
163 | }
|
---|
164 | }
|
---|
165 |
|
---|
166 | // the GUI: control panel, summary tab
|
---|
167 | make_gui();
|
---|
168 |
|
---|
169 | //ready...
|
---|
170 | fStatusBar_->SetText("Ready.", 1);
|
---|
171 | gSystem->ProcessEvents();
|
---|
172 | load_event();
|
---|
173 | gEve->Redraw3D(kTRUE);
|
---|
174 |
|
---|
175 | }
|
---|
176 |
|
---|
177 | // function that parses the config to extract the branches of interest and prepare containers
|
---|
178 | void DelphesEventDisplay::readConfig(const char *configFile, std::vector<DelphesBranchBase*>& elements) {
|
---|
179 | ExRootConfReader *confReader = new ExRootConfReader;
|
---|
180 | confReader->ReadFile(configFile);
|
---|
181 | ExRootConfParam branches = confReader->GetParam("TreeWriter::Branch");
|
---|
182 | Int_t nBranches = branches.GetSize()/3;
|
---|
183 | DelphesBranchElement<TEveTrackList>* tlist;
|
---|
184 | DelphesBranchElement<DelphesCaloData>* clist;
|
---|
185 | DelphesBranchElement<TEveElementList>* elist;
|
---|
186 | // first loop with all but tracks
|
---|
187 | for(Int_t b = 0; b<nBranches; ++b) {
|
---|
188 | TString input = branches[b*3].GetString();
|
---|
189 | TString name = branches[b*3+1].GetString();
|
---|
190 | TString className = branches[b*3+2].GetString();
|
---|
191 | if(className=="Tower") {
|
---|
192 | if(input.Contains("eflow",TString::kIgnoreCase) || name.Contains("eflow",TString::kIgnoreCase)) continue; //no eflow
|
---|
193 | clist = new DelphesBranchElement<DelphesCaloData>(name,treeReader_->UseBranch(name),kBlack);
|
---|
194 | clist->GetContainer()->SetEtaBins(etaAxis_);
|
---|
195 | clist->GetContainer()->SetPhiBins(phiAxis_);
|
---|
196 | elements.push_back(clist);
|
---|
197 | } else if(className=="Jet") {
|
---|
198 | if(input.Contains("GenJetFinder")) {
|
---|
199 | elist = new DelphesBranchElement<TEveElementList>(name,treeReader_->UseBranch(name),kCyan);
|
---|
200 | elist->GetContainer()->SetRnrSelf(false);
|
---|
201 | elist->GetContainer()->SetRnrChildren(false);
|
---|
202 | elist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
|
---|
203 | elements.push_back(elist);
|
---|
204 | } else {
|
---|
205 | elist = new DelphesBranchElement<TEveElementList>(name,treeReader_->UseBranch(name),kYellow);
|
---|
206 | elist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
|
---|
207 | elements.push_back(elist);
|
---|
208 | }
|
---|
209 | } else if(className=="Electron") {
|
---|
210 | tlist = new DelphesBranchElement<TEveTrackList>(name,treeReader_->UseBranch(name),kRed);
|
---|
211 | tlist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
|
---|
212 | elements.push_back(tlist);
|
---|
213 | } else if(className=="Photon") {
|
---|
214 | tlist = new DelphesBranchElement<TEveTrackList>(name,treeReader_->UseBranch(name),kYellow);
|
---|
215 | tlist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
|
---|
216 | elements.push_back(tlist);
|
---|
217 | } else if(className=="Muon") {
|
---|
218 | tlist = new DelphesBranchElement<TEveTrackList>(name,treeReader_->UseBranch(name),kGreen);
|
---|
219 | tlist->SetTrackingVolume(totRadius_, muHalfLength_, bz_);
|
---|
220 | elements.push_back(tlist);
|
---|
221 | } else if(className=="MissingET") {
|
---|
222 | elist = new DelphesBranchElement<TEveElementList>(name,treeReader_->UseBranch(name),kViolet);
|
---|
223 | elist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
|
---|
224 | elements.push_back(elist);
|
---|
225 | } else if(className=="GenParticle") {
|
---|
226 | tlist = new DelphesBranchElement<TEveTrackList>(name,treeReader_->UseBranch(name),kCyan);
|
---|
227 | tlist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
|
---|
228 | tlist->GetContainer()->SetRnrSelf(false);
|
---|
229 | tlist->GetContainer()->SetRnrChildren(false);
|
---|
230 | elements.push_back(tlist);
|
---|
231 | } else {
|
---|
232 | continue;
|
---|
233 | }
|
---|
234 | }
|
---|
235 | // second loop for tracks
|
---|
236 | for(Int_t b = 0; b<nBranches; ++b) {
|
---|
237 | TString input = branches[b*3].GetString();
|
---|
238 | TString name = branches[b*3+1].GetString();
|
---|
239 | TString className = branches[b*3+2].GetString();
|
---|
240 | if(className=="Track") {
|
---|
241 | if(input.Contains("eflow",TString::kIgnoreCase) || name.Contains("eflow",TString::kIgnoreCase)) continue; //no eflow
|
---|
242 | tlist = new DelphesBranchElement<TEveTrackList>(name,treeReader_->UseBranch(name),kBlue);
|
---|
243 | tlist->SetTrackingVolume(tkRadius_, tkHalfLength_, bz_);
|
---|
244 | elements.push_back(tlist);
|
---|
245 | }
|
---|
246 | }
|
---|
247 | }
|
---|
248 |
|
---|
249 | void DelphesEventDisplay::load_event()
|
---|
250 | {
|
---|
251 | // Load event specified in global event_id_.
|
---|
252 | // The contents of previous event are removed.
|
---|
253 |
|
---|
254 | // safety
|
---|
255 | if(event_id_ >= treeReader_->GetEntries() || event_id_<0 ) return;
|
---|
256 |
|
---|
257 | // message
|
---|
258 | fStatusBar_->SetText(Form("Loading event %d.", event_id_), 1);
|
---|
259 | gSystem->ProcessEvents();
|
---|
260 |
|
---|
261 | // clear the previous event
|
---|
262 | gEve->GetViewers()->DeleteAnnotations();
|
---|
263 | for(std::vector<DelphesBranchBase*>::iterator data=elements_.begin();data<elements_.end();++data) {
|
---|
264 | (*data)->Reset();
|
---|
265 | }
|
---|
266 |
|
---|
267 | // Load selected branches with data from specified event
|
---|
268 | treeReader_->ReadEntry(event_id_);
|
---|
269 | for(std::vector<DelphesBranchBase*>::iterator data=elements_.begin();data<elements_.end();++data) {
|
---|
270 | (*data)->ReadBranch();
|
---|
271 | }
|
---|
272 |
|
---|
273 | // update display
|
---|
274 | TEveElement* top = (TEveElement*)gEve->GetCurrentEvent();
|
---|
275 | delphesDisplay_->DestroyEventRPhi();
|
---|
276 | delphesDisplay_->ImportEventRPhi(top);
|
---|
277 | delphesDisplay_->DestroyEventRhoZ();
|
---|
278 | delphesDisplay_->ImportEventRhoZ(top);
|
---|
279 | update_html_summary();
|
---|
280 | plotSummary_->FillEvent();
|
---|
281 | plotSummary_->Draw();
|
---|
282 |
|
---|
283 | gEve->Redraw3D(kFALSE, kTRUE);
|
---|
284 | fStatusBar_->SetText(Form("Loaded event %d.", event_id_), 1);
|
---|
285 | gSystem->ProcessEvents();
|
---|
286 | }
|
---|
287 |
|
---|
288 | void DelphesEventDisplay::update_html_summary()
|
---|
289 | {
|
---|
290 | // Update summary of current event.
|
---|
291 |
|
---|
292 | TEveElement::List_i i;
|
---|
293 | TEveElement::List_i j;
|
---|
294 | Int_t k;
|
---|
295 | TEveElement *el;
|
---|
296 | DelphesHtmlObjTable *table;
|
---|
297 | TEveEventManager *mgr = gEve ? gEve->GetCurrentEvent() : 0;
|
---|
298 | if (mgr) {
|
---|
299 | htmlSummary_->Clear("D");
|
---|
300 | for (i=mgr->BeginChildren(); i!=mgr->EndChildren(); ++i) {
|
---|
301 | el = ((TEveElement*)(*i));
|
---|
302 | if (el->IsA() == TEvePointSet::Class()) {
|
---|
303 | TEvePointSet *ps = (TEvePointSet *)el;
|
---|
304 | TString ename = ps->GetElementName();
|
---|
305 | TString etitle = ps->GetElementTitle();
|
---|
306 | if (ename.First('\'') != kNPOS)
|
---|
307 | ename.Remove(ename.First('\''));
|
---|
308 | etitle.Remove(0, 2);
|
---|
309 | Int_t nel = atoi(etitle.Data());
|
---|
310 | table = htmlSummary_->AddTable(ename, 0, nel);
|
---|
311 | }
|
---|
312 | else if (el->IsA() == TEveTrackList::Class()) {
|
---|
313 | TEveTrackList *tracks = (TEveTrackList *)el;
|
---|
314 | TString ename = tracks->GetElementName();
|
---|
315 | if (ename.First('\'') != kNPOS)
|
---|
316 | ename.Remove(ename.First('\''));
|
---|
317 | table = htmlSummary_->AddTable(ename.Data(), 5,
|
---|
318 | tracks->NumChildren(), kTRUE, "first");
|
---|
319 | table->SetLabel(0, "Momentum");
|
---|
320 | table->SetLabel(1, "P_t");
|
---|
321 | table->SetLabel(2, "Phi");
|
---|
322 | table->SetLabel(3, "Theta");
|
---|
323 | table->SetLabel(4, "Eta");
|
---|
324 | k=0;
|
---|
325 | for (j=tracks->BeginChildren(); j!=tracks->EndChildren(); ++j) {
|
---|
326 | Float_t p = ((TEveTrack*)(*j))->GetMomentum().Mag();
|
---|
327 | table->SetValue(0, k, p);
|
---|
328 | Float_t pt = ((TEveTrack*)(*j))->GetMomentum().Perp();
|
---|
329 | table->SetValue(1, k, pt);
|
---|
330 | Float_t phi = ((TEveTrack*)(*j))->GetMomentum().Phi();
|
---|
331 | table->SetValue(2, k, phi);
|
---|
332 | Float_t theta = ((TEveTrack*)(*j))->GetMomentum().Theta();
|
---|
333 | table->SetValue(3, k, theta);
|
---|
334 | Float_t eta = theta>0.0005 && theta<3.1413 ? ((TEveTrack*)(*j))->GetMomentum().Eta() : 1e10;
|
---|
335 | table->SetValue(4, k, eta);
|
---|
336 | ++k;
|
---|
337 | }
|
---|
338 | }
|
---|
339 | }
|
---|
340 | htmlSummary_->Build();
|
---|
341 | gHtml_->Clear();
|
---|
342 | gHtml_->ParseText((char*)htmlSummary_->Html().Data());
|
---|
343 | gHtml_->Layout();
|
---|
344 | }
|
---|
345 |
|
---|
346 | }
|
---|
347 |
|
---|
348 | /******************************************************************************/
|
---|
349 | // GUI
|
---|
350 | /******************************************************************************/
|
---|
351 |
|
---|
352 | void DelphesEventDisplay::make_gui()
|
---|
353 | {
|
---|
354 | // Create minimal GUI for event navigation.
|
---|
355 |
|
---|
356 | // add a tab on the left
|
---|
357 | TEveBrowser* browser = gEve->GetBrowser();
|
---|
358 | browser->SetWindowName("Delphes Event Display");
|
---|
359 | browser->StartEmbedding(TRootBrowser::kLeft);
|
---|
360 |
|
---|
361 | // set the main title
|
---|
362 | TGMainFrame* frmMain = new TGMainFrame(gClient->GetRoot(), 1000, 600);
|
---|
363 | frmMain->SetWindowName("Delphes Event Display");
|
---|
364 | frmMain->SetCleanup(kDeepCleanup);
|
---|
365 |
|
---|
366 | // build the navigation menu
|
---|
367 | TString icondir;
|
---|
368 | if(gSystem->Getenv("ROOTSYS"))
|
---|
369 | icondir = Form("%s/icons/", gSystem->Getenv("ROOTSYS"));
|
---|
370 | if(!gSystem->OpenDirectory(icondir))
|
---|
371 | icondir = Form("%s/icons/", (const char*)gSystem->GetFromPipe("root-config --etcdir") );
|
---|
372 | TGGroupFrame* vf = new TGGroupFrame(frmMain,"Event navigation",kVerticalFrame | kFitWidth );
|
---|
373 | {
|
---|
374 | TGHorizontalFrame* hf = new TGHorizontalFrame(frmMain);
|
---|
375 | {
|
---|
376 | TGPictureButton* b = 0;
|
---|
377 |
|
---|
378 | b = new TGPictureButton(hf, gClient->GetPicture(icondir+"GoBack.gif"));
|
---|
379 | hf->AddFrame(b, new TGLayoutHints(kLHintsLeft | kLHintsCenterY , 10, 2, 10, 10));
|
---|
380 | b->Connect("Clicked()", "DelphesEventDisplay", this, "Bck()");
|
---|
381 |
|
---|
382 | TGNumberEntry* numberEntry = new TGNumberEntry(hf,0,9,-1,TGNumberFormat::kNESInteger, TGNumberFormat::kNEANonNegative, TGNumberFormat::kNELLimitMinMax, 0, treeReader_->GetEntries());
|
---|
383 | hf->AddFrame(numberEntry, new TGLayoutHints(kLHintsCenterX | kLHintsCenterY , 2, 0, 10, 10));
|
---|
384 | this->Connect("EventChanged(Int_t)","TGNumberEntry",numberEntry,"SetIntNumber(Long_t)");
|
---|
385 | numberEntry->GetNumberEntry()->Connect("TextChanged(char*)", "DelphesEventDisplay", this, "PreSetEv(char*)");
|
---|
386 | numberEntry->GetNumberEntry()->Connect("ReturnPressed()", "DelphesEventDisplay", this, "GoTo()");
|
---|
387 |
|
---|
388 | b = new TGPictureButton(hf, gClient->GetPicture(icondir+"GoForward.gif"));
|
---|
389 | hf->AddFrame(b, new TGLayoutHints(kLHintsRight | kLHintsCenterY , 2, 10, 10, 10));
|
---|
390 | b->Connect("Clicked()", "DelphesEventDisplay", this, "Fwd()");
|
---|
391 |
|
---|
392 | }
|
---|
393 | vf->AddFrame(hf, new TGLayoutHints(kLHintsExpandX , 2, 2, 2, 2));
|
---|
394 |
|
---|
395 | TGHProgressBar* progress = new TGHProgressBar(frmMain, TGProgressBar::kFancy, 100);
|
---|
396 | progress->SetMax( treeReader_->GetEntries());
|
---|
397 | progress->ShowPosition(kTRUE, kFALSE, "Event %.0f");
|
---|
398 | progress->SetBarColor("green");
|
---|
399 | vf->AddFrame(progress, new TGLayoutHints(kLHintsExpandX, 10, 10, 5, 5));
|
---|
400 | this->Connect("EventChanged(Int_t)","TGHProgressBar",progress,"SetPosition(Float_t)");
|
---|
401 | }
|
---|
402 | frmMain->AddFrame(vf, new TGLayoutHints(kLHintsExpandX , 5, 5, 5, 5));
|
---|
403 | vf = new TGGroupFrame(frmMain,"Batch operations",kVerticalFrame | kFitWidth );
|
---|
404 | {
|
---|
405 | TGTextButton *b = new TGTextButton(vf, "Initialize Summary Plots");
|
---|
406 | vf->AddFrame(b, new TGLayoutHints(kLHintsCenterX | kLHintsCenterY | kLHintsExpandX, 10, 10, 10, 10));
|
---|
407 | b->Connect("Clicked()", "DelphesEventDisplay", this, "InitSummaryPlots()");
|
---|
408 | }
|
---|
409 | frmMain->AddFrame(vf, new TGLayoutHints(kLHintsExpandX , 5, 5, 5, 5));
|
---|
410 |
|
---|
411 | frmMain->MapSubwindows();
|
---|
412 | frmMain->Resize();
|
---|
413 | frmMain->MapWindow();
|
---|
414 | browser->StopEmbedding();
|
---|
415 | browser->SetTabTitle("Event Control", 0);
|
---|
416 |
|
---|
417 | // the summary tab
|
---|
418 | htmlSummary_ = new DelphesHtmlSummary("Delphes Event Display Summary Table");
|
---|
419 | TEveWindowSlot* slot = TEveWindow::CreateWindowInTab(gEve->GetBrowser()->GetTabRight());
|
---|
420 | gHtml_ = new TGHtml(0, 100, 100);
|
---|
421 | TEveWindowFrame *wf = slot->MakeFrame(gHtml_);
|
---|
422 | gHtml_->MapSubwindows();
|
---|
423 | wf->SetElementName("Summary tables");
|
---|
424 |
|
---|
425 | // plot tab
|
---|
426 | slot = TEveWindow::CreateWindowInTab(gEve->GetBrowser()->GetTabRight());
|
---|
427 | TEveWindowTab* tab = slot->MakeTab();
|
---|
428 | tab->SetElementName("Summary plots");
|
---|
429 | tab->SetShowTitleBar(kFALSE);
|
---|
430 | plotSummary_ = new DelphesPlotSummary(tab);
|
---|
431 | plotSummary_->Init(elements_);
|
---|
432 | }
|
---|
433 |
|
---|