1 | #include <cassert>
|
---|
2 | #include <utility>
|
---|
3 | #include <algorithm>
|
---|
4 | #include "TGeoManager.h"
|
---|
5 | #include "TGeoVolume.h"
|
---|
6 | #include "external/ExRootAnalysis/ExRootConfReader.h"
|
---|
7 | #include "external/ExRootAnalysis/ExRootTreeReader.h"
|
---|
8 | #include "display/DelphesCaloData.h"
|
---|
9 | #include "display/DelphesBranchElement.h"
|
---|
10 | #include "display/Delphes3DGeometry.h"
|
---|
11 | #include "display/DelphesEventDisplay.h"
|
---|
12 | #include "classes/DelphesClasses.h"
|
---|
13 | #include "TEveElement.h"
|
---|
14 | #include "TEveJetCone.h"
|
---|
15 | #include "TEveTrack.h"
|
---|
16 | #include "TEveTrackPropagator.h"
|
---|
17 | #include "TEveCalo.h"
|
---|
18 | #include "TEveManager.h"
|
---|
19 | #include "TEveGeoNode.h"
|
---|
20 | #include "TEveTrans.h"
|
---|
21 | #include "TEveViewer.h"
|
---|
22 | #include "TEveBrowser.h"
|
---|
23 | #include "TEveArrow.h"
|
---|
24 | #include "TMath.h"
|
---|
25 | #include "TSystem.h"
|
---|
26 | #include "TRootBrowser.h"
|
---|
27 | #include "TGButton.h"
|
---|
28 | #include "TClonesArray.h"
|
---|
29 |
|
---|
30 | DelphesEventDisplay::DelphesEventDisplay()
|
---|
31 | {
|
---|
32 | event_id = 0;
|
---|
33 | tkRadius_ = 1.29;
|
---|
34 | totRadius_ = 2.0;
|
---|
35 | tkHalfLength_ = 3.0;
|
---|
36 | bz_ = 3.8;
|
---|
37 | chain_ = new TChain("Delphes");
|
---|
38 | gTreeReader = 0;
|
---|
39 | gDelphesDisplay = 0;
|
---|
40 | }
|
---|
41 |
|
---|
42 | DelphesEventDisplay::~DelphesEventDisplay()
|
---|
43 | {
|
---|
44 | delete chain_;
|
---|
45 | }
|
---|
46 |
|
---|
47 | DelphesEventDisplay::DelphesEventDisplay(const char *configFile, const char *inputFile, Delphes3DGeometry& det3D)
|
---|
48 | {
|
---|
49 | event_id = 0;
|
---|
50 | tkRadius_ = 1.29;
|
---|
51 | totRadius_ = 2.0;
|
---|
52 | tkHalfLength_ = 3.0;
|
---|
53 | bz_ = 3.8;
|
---|
54 | chain_ = new TChain("Delphes");
|
---|
55 | gTreeReader = 0;
|
---|
56 | gDelphesDisplay = 0;
|
---|
57 |
|
---|
58 | // initialize the application
|
---|
59 | TEveManager::Create(kTRUE, "IV");
|
---|
60 | TGeoManager* geom = gGeoManager;
|
---|
61 |
|
---|
62 | // build the detector
|
---|
63 | tkRadius_ = det3D.getTrackerRadius();
|
---|
64 | totRadius_ = det3D.getDetectorRadius();
|
---|
65 | tkHalfLength_ = det3D.getTrackerHalfLength();
|
---|
66 | bz_ = det3D.getBField();
|
---|
67 |
|
---|
68 | //TODO specific to some classical detector... could use better the det3D
|
---|
69 | TGeoVolume* top = det3D.getDetector(false);
|
---|
70 | geom->SetTopVolume(top);
|
---|
71 | TEveElementList *geometry = new TEveElementList("Geometry");
|
---|
72 | TEveGeoTopNode* trk = new TEveGeoTopNode(gGeoManager, top->FindNode("tracker_1"));
|
---|
73 | trk->SetVisLevel(6);
|
---|
74 | geometry->AddElement(trk);
|
---|
75 | TEveGeoTopNode* calo = new TEveGeoTopNode(gGeoManager, top->FindNode("Calorimeter_barrel_1"));
|
---|
76 | calo->SetVisLevel(3);
|
---|
77 | geometry->AddElement(calo);
|
---|
78 | calo = new TEveGeoTopNode(gGeoManager, top->FindNode("Calorimeter_endcap_1"));
|
---|
79 | calo->SetVisLevel(3);
|
---|
80 | calo->UseNodeTrans();
|
---|
81 | geometry->AddElement(calo);
|
---|
82 | calo = new TEveGeoTopNode(gGeoManager, top->FindNode("Calorimeter_endcap_2"));
|
---|
83 | calo->SetVisLevel(3);
|
---|
84 | calo->UseNodeTrans();
|
---|
85 | geometry->AddElement(calo);
|
---|
86 | TEveGeoTopNode* muon = new TEveGeoTopNode(gGeoManager, top->FindNode("muons_barrel_1"));
|
---|
87 | muon->SetVisLevel(4);
|
---|
88 | geometry->AddElement(muon);
|
---|
89 | muon = new TEveGeoTopNode(gGeoManager, top->FindNode("muons_endcap_1"));
|
---|
90 | muon->SetVisLevel(4);
|
---|
91 | muon->UseNodeTrans();
|
---|
92 | geometry->AddElement(muon);
|
---|
93 | muon = new TEveGeoTopNode(gGeoManager, top->FindNode("muons_endcap_2"));
|
---|
94 | muon->SetVisLevel(4);
|
---|
95 | muon->UseNodeTrans();
|
---|
96 | geometry->AddElement(muon);
|
---|
97 | //gGeoManager->DefaultColors();
|
---|
98 |
|
---|
99 | // Create chain of root trees
|
---|
100 | chain_->Add(inputFile);
|
---|
101 |
|
---|
102 | // Create object of class ExRootTreeReader
|
---|
103 | printf("*** Opening Delphes data file ***\n");
|
---|
104 | gTreeReader = new ExRootTreeReader(chain_);
|
---|
105 |
|
---|
106 | // prepare data collections
|
---|
107 | readConfig(configFile, det3D, gElements, gArrays);
|
---|
108 |
|
---|
109 | // viewers and scenes
|
---|
110 | gDelphesDisplay = new DelphesDisplay;
|
---|
111 | gEve->AddGlobalElement(geometry);
|
---|
112 | gDelphesDisplay->ImportGeomRPhi(geometry);
|
---|
113 | gDelphesDisplay->ImportGeomRhoZ(geometry);
|
---|
114 | // find the first calo data and use that to initialize the calo display
|
---|
115 | for(std::vector<DelphesBranchBase*>::iterator data=gElements.begin();data<gElements.end();++data) {
|
---|
116 | if(TString((*data)->GetType())=="tower") {
|
---|
117 | DelphesCaloData* container = dynamic_cast<DelphesBranchElement<DelphesCaloData>*>((*data))->GetContainer();
|
---|
118 | assert(container);
|
---|
119 | TEveCalo3D *calo3d = new TEveCalo3D(container);
|
---|
120 | calo3d->SetBarrelRadius(tkRadius_);
|
---|
121 | calo3d->SetEndCapPos(tkHalfLength_);
|
---|
122 | gEve->AddGlobalElement(calo3d);
|
---|
123 | gDelphesDisplay->ImportCaloRPhi(calo3d);
|
---|
124 | gDelphesDisplay->ImportCaloRhoZ(calo3d);
|
---|
125 | TEveCaloLego *lego = new TEveCaloLego(container);
|
---|
126 | lego->InitMainTrans();
|
---|
127 | lego->RefMainTrans().SetScale(TMath::TwoPi(), TMath::TwoPi(), TMath::Pi());
|
---|
128 | lego->SetAutoRebin(kFALSE);
|
---|
129 | lego->Set2DMode(TEveCaloLego::kValSizeOutline);
|
---|
130 | gDelphesDisplay->ImportCaloLego(lego);
|
---|
131 | break;
|
---|
132 | }
|
---|
133 | }
|
---|
134 |
|
---|
135 | make_gui();
|
---|
136 | load_event();
|
---|
137 | gEve->Redraw3D(kTRUE);
|
---|
138 |
|
---|
139 | }
|
---|
140 | // function that parses the config to extract the branches of interest and prepare containers
|
---|
141 | void DelphesEventDisplay::readConfig(const char *configFile, Delphes3DGeometry& det3D, std::vector<DelphesBranchBase*>& elements, std::vector<TClonesArray*>& arrays) {
|
---|
142 | ExRootConfReader *confReader = new ExRootConfReader;
|
---|
143 | confReader->ReadFile(configFile);
|
---|
144 | Double_t tk_radius = det3D.getTrackerRadius();
|
---|
145 | Double_t tk_length = det3D.getTrackerHalfLength();
|
---|
146 | Double_t tk_Bz = det3D.getBField();
|
---|
147 | Double_t mu_radius = det3D.getDetectorRadius();
|
---|
148 | Double_t mu_length = det3D.getDetectorHalfLength();
|
---|
149 | TAxis* etaAxis = det3D.getCaloAxes().first;
|
---|
150 | TAxis* phiAxis = det3D.getCaloAxes().second;
|
---|
151 | ExRootConfParam branches = confReader->GetParam("TreeWriter::Branch");
|
---|
152 | Int_t nBranches = branches.GetSize()/3;
|
---|
153 | DelphesBranchElement<TEveTrackList>* tlist;
|
---|
154 | DelphesBranchElement<DelphesCaloData>* clist;
|
---|
155 | DelphesBranchElement<TEveElementList>* elist;
|
---|
156 | for(Int_t b = 0; b<nBranches; ++b) {
|
---|
157 | TString input = branches[b*3].GetString();
|
---|
158 | TString name = branches[b*3+1].GetString();
|
---|
159 | TString className = branches[b*3+2].GetString();
|
---|
160 | if(className=="Track") {
|
---|
161 | if(input.Contains("eflow",TString::kIgnoreCase) || name.Contains("eflow",TString::kIgnoreCase)) continue; //no eflow
|
---|
162 | tlist = new DelphesBranchElement<TEveTrackList>(name,"track",kBlue);
|
---|
163 | elements.push_back(tlist);
|
---|
164 | TEveTrackPropagator *trkProp = tlist->GetContainer()->GetPropagator();
|
---|
165 | trkProp->SetMagField(0., 0., -tk_Bz);
|
---|
166 | trkProp->SetMaxR(tk_radius);
|
---|
167 | trkProp->SetMaxZ(tk_length);
|
---|
168 | } else 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,"track",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,"track",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,"track",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 | }
|
---|
216 | //TODO one possible simplification could be to add the array to the element class.
|
---|
217 | arrays.push_back(gTreeReader->UseBranch(name));
|
---|
218 | }
|
---|
219 | }
|
---|
220 |
|
---|
221 |
|
---|
222 | //______________________________________________________________________________
|
---|
223 | void DelphesEventDisplay::load_event()
|
---|
224 | {
|
---|
225 | // Load event specified in global event_id.
|
---|
226 | // The contents of previous event are removed.
|
---|
227 |
|
---|
228 | //TODO move this to the status bar ???
|
---|
229 | printf("Loading event %d.\n", event_id);
|
---|
230 |
|
---|
231 | // clear the previous event
|
---|
232 | gEve->GetViewers()->DeleteAnnotations();
|
---|
233 | for(std::vector<DelphesBranchBase*>::iterator data=gElements.begin();data<gElements.end();++data) {
|
---|
234 | (*data)->Reset();
|
---|
235 | }
|
---|
236 |
|
---|
237 | // read the new event
|
---|
238 | delphes_read();
|
---|
239 |
|
---|
240 | // update display
|
---|
241 | TEveElement* top = (TEveElement*)gEve->GetCurrentEvent();
|
---|
242 | gDelphesDisplay->DestroyEventRPhi();
|
---|
243 | gDelphesDisplay->ImportEventRPhi(top);
|
---|
244 | gDelphesDisplay->DestroyEventRhoZ();
|
---|
245 | gDelphesDisplay->ImportEventRhoZ(top);
|
---|
246 | //update_html_summary();
|
---|
247 | gEve->Redraw3D(kFALSE, kTRUE);
|
---|
248 | }
|
---|
249 |
|
---|
250 | void DelphesEventDisplay::delphes_read()
|
---|
251 | {
|
---|
252 |
|
---|
253 | // safety
|
---|
254 | if(event_id >= gTreeReader->GetEntries() || event_id<0 ) return;
|
---|
255 |
|
---|
256 | // Load selected branches with data from specified event
|
---|
257 | gTreeReader->ReadEntry(event_id);
|
---|
258 |
|
---|
259 | // loop over selected branches, and apply the proper recipe to fill the collections.
|
---|
260 | // this is basically to loop on gArrays to fill gElements.
|
---|
261 |
|
---|
262 | //TODO: one option would be to have templated methods in the element classes. We could simply call "element.fill()"
|
---|
263 | std::vector<TClonesArray*>::iterator data = gArrays.begin();
|
---|
264 | std::vector<DelphesBranchBase*>::iterator element = gElements.begin();
|
---|
265 | std::vector<TClonesArray*>::iterator data_tracks = gArrays.begin();
|
---|
266 | std::vector<DelphesBranchBase*>::iterator element_tracks = gElements.begin();
|
---|
267 | Int_t nTracks = 0;
|
---|
268 | for(; data<gArrays.end() && element<gElements.end(); ++data, ++element) {
|
---|
269 | TString type = (*element)->GetType();
|
---|
270 | // keep the most generic track collection for the end
|
---|
271 | if(type=="track" && TString((*element)->GetClassName())=="Track" && nTracks==0) {
|
---|
272 | data_tracks = data;
|
---|
273 | element_tracks = element;
|
---|
274 | nTracks = (*data_tracks)->GetEntries();
|
---|
275 | continue;
|
---|
276 | }
|
---|
277 | // branch on the element type
|
---|
278 | if(type=="tower") delphes_read_towers(*data,*element);
|
---|
279 | else if(type=="track" || type=="photon") delphes_read_tracks(*data,*element);
|
---|
280 | else if(type=="jet") delphes_read_jets(*data,*element);
|
---|
281 | else if(type=="vector") delphes_read_vectors(*data,*element);
|
---|
282 | }
|
---|
283 | // finish whith what we consider to be the main track collection
|
---|
284 | if(nTracks>0) delphes_read_tracks(*data,*element);
|
---|
285 | }
|
---|
286 |
|
---|
287 | void DelphesEventDisplay::delphes_read_towers(TClonesArray* data, DelphesBranchBase* element) {
|
---|
288 | DelphesCaloData* container = dynamic_cast<DelphesBranchElement<DelphesCaloData>*>(element)->GetContainer();
|
---|
289 | assert(container);
|
---|
290 | // Loop over all towers
|
---|
291 | TIter itTower(data);
|
---|
292 | Tower *tower;
|
---|
293 | while((tower = (Tower *) itTower.Next()))
|
---|
294 | {
|
---|
295 | container->AddTower(tower->Edges[0], tower->Edges[1], tower->Edges[2], tower->Edges[3]);
|
---|
296 | container->FillSlice(0, tower->Eem);
|
---|
297 | container->FillSlice(1, tower->Ehad);
|
---|
298 | }
|
---|
299 | container->DataChanged();
|
---|
300 | }
|
---|
301 |
|
---|
302 | void DelphesEventDisplay::delphes_read_tracks(TClonesArray* data, DelphesBranchBase* element) {
|
---|
303 | TEveTrackList* container = dynamic_cast<DelphesBranchElement<TEveTrackList>*>(element)->GetContainer();
|
---|
304 | assert(container);
|
---|
305 | TString className = element->GetClassName();
|
---|
306 | TIter itTrack(data);
|
---|
307 | Int_t counter = 0;
|
---|
308 | TEveTrack *eveTrack;
|
---|
309 | TEveTrackPropagator *trkProp = container->GetPropagator();
|
---|
310 | if(className=="Track") {
|
---|
311 | // Loop over all tracks
|
---|
312 | Track *track;
|
---|
313 | while((track = (Track *) itTrack.Next())) {
|
---|
314 | TParticle pb(track->PID, 1, 0, 0, 0, 0,
|
---|
315 | track->P4().Px(), track->P4().Py(),
|
---|
316 | track->P4().Pz(), track->P4().E(),
|
---|
317 | track->X, track->Y, track->Z, 0.0);
|
---|
318 |
|
---|
319 | eveTrack = new TEveTrack(&pb, counter, trkProp);
|
---|
320 | eveTrack->SetName(Form("%s [%d]", pb.GetName(), counter++));
|
---|
321 | eveTrack->SetStdTitle();
|
---|
322 | eveTrack->SetAttLineAttMarker(container);
|
---|
323 | container->AddElement(eveTrack);
|
---|
324 | eveTrack->SetLineColor(element->GetColor());
|
---|
325 | eveTrack->MakeTrack();
|
---|
326 | }
|
---|
327 | } else if(className=="Electron") {
|
---|
328 | // Loop over all electrons
|
---|
329 | Electron *electron;
|
---|
330 | while((electron = (Electron *) itTrack.Next())) {
|
---|
331 | TParticle pb(electron->Charge<0?11:-11, 1, 0, 0, 0, 0,
|
---|
332 | electron->P4().Px(), electron->P4().Py(),
|
---|
333 | electron->P4().Pz(), electron->P4().E(),
|
---|
334 | 0., 0., 0., 0.);
|
---|
335 |
|
---|
336 | eveTrack = new TEveTrack(&pb, counter, trkProp);
|
---|
337 | eveTrack->SetName(Form("%s [%d]", pb.GetName(), counter++));
|
---|
338 | eveTrack->SetStdTitle();
|
---|
339 | eveTrack->SetAttLineAttMarker(container);
|
---|
340 | container->AddElement(eveTrack);
|
---|
341 | eveTrack->SetLineColor(element->GetColor());
|
---|
342 | eveTrack->MakeTrack();
|
---|
343 | }
|
---|
344 | } else if(className=="Muon") {
|
---|
345 | // Loop over all muons
|
---|
346 | Muon *muon;
|
---|
347 | while((muon = (Muon *) itTrack.Next())) {
|
---|
348 | TParticle pb(muon->Charge<0?13:-13, 1, 0, 0, 0, 0,
|
---|
349 | muon->P4().Px(), muon->P4().Py(),
|
---|
350 | muon->P4().Pz(), muon->P4().E(),
|
---|
351 | 0., 0., 0., 0.);
|
---|
352 |
|
---|
353 | eveTrack = new TEveTrack(&pb, counter, trkProp);
|
---|
354 | eveTrack->SetName(Form("%s [%d]", pb.GetName(), counter++));
|
---|
355 | eveTrack->SetStdTitle();
|
---|
356 | eveTrack->SetAttLineAttMarker(container);
|
---|
357 | container->AddElement(eveTrack);
|
---|
358 | eveTrack->SetLineColor(element->GetColor());
|
---|
359 | eveTrack->MakeTrack();
|
---|
360 | }
|
---|
361 | } else if(className=="Photon") {
|
---|
362 | // Loop over all photons
|
---|
363 | Photon *photon;
|
---|
364 | while((photon = (Photon *) itTrack.Next())) {
|
---|
365 | TParticle pb(22, 1, 0, 0, 0, 0,
|
---|
366 | photon->P4().Px(), photon->P4().Py(),
|
---|
367 | photon->P4().Pz(), photon->P4().E(),
|
---|
368 | 0., 0., 0., 0.);
|
---|
369 |
|
---|
370 | eveTrack = new TEveTrack(&pb, counter, trkProp);
|
---|
371 | eveTrack->SetName(Form("%s [%d]", pb.GetName(), counter++));
|
---|
372 | eveTrack->SetStdTitle();
|
---|
373 | eveTrack->SetAttLineAttMarker(container);
|
---|
374 | container->AddElement(eveTrack);
|
---|
375 | eveTrack->SetLineColor(element->GetColor());
|
---|
376 | eveTrack->MakeTrack();
|
---|
377 | }
|
---|
378 | }
|
---|
379 | }
|
---|
380 |
|
---|
381 | void DelphesEventDisplay::delphes_read_jets(TClonesArray* data, DelphesBranchBase* element) {
|
---|
382 | TEveElementList* container = dynamic_cast<DelphesBranchElement<TEveElementList>*>(element)->GetContainer();
|
---|
383 | assert(container);
|
---|
384 | TIter itJet(data);
|
---|
385 | Jet *jet;
|
---|
386 | TEveJetCone *eveJetCone;
|
---|
387 | // Loop over all jets
|
---|
388 | Int_t counter = 0;
|
---|
389 | while((jet = (Jet *) itJet.Next()))
|
---|
390 | {
|
---|
391 | eveJetCone = new TEveJetCone();
|
---|
392 | eveJetCone->SetTitle(Form("jet [%d]: Pt=%f, Eta=%f, \nPhi=%f, M=%f",counter,jet->PT, jet->Eta, jet->Phi, jet->Mass));
|
---|
393 | eveJetCone->SetName(Form("jet [%d]", counter++));
|
---|
394 | eveJetCone->SetMainTransparency(60);
|
---|
395 | eveJetCone->SetLineColor(element->GetColor());
|
---|
396 | eveJetCone->SetFillColor(element->GetColor());
|
---|
397 | eveJetCone->SetCylinder(tkRadius_ - 10, tkHalfLength_ - 10);
|
---|
398 | eveJetCone->SetPickable(kTRUE);
|
---|
399 | eveJetCone->AddEllipticCone(jet->Eta, jet->Phi, jet->DeltaEta, jet->DeltaPhi);
|
---|
400 | container->AddElement(eveJetCone);
|
---|
401 | }
|
---|
402 | }
|
---|
403 |
|
---|
404 | void DelphesEventDisplay::delphes_read_vectors(TClonesArray* data, DelphesBranchBase* element) {
|
---|
405 | TEveElementList* container = dynamic_cast<DelphesBranchElement<TEveElementList>*>(element)->GetContainer();
|
---|
406 | assert(container);
|
---|
407 | TIter itMet(data);
|
---|
408 | MissingET *MET;
|
---|
409 | TEveArrow *eveMet;
|
---|
410 | // Missing Et
|
---|
411 | Double_t maxPt = 50.;
|
---|
412 | // TODO to be changed as we don't have access to maxPt anymore. MET scale could be a general parameter set in GUI
|
---|
413 | while((MET = (MissingET*) itMet.Next())) {
|
---|
414 | eveMet = new TEveArrow((tkRadius_ * MET->MET/maxPt)*cos(MET->Phi), (tkRadius_ * MET->MET/maxPt)*sin(MET->Phi), 0., 0., 0., 0.);
|
---|
415 | eveMet->SetMainColor(element->GetColor());
|
---|
416 | eveMet->SetTubeR(0.04);
|
---|
417 | eveMet->SetConeR(0.08);
|
---|
418 | eveMet->SetConeL(0.10);
|
---|
419 | eveMet->SetPickable(kTRUE);
|
---|
420 | eveMet->SetName("Missing Et");
|
---|
421 | eveMet->SetTitle(Form("Missing Et (%.1f GeV)",MET->MET));
|
---|
422 | container->AddElement(eveMet);
|
---|
423 | }
|
---|
424 | }
|
---|
425 |
|
---|
426 | /******************************************************************************/
|
---|
427 | // GUI
|
---|
428 | /******************************************************************************/
|
---|
429 |
|
---|
430 | void DelphesEventDisplay::make_gui()
|
---|
431 | {
|
---|
432 | // Create minimal GUI for event navigation.
|
---|
433 | // TODO: better GUI could be made based on the ch15 of the manual (Writing a GUI)
|
---|
434 |
|
---|
435 | // add a tab on the left
|
---|
436 | TEveBrowser* browser = gEve->GetBrowser();
|
---|
437 | browser->StartEmbedding(TRootBrowser::kLeft);
|
---|
438 |
|
---|
439 | // set the main title
|
---|
440 | TGMainFrame* frmMain = new TGMainFrame(gClient->GetRoot(), 1000, 600);
|
---|
441 | frmMain->SetWindowName("Delphes Event Display");
|
---|
442 | frmMain->SetCleanup(kDeepCleanup);
|
---|
443 |
|
---|
444 | // build the navigation menu
|
---|
445 | TGHorizontalFrame* hf = new TGHorizontalFrame(frmMain);
|
---|
446 | {
|
---|
447 | TString icondir;
|
---|
448 | if(gSystem->Getenv("ROOTSYS"))
|
---|
449 | icondir = Form("%s/icons/", gSystem->Getenv("ROOTSYS"));
|
---|
450 | if(!gSystem->OpenDirectory(icondir))
|
---|
451 | icondir = Form("%s/icons/", (const char*)gSystem->GetFromPipe("root-config --etcdir") );
|
---|
452 | TGPictureButton* b = 0;
|
---|
453 | EvNavHandler *fh = new EvNavHandler(this);
|
---|
454 |
|
---|
455 | b = new TGPictureButton(hf, gClient->GetPicture(icondir+"GoBack.gif"));
|
---|
456 | hf->AddFrame(b);
|
---|
457 | b->Connect("Clicked()", "EvNavHandler", fh, "Bck()");
|
---|
458 |
|
---|
459 | b = new TGPictureButton(hf, gClient->GetPicture(icondir+"GoForward.gif"));
|
---|
460 | hf->AddFrame(b);
|
---|
461 | b->Connect("Clicked()", "EvNavHandler", fh, "Fwd()");
|
---|
462 | }
|
---|
463 | frmMain->AddFrame(hf);
|
---|
464 | frmMain->MapSubwindows();
|
---|
465 | frmMain->Resize();
|
---|
466 | frmMain->MapWindow();
|
---|
467 | browser->StopEmbedding();
|
---|
468 | browser->SetTabTitle("Event Control", 0);
|
---|
469 | }
|
---|
470 |
|
---|