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