]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveBase/AliEveLego.cxx
add a better alpha cut to default cut
[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 #include "TStopwatch.h"
35
36 //______________________________________________________________________________
37 // Allow histograms visualization in 2D and 3D.
38 //
39
40 ClassImp(AliEveLego)
41 Double_t kPi = 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),
52   fMaxPtAE(10000),
53   fEsd(0),
54   fPhysicsSelection(0),
55   fHistopos(0),
56   fHistoposAllEvents(0),
57   fHistoneg(0),
58   fHistonegAllEvents(0),
59   fData(0),
60   fDataAllEvents(0),
61   fLego(0),
62   fLegoAllEvents(0),
63   fCalo3d(0),
64   fCalo3dAllEvents(0),
65   fGlv(0),
66   fHisto2dv(0),
67   fHisto2ds(0),
68   fHisto2ds2(0),
69   fHisto2dAllEventsv0(0),
70   fHisto2dAllEventsv1(0),
71   fHisto2dAllEventsv2(0),
72   fHisto2dAllEventsv3(0),
73   fHisto2dAllEventss0(0),
74   fHisto2dAllEventss1(0),
75   fHisto2dAllEventss2(0),
76   fHisto2dAllEventss3(0),
77   fAl(0),
78   fHisto2dLegoOverlay(0),
79   fHisto2dAllEventsLegoOverlay(0),
80   fHisto2dAllEventsSlot(0),
81   fEventSelector(0),
82   fShowEventsInfo(0),
83   fGButton(0),
84   fB1(0),
85   fB2(0)
86 {
87   // Constructor.
88   gEve->AddToListTree(this,0);
89
90   // Get Current ESD event
91   fEsd = AliEveEventManager::AssertESD();
92
93   fPhysicsSelection = new AliPhysicsSelection();
94   fPhysicsSelection->Initialize(fEsd->GetRunNumber());
95
96   fEventSelector = AliEveEventManager::GetMaster()->GetEventSelector();
97
98   fHistopos = new TH2F("histopos","Histo 2d positive", 100, -1.5, 1.5, 80, -kPi, kPi);
99   fHistoneg = new TH2F("histoneg","Histo 2d negative", 100, -1.5, 1.5, 80, -kPi, kPi);
100
101   fHistopos->SetDirectory(0);
102   fHistoneg->SetDirectory(0);
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    // deleting variables
143    delete fEsd;
144    delete fPhysicsSelection;
145    delete fHistopos;
146    delete fHistoposAllEvents;
147    delete fHistoneg;
148    delete fHistonegAllEvents;
149
150    delete fData;
151    delete fDataAllEvents;
152    delete fLego;
153    delete fLegoAllEvents;
154    delete fCalo3d;
155    delete fCalo3dAllEvents;
156    delete fGlv;
157
158    delete fHisto2dv;
159    delete fHisto2ds;
160    delete fHisto2ds2;
161    delete fHisto2dAllEventsv0;
162    delete fHisto2dAllEventsv1;
163    delete fHisto2dAllEventsv2;
164    delete fHisto2dAllEventsv3;
165    delete fHisto2dAllEventss0;
166    delete fHisto2dAllEventss1;
167    delete fHisto2dAllEventss2;
168    delete fHisto2dAllEventss3;
169
170    delete fAl;
171    delete fHisto2dLegoOverlay;
172    delete fHisto2dAllEventsLegoOverlay;
173    delete fHisto2dAllEventsSlot;
174
175    delete fEventSelector;
176    delete fGButton;
177    delete fB1;
178    delete fB2;
179 }
180
181 namespace
182 {
183   //____________________________________________________________________________
184   Double_t getphi(Double_t phi)
185   {
186     // phi correction for alice
187
188     if (phi > TMath::Pi()) {
189       phi -= TMath::TwoPi();
190     }
191     return phi;
192   }
193 }
194
195 //______________________________________________________________________________
196 TEveCaloDataHist* AliEveLego::LoadData()
197 {
198    // Load data from ESD tree
199    fHistopos->Reset();
200    fHistoneg->Reset();
201
202    // Getting current tracks, filling histograms
203    for (int n = 0; n < fEsd->GetNumberOfTracks(); ++n) {
204
205       if (fEsd->GetTrack(n)->GetSign() > 0) {
206          fHistopos->Fill(fEsd->GetTrack(n)->Eta(),
207                          getphi(fEsd->GetTrack(n)->Phi()),
208                          fabs(fEsd->GetTrack(n)->Pt()));
209       }
210
211       if (fEsd->GetTrack(n)->GetSign() < 0) {
212          fHistoneg->Fill(fEsd->GetTrack(n)->Eta(),
213                          getphi(fEsd->GetTrack(n)->Phi()),
214                          fabs(fEsd->GetTrack(n)->Pt()));
215       }
216    }
217
218    fData->DataChanged();
219
220    FilterData();
221
222    return fData;
223 }
224
225 //______________________________________________________________________________
226 TEveCaloDataHist* AliEveLego::LoadAllData()
227 {
228    // load data from all events ESD
229    fHistoposAllEvents->Reset();
230    fHistonegAllEvents->Reset();
231
232    TTree* t = AliEveEventManager::GetMaster()->GetESDTree();
233
234    // Getting current tracks for each event, filling histograms
235    for (int event = 0; event < t->GetEntries(); event++) {
236       t->GetEntry(event);
237          for (int n = 0; n < fEsd->GetNumberOfTracks(); ++n) {
238
239             if (fEsd->GetTrack(n)->GetSign() > 0) {
240                fHistoposAllEvents->Fill(fEsd->GetTrack(n)->Eta(),
241                                           getphi(fEsd->GetTrack(n)->Phi()),
242                                           fabs(fEsd->GetTrack(n)->Pt()));
243             } else {
244                fHistonegAllEvents->Fill(fEsd->GetTrack(n)->Eta(),
245                                           getphi(fEsd->GetTrack(n)->Phi()),
246                                           fabs(fEsd->GetTrack(n)->Pt()));
247             }
248          }
249    }
250
251    fDataAllEvents->DataChanged();
252
253    return fDataAllEvents;
254 }
255
256 //______________________________________________________________________________
257 TEveCaloDataHist* AliEveLego::FilterData()
258 {
259    // Tracks selection
260    if ( fTracksId == 2 )
261    {
262       fHistopos->Reset();
263       fHistoneg->Reset();
264
265       const AliESDVertex *pv  = fEsd->GetPrimaryVertex();
266
267       for (Int_t n = 0; n < pv->GetNIndices(); n++ )
268       {
269          AliESDtrack *at = fEsd->GetTrack(pv->GetIndices()[n]);
270          if (at->GetSign() > 0) {
271             fHistopos->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
272          }
273          if (at->GetSign() < 0) {
274             fHistoneg->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
275          }
276       }
277    }
278
279    fData->DataChanged();
280
281    // Max Pt threshold
282    if (GetPtMax() >= fMaxPt){
283       for (Int_t binx = 1; binx <= 100; binx++) {
284          for (Int_t biny = 1; biny <= 80; biny++) {
285             if (fHistopos->GetBinContent(binx, biny) >= fMaxPt)
286             {
287                fHistopos->SetBinContent(binx, biny, fMaxPt);
288             }
289             if (fHistoneg->GetBinContent(binx, biny) >= fMaxPt)
290             {
291                fHistoneg->SetBinContent(binx, biny, fMaxPt);
292             }           
293          }
294       }
295    }
296
297    // Positive only
298    if ( fChargeId == 2 ) fHistoneg->Reset();
299
300    // Negative only
301    if ( fChargeId == 3 ) fHistopos->Reset();
302
303    fData->DataChanged();
304
305    return fData;
306 }
307
308
309 //______________________________________________________________________________
310 TEveCaloDataHist* AliEveLego::FilterAllData()
311 {
312    // Tracks selection
313    if ( fTracksIdAE == 2 )
314    {
315       fHistoposAllEvents->Reset();
316       fHistonegAllEvents->Reset();
317
318       TTree* t = AliEveEventManager::GetMaster()->GetESDTree();
319
320       // Getting current tracks for each event, filling histograms
321       for (int event = 0; event < t->GetEntries(); event++) {
322
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
331          if (at->GetSign() > 0) {
332             fHistoposAllEvents->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
333          }
334          if (at->GetSign() < 0) {
335             fHistonegAllEvents->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
336          }
337       }
338       }
339    } else {
340       LoadAllData();
341    }
342    
343    fDataAllEvents->DataChanged();
344
345    // Max Pt threshold
346    if (GetPtMaxAE() >= fMaxPtAE){
347       for (Int_t binx = 1; binx <= 100; binx++) {
348          for (Int_t biny = 1; biny <= 80; biny++) {
349             if (fHistoposAllEvents->GetBinContent(binx, biny) >= fMaxPtAE)
350             {
351                fHistoposAllEvents->SetBinContent(binx, biny, fMaxPtAE);
352             }
353             if (fHistonegAllEvents->GetBinContent(binx, biny) >= fMaxPtAE)
354             {
355                fHistonegAllEvents->SetBinContent(binx, biny, fMaxPtAE);
356             }           
357          }
358       }
359    }
360
361    // Positive only
362    if ( fChargeIdAE == 2 ) fHistonegAllEvents->Reset();
363
364    // Negative only
365    if ( fChargeIdAE == 3 ) fHistoposAllEvents->Reset();
366
367    fDataAllEvents->DataChanged();
368
369    gEve->Redraw3D(kTRUE);
370
371    return fDataAllEvents;
372 }
373
374 //______________________________________________________________________________
375 void AliEveLego::Update()
376 {
377    // Load/Reload data
378    LoadData();
379
380    // Create new histo2d
381    CreateHistoLego();
382
383    // Create 3d view
384    Create3DView();
385
386    // Show information about event;
387    ShowEventSeletion(fShowEventsInfo, kTRUE);
388
389    // Update the viewers
390    gEve->Redraw3D(kTRUE);
391 }
392
393 //______________________________________________________________________________
394 TEveCaloLego* AliEveLego::CreateHistoLego()
395 {
396    // Viewer initialization, tab creation
397    if (fHisto2dv == 0) {
398       TEveWindowSlot *fslot    = 0;
399       TEveBrowser    *fbrowser = gEve->GetBrowser();
400
401       fslot = TEveWindow::CreateWindowInTab(fbrowser->GetTabRight());
402       fslot->MakeCurrent();
403       fHisto2dv = gEve->SpawnNewViewer("2D Lego Histogram", "2D Lego Histogram");
404       fHisto2ds = gEve->SpawnNewScene("2D Lego Histogram", "2D Lego Histogram");
405       fHisto2dv->AddScene(fHisto2ds);
406       fHisto2dv->SetElementName("2D Lego Viewer");
407       fHisto2ds->SetElementName("2D Lego Scene");
408
409       fGlv = fHisto2dv->GetGLViewer();
410       fHisto2dLegoOverlay = new TEveCaloLegoOverlay();
411       fGlv->AddOverlayElement(fHisto2dLegoOverlay);
412       fGlv->SetCurrentCamera(TGLViewer::kCameraPerspXOY);
413
414       fHisto2ds->AddElement(fLego);
415
416       // move to real world coordinates
417       fLego->InitMainTrans();
418       Float_t sc = TMath::Min(fLego->GetEtaRng(), fLego->GetPhiRng());
419       fLego->RefMainTrans().SetScale(sc, sc, sc);
420
421       // set event handler to move from perspective to orthographic view.
422       fGlv->SetEventHandler(new TEveLegoEventHandler(fGlv->GetGLWidget(), fGlv, fLego));
423
424       fHisto2dLegoOverlay->SetCaloLego(fLego);
425    }
426
427    return fLego;
428 }
429
430 //______________________________________________________________________________
431 TEveCaloLego* AliEveLego::CreateHistoLego(TEveWindowSlot *slot)
432 {
433    // Viewer initialization, tab creation
434    if (fHisto2dAllEventsv0 == 0) {
435
436       slot->MakeCurrent();
437       fHisto2dAllEventsv0 = gEve->SpawnNewViewer("2D Lego Histogram", "2D Lego Histogram");
438       fHisto2dAllEventss0 = gEve->SpawnNewScene("2D Lego Histogram", "2D Lego Histogram");
439       fHisto2dAllEventsv0->AddScene(fHisto2dAllEventss0);
440       fHisto2dAllEventsv0->SetElementName("2D Lego Viewer");
441       fHisto2dAllEventss0->SetElementName("2D Lego Scene");
442
443       TGLViewer* glv = fHisto2dAllEventsv0->GetGLViewer();
444       fHisto2dAllEventsLegoOverlay = new TEveCaloLegoOverlay();
445       glv->AddOverlayElement(fHisto2dAllEventsLegoOverlay);
446       glv->SetCurrentCamera(TGLViewer::kCameraPerspXOY);
447
448       // Plotting histogram lego
449       fLegoAllEvents = new TEveCaloLego(fDataAllEvents);
450       fHisto2dAllEventss0->AddElement(fLegoAllEvents);
451
452       // Move to real world coordinates
453       fLegoAllEvents->InitMainTrans();
454       Float_t sc = TMath::Min(fLegoAllEvents->GetEtaRng(), fLegoAllEvents->GetPhiRng());
455       fLegoAllEvents->RefMainTrans().SetScale(sc, sc, sc);
456
457       // Set event handler to move from perspective to orthographic view.
458       glv->SetEventHandler(new TEveLegoEventHandler(glv->GetGLWidget(), glv, fLegoAllEvents));
459
460       fHisto2dAllEventsLegoOverlay->SetCaloLego(fLegoAllEvents);
461    }
462
463    return fLegoAllEvents;
464 }
465
466 //______________________________________________________________________________
467 TEveCalo3D* AliEveLego::Create3DView()
468 {
469    //initialization
470    if (fHisto2ds2 == 0) {
471       fHisto2ds2 = gEve->SpawnNewScene("3D Histogram", "3D Histogram");
472       gEve->GetDefaultViewer()->AddScene(fHisto2ds2);
473       fHisto2ds2->SetElementName("3D Histogram Scene");
474       fHisto2ds2->AddElement(fCalo3d);
475    }
476
477    return fCalo3d;
478 }
479
480 //______________________________________________________________________________
481 TEveCalo3D* AliEveLego::Create3DView(TEveWindowSlot *slot)
482 {
483    // creates a 3d view for the 3d histogram
484    if ( fHisto2dAllEventsv1 == 0 ) {
485
486       slot->MakeCurrent();
487       fHisto2dAllEventsv1 = gEve->SpawnNewViewer("3D Histogram", "3D Histogram");
488       fHisto2dAllEventss1 = gEve->SpawnNewScene("3D Histogram", "3D Histogram");
489       fHisto2dAllEventsv1->AddScene(fHisto2dAllEventss1);
490       fHisto2dAllEventsv1->SetElementName("3D Histogram Viewer");
491       fHisto2dAllEventss1->SetElementName("3D Histogram Scene");
492
493       fCalo3dAllEvents = new TEveCalo3D(fDataAllEvents);
494
495       fCalo3dAllEvents->SetBarrelRadius(550);
496       fCalo3dAllEvents->SetEndCapPos(550);
497       fHisto2dAllEventss1->AddElement(fCalo3dAllEvents);
498    }
499
500    return fCalo3dAllEvents;
501 }
502
503 //______________________________________________________________________________
504 void AliEveLego::CreateProjections(TEveWindowSlot* slot1, TEveWindowSlot* slot2)
505 {
506    // create projections
507    if (fHisto2dAllEventsv2 == 0) {
508
509       slot1->MakeCurrent();
510       fHisto2dAllEventsv2 = gEve->SpawnNewViewer("RPhi projection", "RPhi projection");
511       fHisto2dAllEventss2 = gEve->SpawnNewScene("RPhi projection", "RPhi projection");
512       fHisto2dAllEventsv2->AddScene(fHisto2dAllEventss2);
513       fHisto2dAllEventsv2->SetElementName("RPhi Projection Viewer");
514       fHisto2dAllEventss2->SetElementName("RPhi Projection Scene");
515
516       TEveProjectionManager* mng1 = new TEveProjectionManager();
517       mng1->SetProjection(TEveProjection::kPT_RPhi);
518
519       TEveProjectionAxes* axeghisto2dAllEventss1 = new TEveProjectionAxes(mng1);
520       fHisto2dAllEventss2->AddElement(axeghisto2dAllEventss1);
521       TEveCalo2D* fcalo2d1 = (TEveCalo2D*) mng1->ImportElements(fCalo3dAllEvents);
522       fHisto2dAllEventss2->AddElement(fcalo2d1);
523
524       fHisto2dAllEventsv2->GetGLViewer()->SetCurrentCamera(TGLViewer::kCameraOrthoXOY);
525    }
526
527    if (fHisto2dAllEventsv3 == 0) {
528
529       slot2->MakeCurrent();
530       fHisto2dAllEventsv3 = gEve->SpawnNewViewer("RhoZ projection", "RhoZ projection");
531       fHisto2dAllEventss3 = gEve->SpawnNewScene("RhoZ projection", "RhoZ projection");
532       fHisto2dAllEventsv3->AddScene(fHisto2dAllEventss3);
533       fHisto2dAllEventsv3->SetElementName("RhoZ Projection Viewer");
534       fHisto2dAllEventss3->SetElementName("RhoZ Projection Viewer");
535
536       TEveProjectionManager* mng2 = new TEveProjectionManager();
537       mng2->SetProjection(TEveProjection::kPT_RhoZ);
538
539       TEveProjectionAxes* axeghisto2dAllEventss2 = new TEveProjectionAxes(mng2);
540       fHisto2dAllEventss3->AddElement(axeghisto2dAllEventss2);
541       TEveCalo2D* fcalo2d2 = (TEveCalo2D*) mng2->ImportElements(fCalo3dAllEvents);
542       fHisto2dAllEventss3->AddElement(fcalo2d2);
543
544       fHisto2dAllEventsv3->GetGLViewer()->SetCurrentCamera(TGLViewer::kCameraOrthoXOY);
545    }
546
547    return;
548 }
549
550 //______________________________________________________________________________
551 TEveCaloDataHist* AliEveLego::LoadAllEvents()
552 {
553    // load all events data from ESD
554    if ( fHisto2dAllEventsSlot == 0 ) {
555
556       printf("Filling histogram...\n");
557       TStopwatch timer;
558       timer.Start();
559
560       // Creating 2D histograms
561       fHistoposAllEvents = new TH2F("fHistoposAllEvents","Histo 2d positive",
562                                  100,-1.5,1.5,80,-kPi,kPi);
563       fHistonegAllEvents = new TH2F("fHistonegAllEvents","Histo 2d negative",
564                                  100,-1.5,1.5,80,-kPi,kPi);
565
566       fHistoposAllEvents->SetDirectory(0);
567       fHistonegAllEvents->SetDirectory(0);
568
569       fDataAllEvents = new TEveCaloDataHist();
570       fDataAllEvents->AddHistogram(fHistonegAllEvents);
571       fDataAllEvents->RefSliceInfo(0).Setup("NegCg:", 0, kBlue);
572       fDataAllEvents->AddHistogram(fHistoposAllEvents);
573       fDataAllEvents->RefSliceInfo(1).Setup("PosCg:", 0, kRed);
574       fDataAllEvents->GetEtaBins()->SetTitleFont(120);
575       fDataAllEvents->GetEtaBins()->SetTitle("h");
576       fDataAllEvents->GetPhiBins()->SetTitleFont(120);
577       fDataAllEvents->GetPhiBins()->SetTitle("f");
578       fDataAllEvents->IncDenyDestroy();
579
580       // Creating frames
581       fHisto2dAllEventsSlot = TEveWindow::CreateWindowInTab(gEve->GetBrowser()->GetTabRight());
582       TEveWindowPack* packH = fHisto2dAllEventsSlot->MakePack();
583       packH->SetElementName("Projections");
584       packH->SetHorizontal();
585       packH->SetShowTitleBar(kFALSE);
586
587       fHisto2dAllEventsSlot = packH->NewSlot();
588       TEveWindowPack* pack0 = fHisto2dAllEventsSlot->MakePack();
589       pack0->SetShowTitleBar(kFALSE);
590       TEveWindowSlot*  slotLeftTop   = pack0->NewSlot();
591       TEveWindowSlot* slotLeftBottom = pack0->NewSlot();
592
593       fHisto2dAllEventsSlot = packH->NewSlot();
594       TEveWindowPack* pack1 = fHisto2dAllEventsSlot->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");
609       timer.Stop();
610       timer.Print();
611
612    }
613
614    return fDataAllEvents;
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 fDataAllEvents->GetMaxVal(fLegoAllEvents->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    fDataAllEvents->SetSliceThreshold(0,val);
661    fDataAllEvents->SetSliceThreshold(1,val);
662    fDataAllEvents->DataChanged();
663
664    gEve->Redraw3D(kTRUE);
665 }
666
667 //______________________________________________________________________________
668 void AliEveLego::SetEventSelection()
669 {
670    // activate/deactivate info box
671    if (fShowEventsInfo == 0)
672    {
673       fShowEventsInfo = 1;
674    } else {
675       fShowEventsInfo = 0;
676    }
677
678    ShowEventSeletion(fShowEventsInfo);
679 }
680
681 //______________________________________________________________________________
682 void AliEveLego::ShowEventSeletion(Bool_t show, Bool_t updateonly)
683 {
684    // activate/deactivate info box
685    if (show == 0)
686    {
687       gEve->GetDefaultGLViewer()->RemoveOverlayElement(fGButton);
688       fAl->Get3DView()->GetGLViewer()->RemoveOverlayElement(fGButton);
689       fHisto2dv->GetGLViewer()->RemoveOverlayElement(fGButton);
690
691       gEve->GetDefaultGLViewer()->RemoveOverlayElement(fB1);
692       fAl->Get3DView()->GetGLViewer()->RemoveOverlayElement(fB1);
693       fHisto2dv->GetGLViewer()->RemoveOverlayElement(fB1);
694
695       gEve->GetDefaultGLViewer()->RemoveOverlayElement(fB2);
696       fAl->Get3DView()->GetGLViewer()->RemoveOverlayElement(fB2);
697       fHisto2dv->GetGLViewer()->RemoveOverlayElement(fB2);
698
699    } else {
700
701       //Collision candidate
702       if (updateonly == kFALSE) {
703          gEve->GetDefaultGLViewer()->AddOverlayElement(fGButton);
704          fAl->Get3DView()->GetGLViewer()->AddOverlayElement(fGButton);
705          fHisto2dv->GetGLViewer()->AddOverlayElement(fGButton);
706       }
707
708       Bool_t ev = fPhysicsSelection->IsCollisionCandidate(fEsd);
709
710       if (ev == 1)
711       {
712          fGButton->SetText("Collision candidate: YES");
713       } else {
714          fGButton->SetText("Collision candidate: NO ");
715       }
716
717       // Beam 1 & 2 setup: method 1
718       if (updateonly == kFALSE) {
719          gEve->GetDefaultGLViewer()->AddOverlayElement(fB1);
720          fAl->Get3DView()->GetGLViewer()->AddOverlayElement(fB1);
721          fHisto2dv->GetGLViewer()->AddOverlayElement(fB1);
722
723          gEve->GetDefaultGLViewer()->AddOverlayElement(fB2);
724          fAl->Get3DView()->GetGLViewer()->AddOverlayElement(fB2);
725          fHisto2dv->GetGLViewer()->AddOverlayElement(fB2);
726       }
727
728       Bool_t b1  = fEsd->IsTriggerClassFired("CINT1A-ABCE-NOPF-ALL");
729       Bool_t b2  = fEsd->IsTriggerClassFired("CINT1C-ABCE-NOPF-ALL");
730       Bool_t b12 = fEsd->IsTriggerClassFired("CINT1B-ABCE-NOPF-ALL");
731
732       if (b1 == 1 || b12 == 1)
733       {
734          fB1->SetText("Beam 1: YES");
735          fB1->SetBackColor(0x00ff00);
736       } else {
737          fB1->SetText("Beam 1: NO");
738          fB1->SetBackColor(0xff0000);
739       }
740
741       if (b2 == 1 || b12 == 1)
742       {
743          fB2->SetText("Beam 2: YES");
744          fB2->SetBackColor(0x00ff00);
745       } else {
746          fB2->SetText("Beam 2: NO");
747          fB2->SetBackColor(0xff0000);
748       }
749    }
750
751    gEve->Redraw3D(kTRUE);
752
753 }
754
755 //______________________________________________________________________________
756 void AliEveLego::SelectEventSelection(Int_t id)
757 {
758    // show trigger information
759    if (id == 0)
760    {
761       fEventSelector->SetSelectOnTriggerType(kFALSE);
762    } else {
763       if (id == 1) fEventSelector->SetTriggerType("CINT1A-ABCE-NOPF-ALL");
764       if (id == 2) fEventSelector->SetTriggerType("CINT1C-ABCE-NOPF-ALL");
765       if (id == 3) fEventSelector->SetTriggerType("CINT1B-ABCE-NOPF-ALL");
766       fEventSelector->SetSelectOnTriggerType(kTRUE);
767    }
768 }
769
770 //______________________________________________________________________________
771 void AliEveLego::ShowPrevEvent()
772 {
773    AliEveEventManager::GetMaster()->PrevEvent();
774 }
775
776 //______________________________________________________________________________
777 void AliEveLego::ShowNextEvent()
778 {
779    AliEveEventManager::GetMaster()->NextEvent();
780 }
781 /******************************************************************************/
782
783