3d94b490 |
1 | // $Id$ |
2 | // Author: Stefano Carrazza 2010 |
3 | |
4 | /************************************************************************** |
5 | * Copyright(c) 1998-2009, ALICE Experiment at CERN, all rights reserved. * |
6 | * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for * |
7 | * full copyright notice. * |
8 | **************************************************************************/ |
9 | |
10 | #include "AliESDEvent.h" |
11 | #include "AliEveLego.h" |
12 | #include "AliEveEventManager.h" |
13 | #include "AliEveMultiView.h" |
14 | #include "AliPhysicsSelection.h" |
15 | #include "AliEveEventSelector.h" |
16 | |
17 | #include "TH2F.h" |
18 | #include "TMath.h" |
19 | #include "TGLViewer.h" |
20 | #include "TEveWindow.h" |
21 | #include "TEveManager.h" |
22 | #include "TEveBrowser.h" |
23 | #include "TEveViewer.h" |
24 | #include "TEveScene.h" |
25 | #include "TEveCaloLegoOverlay.h" |
26 | #include "TEveCalo.h" |
27 | #include "TEveCaloData.h" |
28 | #include "TEveLegoEventHandler.h" |
29 | #include "TEveTrans.h" |
30 | #include "TEveProjectionManager.h" |
31 | #include "TEveProjectionAxes.h" |
32 | #include "TGLWidget.h" |
33 | #include "TGLOverlayButton.h" |
a096a7bf |
34 | #include "TStopwatch.h" |
3d94b490 |
35 | |
36 | //______________________________________________________________________________ |
37 | // Allow histograms visualization in 2D and 3D. |
38 | // |
39 | |
40 | ClassImp(AliEveLego) |
41 | Double_t pi = TMath::Pi(); |
42 | |
43 | //______________________________________________________________________________ |
44 | AliEveLego::AliEveLego(const char* name) : |
45 | TEveElementList(name), |
46 | fChargeId(1), |
47 | fTracksId(1), |
48 | fEventsId(1), |
49 | fMaxPt(10000), |
50 | fChargeIdAE(0), |
51 | fTracksIdAE(0), |
a096a7bf |
52 | fMaxPtAE(10000), |
3d94b490 |
53 | fEsd(0), |
54 | fPhysicsSelection(0), |
55 | fHistopos(0), |
56 | fHistoposclone(0), |
57 | fHistopos_all_events(0), |
58 | fHistoneg(0), |
59 | fHistonegclone(0), |
60 | fHistoneg_all_events(0), |
61 | fData(0), |
62 | fData_all_events(0), |
63 | fLego(0), |
64 | fLego_all_events(0), |
65 | fCalo3d(0), |
66 | fCalo3d_all_events(0), |
67 | fGlv(0), |
68 | fHisto2d_v(0), |
69 | fHisto2d_s(0), |
70 | fHisto2d_s2(0), |
71 | fHisto2d_all_events_v0(0), |
72 | fHisto2d_all_events_v1(0), |
73 | fHisto2d_all_events_v2(0), |
74 | fHisto2d_all_events_v3(0), |
75 | fHisto2d_all_events_s0(0), |
76 | fHisto2d_all_events_s1(0), |
77 | fHisto2d_all_events_s2(0), |
78 | fHisto2d_all_events_s3(0), |
79 | fAl(0), |
80 | fHisto2d_lego_overlay(0), |
81 | fHisto2d_all_events_lego_overlay(0), |
82 | fHisto2d_all_events_slot(0), |
83 | fEventSelector(0), |
84 | fShowEventsInfo(0), |
85 | fGButton(0), |
86 | fB1(0), |
87 | fB2(0) |
88 | { |
89 | // Constructor. |
90 | gEve->AddToListTree(this,0); |
91 | |
92 | // Get Current ESD event |
93 | fEsd = AliEveEventManager::AssertESD(); |
94 | |
95 | fPhysicsSelection = new AliPhysicsSelection(); |
96 | fPhysicsSelection->Initialize(fEsd->GetRunNumber()); |
97 | |
98 | fEventSelector = AliEveEventManager::GetMaster()->GetEventSelector(); |
99 | |
100 | fHistopos = new TH2F("histopos","Histo 2d positive", 100, -1.5, 1.5, 80, -pi, pi); |
101 | fHistoneg = new TH2F("histoneg","Histo 2d negative", 100, -1.5, 1.5, 80, -pi, pi); |
a096a7bf |
102 | |
103 | fHistopos->SetDirectory(0); |
104 | fHistoneg->SetDirectory(0); |
3d94b490 |
105 | |
106 | fData = new TEveCaloDataHist(); |
107 | fData->AddHistogram(fHistoneg); |
108 | fData->RefSliceInfo(0).Setup("NegCg:", 0, kBlue); |
109 | fData->AddHistogram(fHistopos); |
110 | fData->RefSliceInfo(1).Setup("PosCg:", 0, kRed); |
111 | fData->GetEtaBins()->SetTitleFont(120); |
112 | fData->GetEtaBins()->SetTitle("h"); |
113 | fData->GetPhiBins()->SetTitleFont(120); |
114 | fData->GetPhiBins()->SetTitle("f"); |
115 | fData->IncDenyDestroy(); |
116 | |
117 | fCalo3d = new TEveCalo3D(fData); |
118 | fCalo3d->SetBarrelRadius(550); |
119 | fCalo3d->SetEndCapPos(550); |
120 | |
121 | // plotting histo |
122 | fLego = new TEveCaloLego(fData); |
123 | |
124 | // projections |
125 | fAl = AliEveMultiView::Instance(); |
126 | fAl->ImportEventRPhi(fCalo3d); |
127 | fAl->ImportEventRhoZ(fCalo3d); |
128 | |
129 | // glbutton |
130 | fGButton = new TGLOverlayButton(0, "", 10.0, -10.0, 190.0, 20.0); |
131 | fGButton->SetAlphaValues(1.5,1.5); |
132 | fB1 = new TGLOverlayButton(0, "", 10.0, -30.0, 190.0, 20.0); |
133 | fB1->SetAlphaValues(1.5,1.5); |
134 | fB2 = new TGLOverlayButton(0, "", 10.0, -50.0, 190.0, 20.0); |
135 | fB2->SetAlphaValues(1.5,1.5); |
136 | |
137 | // Update |
138 | Update(); |
139 | } |
140 | |
141 | //______________________________________________________________________________ |
142 | AliEveLego::~AliEveLego() |
143 | { |
144 | delete fEsd; |
145 | delete fPhysicsSelection; |
146 | delete fHistopos; |
147 | delete fHistopos_all_events; |
148 | delete fHistoneg; |
149 | delete fHistoneg_all_events; |
150 | |
151 | delete fData; |
152 | delete fData_all_events; |
153 | delete fLego; |
154 | delete fLego_all_events; |
155 | delete fCalo3d; |
156 | delete fCalo3d_all_events; |
157 | delete fGlv; |
158 | |
159 | delete fHisto2d_v; |
160 | delete fHisto2d_s; |
161 | delete fHisto2d_s2; |
162 | delete fHisto2d_all_events_v0; |
163 | delete fHisto2d_all_events_v1; |
164 | delete fHisto2d_all_events_v2; |
165 | delete fHisto2d_all_events_v3; |
166 | delete fHisto2d_all_events_s0; |
167 | delete fHisto2d_all_events_s1; |
168 | delete fHisto2d_all_events_s2; |
169 | delete fHisto2d_all_events_s3; |
170 | |
171 | delete fAl; |
172 | delete fHisto2d_lego_overlay; |
173 | delete fHisto2d_all_events_lego_overlay; |
174 | delete fHisto2d_all_events_slot; |
175 | |
176 | delete fEventSelector; |
177 | delete fGButton; |
178 | delete fB1; |
179 | delete fB2; |
180 | } |
181 | |
182 | //______________________________________________________________________________ |
183 | Double_t getphi(Double_t phi) |
184 | { |
185 | Double_t pi = TMath::Pi(); |
186 | |
187 | if (phi > pi) { |
188 | phi -= 2*pi; |
189 | } |
190 | return phi; |
191 | } |
192 | |
193 | //______________________________________________________________________________ |
194 | TEveCaloDataHist* AliEveLego::LoadData() |
195 | { |
196 | |
197 | fHistopos->Reset(); |
198 | fHistoneg->Reset(); |
199 | |
200 | // Getting current tracks, filling histograms |
201 | for (int n = 0; n < fEsd->GetNumberOfTracks(); ++n) { |
202 | |
203 | if (fEsd->GetTrack(n)->GetSign() > 0) { |
204 | fHistopos->Fill(fEsd->GetTrack(n)->Eta(), |
205 | getphi(fEsd->GetTrack(n)->Phi()), |
206 | fabs(fEsd->GetTrack(n)->Pt())); |
207 | } |
208 | |
209 | if (fEsd->GetTrack(n)->GetSign() < 0) { |
210 | fHistoneg->Fill(fEsd->GetTrack(n)->Eta(), |
211 | getphi(fEsd->GetTrack(n)->Phi()), |
212 | fabs(fEsd->GetTrack(n)->Pt())); |
213 | } |
214 | } |
215 | |
216 | // fHistoposclone = (TH2F*)fHistopos->Clone("histoposclone"); |
217 | // fHistonegclone = (TH2F*)fHistoneg->Clone("histonegclone"); |
218 | // fHistoposclone->SetName("histoposclone"); |
219 | // fHistonegclone->SetName("histonegclone"); |
220 | |
221 | fData->DataChanged(); |
222 | |
223 | FilterData(); |
224 | |
225 | return fData; |
226 | } |
227 | |
228 | //______________________________________________________________________________ |
229 | TEveCaloDataHist* AliEveLego::LoadAllData() |
230 | { |
231 | |
232 | fHistopos_all_events->Reset(); |
233 | fHistoneg_all_events->Reset(); |
234 | |
235 | TTree* t = AliEveEventManager::GetMaster()->GetESDTree(); |
236 | |
237 | // Getting current tracks for each event, filling histograms |
238 | for (int event = 0; event < t->GetEntries(); event++) { |
239 | t->GetEntry(event); |
240 | for (int n = 0; n < fEsd->GetNumberOfTracks(); ++n) { |
241 | |
242 | if (fEsd->GetTrack(n)->GetSign() > 0) { |
243 | fHistopos_all_events->Fill(fEsd->GetTrack(n)->Eta(), |
244 | getphi(fEsd->GetTrack(n)->Phi()), |
245 | fabs(fEsd->GetTrack(n)->Pt())); |
246 | } else { |
247 | fHistoneg_all_events->Fill(fEsd->GetTrack(n)->Eta(), |
248 | getphi(fEsd->GetTrack(n)->Phi()), |
249 | fabs(fEsd->GetTrack(n)->Pt())); |
250 | } |
251 | } |
252 | } |
253 | |
254 | fData_all_events->DataChanged(); |
255 | |
256 | return fData_all_events; |
257 | } |
258 | |
259 | //______________________________________________________________________________ |
260 | TEveCaloDataHist* AliEveLego::FilterData() |
261 | { |
262 | // Tracks selection |
263 | if ( fTracksId == 2 ) |
264 | { |
265 | fHistopos->Reset(); |
266 | fHistoneg->Reset(); |
267 | |
268 | const AliESDVertex *pv = fEsd->GetPrimaryVertex(); |
269 | |
270 | for (Int_t n = 0; n < pv->GetNIndices(); n++ ) |
271 | { |
272 | AliESDtrack *at = fEsd->GetTrack(pv->GetIndices()[n]); |
273 | if (at->GetSign() > 0) { |
274 | fHistopos->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt())); |
275 | } |
276 | if (at->GetSign() < 0) { |
277 | fHistoneg->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt())); |
278 | } |
279 | } |
280 | } |
281 | |
282 | fData->DataChanged(); |
283 | |
284 | // Max Pt threshold |
285 | if (GetPtMax() >= fMaxPt){ |
286 | for (Int_t binx = 1; binx <= 100; binx++) { |
287 | for (Int_t biny = 1; biny <= 80; biny++) { |
288 | if (fHistopos->GetBinContent(binx, biny) >= fMaxPt) |
289 | { |
290 | fHistopos->SetBinContent(binx, biny, fMaxPt); |
291 | } |
292 | if (fHistoneg->GetBinContent(binx, biny) >= fMaxPt) |
293 | { |
294 | fHistoneg->SetBinContent(binx, biny, fMaxPt); |
295 | } |
296 | } |
297 | } |
298 | } |
299 | |
300 | // Positive only |
301 | if ( fChargeId == 2 ) fHistoneg->Reset(); |
302 | |
303 | // Negative only |
304 | if ( fChargeId == 3 ) fHistopos->Reset(); |
305 | |
306 | fData->DataChanged(); |
307 | |
308 | return fData; |
309 | } |
310 | |
311 | |
312 | //______________________________________________________________________________ |
313 | TEveCaloDataHist* AliEveLego::FilterAllData() |
314 | { |
315 | // Tracks selection |
316 | if ( fTracksIdAE == 2 ) |
317 | { |
318 | fHistopos_all_events->Reset(); |
319 | fHistoneg_all_events->Reset(); |
320 | |
321 | TTree* t = AliEveEventManager::GetMaster()->GetESDTree(); |
322 | |
323 | // Getting current tracks for each event, filling histograms |
324 | for (int event = 0; event < t->GetEntries(); event++) { |
a096a7bf |
325 | |
3d94b490 |
326 | t->GetEntry(event); |
327 | |
328 | const AliESDVertex *pv = fEsd->GetPrimaryVertex(); |
329 | |
330 | for (Int_t n = 0; n < pv->GetNIndices(); n++ ) |
331 | { |
332 | AliESDtrack *at = fEsd->GetTrack(pv->GetIndices()[n]); |
a096a7bf |
333 | |
3d94b490 |
334 | if (at->GetSign() > 0) { |
335 | fHistopos_all_events->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt())); |
336 | } |
337 | if (at->GetSign() < 0) { |
338 | fHistoneg_all_events->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt())); |
339 | } |
340 | } |
a096a7bf |
341 | } |
3d94b490 |
342 | } else { |
343 | LoadAllData(); |
344 | } |
345 | |
346 | fData_all_events->DataChanged(); |
a096a7bf |
347 | |
3d94b490 |
348 | // Max Pt threshold |
349 | if (GetPtMaxAE() >= fMaxPtAE){ |
350 | for (Int_t binx = 1; binx <= 100; binx++) { |
351 | for (Int_t biny = 1; biny <= 80; biny++) { |
352 | if (fHistopos_all_events->GetBinContent(binx, biny) >= fMaxPtAE) |
353 | { |
354 | fHistopos_all_events->SetBinContent(binx, biny, fMaxPtAE); |
355 | } |
356 | if (fHistoneg_all_events->GetBinContent(binx, biny) >= fMaxPtAE) |
357 | { |
358 | fHistoneg_all_events->SetBinContent(binx, biny, fMaxPtAE); |
359 | } |
360 | } |
361 | } |
362 | } |
363 | |
364 | // Positive only |
365 | if ( fChargeIdAE == 2 ) fHistoneg_all_events->Reset(); |
366 | |
367 | // Negative only |
368 | if ( fChargeIdAE == 3 ) fHistopos_all_events->Reset(); |
369 | |
370 | fData_all_events->DataChanged(); |
371 | |
372 | gEve->Redraw3D(kTRUE); |
373 | |
374 | return fData_all_events; |
375 | } |
376 | |
377 | //______________________________________________________________________________ |
378 | void AliEveLego::Update() |
379 | { |
380 | // Load/Reload data |
381 | LoadData(); |
382 | |
383 | // Create new histo2d |
384 | CreateHistoLego(); |
385 | |
386 | // Create 3d view |
387 | Create3DView(); |
388 | |
389 | // Show information about event; |
390 | ShowEventSeletion(fShowEventsInfo, kTRUE); |
391 | |
392 | // Update the viewers |
393 | gEve->Redraw3D(kTRUE); |
394 | } |
395 | |
396 | //______________________________________________________________________________ |
397 | TEveCaloLego* AliEveLego::CreateHistoLego() |
398 | { |
399 | // Viewer initialization, tab creation |
400 | if (fHisto2d_v == 0) { |
401 | TEveWindowSlot *fslot = 0; |
402 | TEveBrowser *fbrowser = gEve->GetBrowser(); |
403 | |
404 | fslot = TEveWindow::CreateWindowInTab(fbrowser->GetTabRight()); |
405 | fslot->MakeCurrent(); |
406 | fHisto2d_v = gEve->SpawnNewViewer("2D Lego Histogram", "2D Lego Histogram"); |
407 | fHisto2d_s = gEve->SpawnNewScene("2D Lego Histogram", "2D Lego Histogram"); |
408 | fHisto2d_v->AddScene(fHisto2d_s); |
409 | fHisto2d_v->SetElementName("2D Lego Viewer"); |
410 | fHisto2d_s->SetElementName("2D Lego Scene"); |
411 | |
412 | fGlv = fHisto2d_v->GetGLViewer(); |
413 | fHisto2d_lego_overlay = new TEveCaloLegoOverlay(); |
414 | fGlv->AddOverlayElement(fHisto2d_lego_overlay); |
415 | fGlv->SetCurrentCamera(TGLViewer::kCameraPerspXOY); |
416 | |
417 | fHisto2d_s->AddElement(fLego); |
418 | |
419 | // move to real world coordinates |
420 | fLego->InitMainTrans(); |
421 | Float_t sc = TMath::Min(fLego->GetEtaRng(), fLego->GetPhiRng()); |
422 | fLego->RefMainTrans().SetScale(sc, sc, sc); |
423 | |
424 | // set event handler to move from perspective to orthographic view. |
425 | fGlv->SetEventHandler(new TEveLegoEventHandler(fGlv->GetGLWidget(), fGlv, fLego)); |
426 | |
427 | fHisto2d_lego_overlay->SetCaloLego(fLego); |
428 | } |
429 | |
430 | return fLego; |
431 | } |
432 | |
433 | //______________________________________________________________________________ |
434 | TEveCaloLego* AliEveLego::CreateHistoLego(TEveWindowSlot *slot) |
435 | { |
436 | // Viewer initialization, tab creation |
437 | if (fHisto2d_all_events_v0 == 0) { |
438 | |
439 | slot->MakeCurrent(); |
440 | fHisto2d_all_events_v0 = gEve->SpawnNewViewer("2D Lego Histogram", "2D Lego Histogram"); |
441 | fHisto2d_all_events_s0 = gEve->SpawnNewScene("2D Lego Histogram", "2D Lego Histogram"); |
442 | fHisto2d_all_events_v0->AddScene(fHisto2d_all_events_s0); |
443 | fHisto2d_all_events_v0->SetElementName("2D Lego Viewer"); |
444 | fHisto2d_all_events_s0->SetElementName("2D Lego Scene"); |
445 | |
446 | TGLViewer* glv = fHisto2d_all_events_v0->GetGLViewer(); |
447 | fHisto2d_all_events_lego_overlay = new TEveCaloLegoOverlay(); |
448 | glv->AddOverlayElement(fHisto2d_all_events_lego_overlay); |
449 | glv->SetCurrentCamera(TGLViewer::kCameraPerspXOY); |
450 | |
451 | // Plotting histogram lego |
452 | fLego_all_events = new TEveCaloLego(fData_all_events); |
453 | fHisto2d_all_events_s0->AddElement(fLego_all_events); |
454 | |
455 | // Move to real world coordinates |
456 | fLego_all_events->InitMainTrans(); |
457 | Float_t sc = TMath::Min(fLego_all_events->GetEtaRng(), fLego_all_events->GetPhiRng()); |
458 | fLego_all_events->RefMainTrans().SetScale(sc, sc, sc); |
459 | |
460 | // Set event handler to move from perspective to orthographic view. |
461 | glv->SetEventHandler(new TEveLegoEventHandler(glv->GetGLWidget(), glv, fLego_all_events)); |
462 | |
463 | fHisto2d_all_events_lego_overlay->SetCaloLego(fLego_all_events); |
464 | } |
465 | |
466 | return fLego_all_events; |
467 | } |
468 | |
469 | //______________________________________________________________________________ |
470 | TEveCalo3D* AliEveLego::Create3DView() |
471 | { |
472 | //initialization |
473 | if (fHisto2d_s2 == 0) { |
474 | fHisto2d_s2 = gEve->SpawnNewScene("3D Histogram", "3D Histogram"); |
475 | gEve->GetDefaultViewer()->AddScene(fHisto2d_s2); |
476 | fHisto2d_s2->SetElementName("3D Histogram Scene"); |
477 | fHisto2d_s2->AddElement(fCalo3d); |
478 | } |
479 | |
480 | return fCalo3d; |
481 | } |
482 | |
483 | //______________________________________________________________________________ |
484 | TEveCalo3D* AliEveLego::Create3DView(TEveWindowSlot *slot) |
485 | { |
486 | if ( fHisto2d_all_events_v1 == 0 ) { |
487 | |
488 | slot->MakeCurrent(); |
489 | fHisto2d_all_events_v1 = gEve->SpawnNewViewer("3D Histogram", "3D Histogram"); |
490 | fHisto2d_all_events_s1 = gEve->SpawnNewScene("3D Histogram", "3D Histogram"); |
491 | fHisto2d_all_events_v1->AddScene(fHisto2d_all_events_s1); |
492 | fHisto2d_all_events_v1->SetElementName("3D Histogram Viewer"); |
493 | fHisto2d_all_events_s1->SetElementName("3D Histogram Scene"); |
494 | |
495 | fCalo3d_all_events = new TEveCalo3D(fData_all_events); |
496 | |
497 | fCalo3d_all_events->SetBarrelRadius(550); |
498 | fCalo3d_all_events->SetEndCapPos(550); |
499 | fHisto2d_all_events_s1->AddElement(fCalo3d_all_events); |
500 | } |
501 | |
502 | return fCalo3d_all_events; |
503 | } |
504 | |
505 | //______________________________________________________________________________ |
506 | void AliEveLego::CreateProjections(TEveWindowSlot* slot1, TEveWindowSlot* slot2){ |
507 | |
508 | if (fHisto2d_all_events_v2 == 0) { |
509 | |
510 | slot1->MakeCurrent(); |
511 | fHisto2d_all_events_v2 = gEve->SpawnNewViewer("RPhi projection", "RPhi projection"); |
512 | fHisto2d_all_events_s2 = gEve->SpawnNewScene("RPhi projection", "RPhi projection"); |
513 | fHisto2d_all_events_v2->AddScene(fHisto2d_all_events_s2); |
514 | fHisto2d_all_events_v2->SetElementName("RPhi Projection Viewer"); |
515 | fHisto2d_all_events_s2->SetElementName("RPhi Projection Scene"); |
516 | |
517 | TEveProjectionManager* mng1 = new TEveProjectionManager(); |
518 | mng1->SetProjection(TEveProjection::kPT_RPhi); |
519 | |
520 | TEveProjectionAxes* axeg_histo2d_all_events_s1 = new TEveProjectionAxes(mng1); |
521 | fHisto2d_all_events_s2->AddElement(axeg_histo2d_all_events_s1); |
522 | TEveCalo2D* fcalo2d1 = (TEveCalo2D*) mng1->ImportElements(fCalo3d_all_events); |
523 | fHisto2d_all_events_s2->AddElement(fcalo2d1); |
524 | |
525 | fHisto2d_all_events_v2->GetGLViewer()->SetCurrentCamera(TGLViewer::kCameraOrthoXOY); |
526 | } |
527 | |
528 | if (fHisto2d_all_events_v3 == 0) { |
529 | |
530 | slot2->MakeCurrent(); |
531 | fHisto2d_all_events_v3 = gEve->SpawnNewViewer("RhoZ projection", "RhoZ projection"); |
532 | fHisto2d_all_events_s3 = gEve->SpawnNewScene("RhoZ projection", "RhoZ projection"); |
533 | fHisto2d_all_events_v3->AddScene(fHisto2d_all_events_s3); |
534 | fHisto2d_all_events_v3->SetElementName("RhoZ Projection Viewer"); |
535 | fHisto2d_all_events_s3->SetElementName("RhoZ Projection Viewer"); |
536 | |
537 | TEveProjectionManager* mng2 = new TEveProjectionManager(); |
538 | mng2->SetProjection(TEveProjection::kPT_RhoZ); |
539 | |
540 | TEveProjectionAxes* axeg_histo2d_all_events_s2 = new TEveProjectionAxes(mng2); |
541 | fHisto2d_all_events_s3->AddElement(axeg_histo2d_all_events_s2); |
542 | TEveCalo2D* fcalo2d2 = (TEveCalo2D*) mng2->ImportElements(fCalo3d_all_events); |
543 | fHisto2d_all_events_s3->AddElement(fcalo2d2); |
544 | |
545 | fHisto2d_all_events_v3->GetGLViewer()->SetCurrentCamera(TGLViewer::kCameraOrthoXOY); |
546 | } |
547 | |
548 | return; |
549 | } |
550 | |
551 | //______________________________________________________________________________ |
552 | TEveCaloDataHist* AliEveLego::LoadAllEvents() |
553 | { |
554 | if ( fHisto2d_all_events_slot == 0 ) { |
555 | |
556 | printf("Filling histogram...\n"); |
a096a7bf |
557 | TStopwatch timer; |
558 | timer.Start(); |
3d94b490 |
559 | |
560 | // Creating 2D histograms |
561 | fHistopos_all_events = new TH2F("fHistopos_all_events","Histo 2d positive", |
562 | 100,-1.5,1.5,80,-pi,pi); |
563 | fHistoneg_all_events = new TH2F("fHistoneg_all_events","Histo 2d negative", |
564 | 100,-1.5,1.5,80,-pi,pi); |
565 | |
a096a7bf |
566 | fHistopos_all_events->SetDirectory(0); |
567 | fHistoneg_all_events->SetDirectory(0); |
568 | |
3d94b490 |
569 | fData_all_events = new TEveCaloDataHist(); |
570 | fData_all_events->AddHistogram(fHistoneg_all_events); |
571 | fData_all_events->RefSliceInfo(0).Setup("NegCg:", 0, kBlue); |
572 | fData_all_events->AddHistogram(fHistopos_all_events); |
573 | fData_all_events->RefSliceInfo(1).Setup("PosCg:", 0, kRed); |
574 | fData_all_events->GetEtaBins()->SetTitleFont(120); |
575 | fData_all_events->GetEtaBins()->SetTitle("h"); |
576 | fData_all_events->GetPhiBins()->SetTitleFont(120); |
577 | fData_all_events->GetPhiBins()->SetTitle("f"); |
578 | fData_all_events->IncDenyDestroy(); |
579 | |
580 | // Creating frames |
581 | fHisto2d_all_events_slot = TEveWindow::CreateWindowInTab(gEve->GetBrowser()->GetTabRight()); |
582 | TEveWindowPack* packH = fHisto2d_all_events_slot->MakePack(); |
583 | packH->SetElementName("Projections"); |
584 | packH->SetHorizontal(); |
585 | packH->SetShowTitleBar(kFALSE); |
586 | |
587 | fHisto2d_all_events_slot = packH->NewSlot(); |
588 | TEveWindowPack* pack0 = fHisto2d_all_events_slot->MakePack(); |
589 | pack0->SetShowTitleBar(kFALSE); |
590 | TEveWindowSlot* slotLeftTop = pack0->NewSlot(); |
591 | TEveWindowSlot* slotLeftBottom = pack0->NewSlot(); |
592 | |
593 | fHisto2d_all_events_slot = packH->NewSlot(); |
594 | TEveWindowPack* pack1 = fHisto2d_all_events_slot->MakePack(); |
595 | pack1->SetShowTitleBar(kFALSE); |
596 | TEveWindowSlot* slotRightTop = pack1->NewSlot(); |
597 | TEveWindowSlot* slotRightBottom = pack1->NewSlot(); |
598 | |
599 | // Creating viewers and scenes |
600 | Create3DView(slotLeftTop); |
601 | CreateHistoLego(slotLeftBottom); |
602 | CreateProjections(slotRightTop, slotRightBottom); |
603 | |
604 | LoadAllData(); |
605 | |
606 | gEve->Redraw3D(kTRUE); |
607 | |
608 | printf("Filling histogram... Finished\n"); |
a096a7bf |
609 | timer.Stop(); |
610 | timer.Print(); |
611 | |
3d94b490 |
612 | } |
613 | |
614 | return fData_all_events; |
615 | } |
616 | |
617 | //______________________________________________________________________________ |
618 | Float_t AliEveLego::GetPtMax() |
619 | { |
620 | return fData->GetMaxVal(fLego->GetPlotEt()); |
621 | } |
622 | |
623 | //______________________________________________________________________________ |
624 | Float_t AliEveLego::GetPtMaxAE() |
625 | { |
626 | return fData_all_events->GetMaxVal(fLego_all_events->GetPlotEt()); |
627 | } |
628 | |
629 | //______________________________________________________________________________ |
630 | void AliEveLego::SetMaxPt(Double_t val) |
631 | { |
632 | // Add new maximum |
633 | fMaxPt = val; |
634 | Update(); |
635 | } |
636 | |
637 | //______________________________________________________________________________ |
638 | void AliEveLego::SetMaxPtAE(Double_t val) |
639 | { |
640 | // Add new maximum |
641 | fMaxPtAE = val; |
642 | FilterAllData(); |
643 | } |
644 | |
645 | //______________________________________________________________________________ |
646 | void AliEveLego::SetThreshold(Double_t val) |
647 | { |
648 | // Setting up the new threshold for all histograms |
649 | fData->SetSliceThreshold(0,val); |
650 | fData->SetSliceThreshold(1,val); |
651 | fData->DataChanged(); |
652 | |
653 | gEve->Redraw3D(kTRUE); |
654 | } |
655 | |
656 | //______________________________________________________________________________ |
657 | void AliEveLego::SetThresholdAE(Double_t val) |
658 | { |
659 | // Setting up the new threshold for all histograms |
660 | fData_all_events->SetSliceThreshold(0,val); |
661 | fData_all_events->SetSliceThreshold(1,val); |
662 | fData_all_events->DataChanged(); |
663 | |
664 | gEve->Redraw3D(kTRUE); |
665 | } |
666 | |
667 | //______________________________________________________________________________ |
668 | void AliEveLego::SetEventSelection() |
669 | { |
670 | if (fShowEventsInfo == 0) |
671 | { |
672 | fShowEventsInfo = 1; |
673 | } else { |
674 | fShowEventsInfo = 0; |
675 | } |
676 | |
677 | ShowEventSeletion(fShowEventsInfo); |
678 | } |
679 | |
680 | //______________________________________________________________________________ |
681 | void AliEveLego::ShowEventSeletion(Bool_t show, Bool_t updateonly) |
682 | { |
683 | if (show == 0) |
684 | { |
685 | gEve->GetDefaultGLViewer()->RemoveOverlayElement(fGButton); |
686 | fAl->Get3DView()->GetGLViewer()->RemoveOverlayElement(fGButton); |
687 | fHisto2d_v->GetGLViewer()->RemoveOverlayElement(fGButton); |
688 | |
689 | gEve->GetDefaultGLViewer()->RemoveOverlayElement(fB1); |
690 | fAl->Get3DView()->GetGLViewer()->RemoveOverlayElement(fB1); |
691 | fHisto2d_v->GetGLViewer()->RemoveOverlayElement(fB1); |
692 | |
693 | gEve->GetDefaultGLViewer()->RemoveOverlayElement(fB2); |
694 | fAl->Get3DView()->GetGLViewer()->RemoveOverlayElement(fB2); |
695 | fHisto2d_v->GetGLViewer()->RemoveOverlayElement(fB2); |
696 | |
697 | } else { |
698 | |
699 | //Collision candidate |
700 | if (updateonly == kFALSE) { |
701 | gEve->GetDefaultGLViewer()->AddOverlayElement(fGButton); |
702 | fAl->Get3DView()->GetGLViewer()->AddOverlayElement(fGButton); |
703 | fHisto2d_v->GetGLViewer()->AddOverlayElement(fGButton); |
704 | } |
705 | |
706 | Bool_t ev = fPhysicsSelection->IsCollisionCandidate(fEsd); |
707 | |
708 | if (ev == 1) |
709 | { |
710 | fGButton->SetText("Collision candidate: YES"); |
711 | } else { |
712 | fGButton->SetText("Collision candidate: NO "); |
713 | } |
714 | |
715 | // Beam 1 & 2 setup: method 1 |
716 | if (updateonly == kFALSE) { |
717 | gEve->GetDefaultGLViewer()->AddOverlayElement(fB1); |
718 | fAl->Get3DView()->GetGLViewer()->AddOverlayElement(fB1); |
719 | fHisto2d_v->GetGLViewer()->AddOverlayElement(fB1); |
720 | |
721 | gEve->GetDefaultGLViewer()->AddOverlayElement(fB2); |
722 | fAl->Get3DView()->GetGLViewer()->AddOverlayElement(fB2); |
723 | fHisto2d_v->GetGLViewer()->AddOverlayElement(fB2); |
724 | } |
725 | |
726 | Bool_t b1 = fEsd->IsTriggerClassFired("CINT1A-ABCE-NOPF-ALL"); |
727 | Bool_t b2 = fEsd->IsTriggerClassFired("CINT1C-ABCE-NOPF-ALL"); |
728 | Bool_t b12 = fEsd->IsTriggerClassFired("CINT1B-ABCE-NOPF-ALL"); |
729 | |
730 | if (b1 == 1 || b12 == 1) |
731 | { |
732 | fB1->SetText("Beam 1: YES"); |
733 | fB1->SetBackColor(0x00ff00); |
734 | } else { |
735 | fB1->SetText("Beam 1: NO"); |
736 | fB1->SetBackColor(0xff0000); |
737 | } |
738 | |
739 | if (b2 == 1 || b12 == 1) |
740 | { |
741 | fB2->SetText("Beam 2: YES"); |
742 | fB2->SetBackColor(0x00ff00); |
743 | } else { |
744 | fB2->SetText("Beam 2: NO"); |
745 | fB2->SetBackColor(0xff0000); |
746 | } |
747 | } |
748 | |
749 | gEve->Redraw3D(kTRUE); |
750 | |
751 | } |
752 | |
753 | //______________________________________________________________________________ |
754 | void AliEveLego::SelectEventSelection(Int_t id) |
755 | { |
756 | if (id == 0) |
757 | { |
758 | fEventSelector->SetSelectOnTriggerType(kFALSE); |
759 | } else { |
760 | if (id == 1) fEventSelector->SetTriggerType("CINT1A-ABCE-NOPF-ALL"); |
761 | if (id == 2) fEventSelector->SetTriggerType("CINT1C-ABCE-NOPF-ALL"); |
762 | if (id == 3) fEventSelector->SetTriggerType("CINT1B-ABCE-NOPF-ALL"); |
763 | fEventSelector->SetSelectOnTriggerType(kTRUE); |
764 | } |
765 | } |
766 | |
767 | //______________________________________________________________________________ |
768 | void AliEveLego::ShowPrevEvent() |
769 | { |
770 | AliEveEventManager::GetMaster()->PrevEvent(); |
771 | } |
772 | |
773 | //______________________________________________________________________________ |
774 | void AliEveLego::ShowNextEvent() |
775 | { |
776 | AliEveEventManager::GetMaster()->NextEvent(); |
777 | } |
778 | /******************************************************************************/ |
779 | |
780 | |