]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveBase/AliEveLego.cxx
Satoshi request: fTrig and fTime
[u/mrichter/AliRoot.git] / EVE / EveBase / AliEveLego.cxx
1 // $Id$
2 // Author: Stefano Carrazza 2010, CERN, stefano.carrazza@cern.ch
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
16 #include "TH2F.h"
17 #include "TMath.h"
18 #include "TGLViewer.h"
19 #include "TEveWindow.h"
20 #include "TEveManager.h"
21 #include "TEveBrowser.h"
22 #include "TEveViewer.h"
23 #include "TEveScene.h"
24 #include "TEveCaloLegoOverlay.h"
25 #include "TEveCalo.h"
26 #include "TEveCaloData.h"
27 #include "TEveLegoEventHandler.h"
28 #include "TEveTrans.h"
29 #include "TEveProjectionManager.h"
30 #include "TEveProjectionAxes.h"
31 #include "TGLWidget.h"
32 #include "TStopwatch.h"
33
34 //______________________________________________________________________________
35 // This class provides the following features:
36 //
37 // 1) 2D and 3D calorimeter like histograms of tracks distribution,
38 // providing particle selection by charge (positive and negative),
39 // and by specie: electrons, muons, pions, kaons and protons.
40 //
41 // 2) Histograms are plotted around detectors, allowing track selection: all tracks,
42 // only primary tracks. It allows pT maximum and minimum threashold.
43 //
44 // 3) 2D and 3D calorimeter like histograms of particles distribution obtained
45 // from all events.
46 //
47 // 4) Possibility to use AliPhysicsSelection during the all events histograms creation.
48 // It is also possible to switch between real data and simulation data (MC).
49 //
50
51 ClassImp(AliEveLego)
52 Double_t kPi = TMath::Pi();
53
54 //______________________________________________________________________________
55 AliEveLego::AliEveLego(const char* name) :
56   TEveElementList(name),
57   fIsMC(kFALSE),
58   fCollisionCandidatesOnly(kFALSE),
59   fParticleTypeId(0),
60   fParticleTypeIdAE(0),
61   fTracksId(1),
62   fMaxPt(10000),
63   fTracksIdAE(0),
64   fMaxPtAE(10000),
65   fEsd(0),
66   fPhysicsSelection(0),
67   fHistopos(0),
68   fHistoposAllEvents(0),
69   fHistoneg(0),
70   fHistonegAllEvents(0),
71   fHistoElectrons(0),
72   fHistoElectronsAllEvents(0),
73   fHistoMuons(0),
74   fHistoMuonsAllEvents(0),
75   fHistoPions(0),
76   fHistoPionsAllEvents(0),
77   fHistoKaons(0),
78   fHistoKaonsAllEvents(0),
79   fHistoProtons(0),
80   fHistoProtonsAllEvents(0),
81   fData(0),
82   fDataAllEvents(0),
83   fLego(0),
84   fLegoAllEvents(0),
85   fCalo3d(0),
86   fCalo3dAllEvents(0),
87   fGlv(0),
88   fHisto2dv(0),
89   fHisto2ds(0),
90   fHisto2ds2(0),
91   fHisto2dAllEventsv0(0),
92   fHisto2dAllEventsv1(0),
93   fHisto2dAllEventsv2(0),
94   fHisto2dAllEventsv3(0),
95   fHisto2dAllEventss0(0),
96   fHisto2dAllEventss1(0),
97   fHisto2dAllEventss2(0),
98   fHisto2dAllEventss3(0),
99   fAl(0),
100   fHisto2dLegoOverlay(0),
101   fHisto2dAllEventsLegoOverlay(0),
102   fHisto2dAllEventsSlot(0)
103 {
104   // Constructor.
105   gEve->AddToListTree(this,0);
106
107   // Get Current ESD event
108   fEsd = AliEveEventManager::AssertESD();
109
110   // Particles types per default
111   fParticleTypeId = new Bool_t[7];
112   fParticleTypeIdAE = new Bool_t[7];
113
114   for (Int_t s = 0; s < 7; s++)
115   {
116     if (s > 1)
117     {
118       fParticleTypeId[s] = kFALSE;
119       fParticleTypeIdAE[s] = kFALSE;
120     } else {
121       fParticleTypeId[s] = kTRUE;
122       fParticleTypeIdAE[s] = kTRUE;
123     }
124   }
125
126   // Loading Physics Selection to determine the collision candidates
127   fPhysicsSelection = new AliPhysicsSelection();
128   fPhysicsSelection->Initialize(fEsd);
129
130   fHistopos       = new TH2F("histopos","Histo 2d positive", 100, -1.5, 1.5, 80, -kPi, kPi);
131   fHistoneg       = new TH2F("histoneg","Histo 2d negative", 100, -1.5, 1.5, 80, -kPi, kPi);
132   fHistoElectrons = new TH2F("histoele","Histo 2d electron", 100, -1.5, 1.5, 80, -kPi, kPi);
133   fHistoMuons     = new TH2F("histomuo","Histo 2d muons   ", 100, -1.5, 1.5, 80, -kPi, kPi);
134   fHistoPions     = new TH2F("histopio","Histo 2d pions   ", 100, -1.5, 1.5, 80, -kPi, kPi);
135   fHistoKaons     = new TH2F("histokao","Histo 2d kaons   ", 100, -1.5, 1.5, 80, -kPi, kPi);
136   fHistoProtons   = new TH2F("histopro","Histo 2d protons ", 100, -1.5, 1.5, 80, -kPi, kPi);
137
138   fHistopos->SetDirectory(0);
139   fHistoneg->SetDirectory(0);
140   fHistoElectrons->SetDirectory(0);
141   fHistoMuons->SetDirectory(0);
142   fHistoPions->SetDirectory(0);
143   fHistoKaons->SetDirectory(0);
144   fHistoProtons->SetDirectory(0);
145
146   // Colors from get_pdg_color() in /alice-macros/kine_tracks.C
147   static Color_t fElectro = 5;
148   static Color_t fMuon = 30;
149   static Color_t fPion = 3;
150   static Color_t fKaon = 38;
151   static Color_t fProton = 10;
152
153   // Adding data to TEveCaloDataHist
154   fData = new TEveCaloDataHist();
155   fData->AddHistogram(fHistoneg);
156   fData->RefSliceInfo(0).Setup("NegCg:", 0, kBlue);
157   fData->AddHistogram(fHistopos);
158   fData->RefSliceInfo(1).Setup("PosCg:", 0, kRed);
159   fData->AddHistogram(fHistoElectrons);
160   fData->RefSliceInfo(2).Setup("Elect.:", 0, fElectro);
161   fData->AddHistogram(fHistoMuons);
162   fData->RefSliceInfo(3).Setup("Muons:", 0, fMuon);
163   fData->AddHistogram(fHistoPions);
164   fData->RefSliceInfo(4).Setup("Pions:", 0, fPion);
165   fData->AddHistogram(fHistoKaons);
166   fData->RefSliceInfo(5).Setup("Kaons:", 0, fKaon);
167   fData->AddHistogram(fHistoProtons);
168   fData->RefSliceInfo(6).Setup("Proto.:", 0, fProton);
169
170   fData->GetEtaBins()->SetTitleFont(120);
171   fData->GetEtaBins()->SetTitle("h");
172   fData->GetPhiBins()->SetTitleFont(120);
173   fData->GetPhiBins()->SetTitle("f");
174   fData->IncDenyDestroy();
175
176   // Setting up the position of the 3D calorimeter histogram view.
177   fCalo3d = new TEveCalo3D(fData);
178   fCalo3d->SetBarrelRadius(550);
179   fCalo3d->SetEndCapPos(550);
180
181   // Adding data to the Lego object
182   fLego = new TEveCaloLego(fData);
183
184   // Creating projections
185   fAl = AliEveMultiView::Instance();
186   fAl->ImportEventRPhi(fCalo3d);
187   fAl->ImportEventRhoZ(fCalo3d);
188
189   // Update viewers and scenes
190   Update();
191 }
192
193 //______________________________________________________________________________
194 AliEveLego::~AliEveLego()
195 {
196    // Deleting variables
197    delete fEsd;
198    delete fPhysicsSelection;
199    delete fHistopos;
200    delete fHistoposAllEvents;
201    delete fHistoneg;
202    delete fHistonegAllEvents;
203    delete fHistoElectrons;
204    delete fHistoElectronsAllEvents;
205    delete fHistoMuons;
206    delete fHistoMuonsAllEvents;
207    delete fHistoPions;
208    delete fHistoPionsAllEvents;
209    delete fHistoKaons;
210    delete fHistoKaonsAllEvents;
211    delete fHistoProtons;
212    delete fHistoProtonsAllEvents;
213    delete[] fParticleTypeId;
214    delete[] fParticleTypeIdAE;
215    delete fData;
216    delete fDataAllEvents;
217    delete fLego;
218    delete fLegoAllEvents;
219    delete fCalo3d;
220    delete fCalo3dAllEvents;
221    delete fGlv;
222    delete fHisto2dv;
223    delete fHisto2ds;
224    delete fHisto2ds2;
225    delete fHisto2dAllEventsv0;
226    delete fHisto2dAllEventsv1;
227    delete fHisto2dAllEventsv2;
228    delete fHisto2dAllEventsv3;
229    delete fHisto2dAllEventss0;
230    delete fHisto2dAllEventss1;
231    delete fHisto2dAllEventss2;
232    delete fHisto2dAllEventss3;
233    delete fAl;
234    delete fHisto2dLegoOverlay;
235    delete fHisto2dAllEventsLegoOverlay;
236    delete fHisto2dAllEventsSlot;
237 }
238
239 namespace
240 {
241   //____________________________________________________________________________
242   Double_t getphi(Double_t phi)
243   {
244     // Phi correction for alice
245
246     if (phi > TMath::Pi()) {
247       phi -= TMath::TwoPi();
248     }
249     return phi;
250   }
251 }
252
253 //______________________________________________________________________________
254 TEveCaloDataHist* AliEveLego::LoadData()
255 {
256    // Load data from ESD tree
257    fHistopos->Reset();
258    fHistoneg->Reset();
259    fHistoElectrons->Reset();
260    fHistoMuons->Reset();
261    fHistoPions->Reset();
262    fHistoKaons->Reset();
263    fHistoProtons->Reset();
264
265    // Getting current tracks, filling histograms
266    for (int n = 0; n < fEsd->GetNumberOfTracks(); ++n) {
267
268      AliESDtrack *track = fEsd->GetTrack(n);
269      const Double_t sign = fEsd->GetTrack(n)->GetSign();
270      const Int_t prob = GetParticleType(track);     
271
272      // Filling histograms
273      if (sign > 0)
274        fHistopos->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
275
276      if (sign < 0)
277        fHistoneg->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
278
279      if (prob == 0)
280        fHistoElectrons->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
281
282      if (prob == 1)
283        fHistoMuons->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
284
285      if (prob == 2)
286        fHistoPions->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
287
288      if (prob == 3)
289        fHistoKaons->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
290
291      if (prob == 4)
292        fHistoProtons->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
293    }
294
295    fData->DataChanged();
296
297    FilterData();
298
299    return fData;
300 }
301
302 //______________________________________________________________________________
303 TEveCaloDataHist* AliEveLego::LoadAllData()
304 {
305    // Load data from all events ESD
306    fHistoposAllEvents->Reset();
307    fHistonegAllEvents->Reset();
308    fHistoElectronsAllEvents->Reset();
309    fHistoMuonsAllEvents->Reset();
310    fHistoPionsAllEvents->Reset();
311    fHistoKaonsAllEvents->Reset();
312    fHistoProtonsAllEvents->Reset();
313
314    TTree* t = AliEveEventManager::GetMaster()->GetESDTree();
315
316    // Getting current tracks for each event, filling histograms
317    Int_t fAcceptedEvents = 0;
318    for (int event = 0; event < t->GetEntries(); event++) {
319       t->GetEntry(event);
320
321       if (fCollisionCandidatesOnly == kTRUE)
322         if (fPhysicsSelection->IsCollisionCandidate(fEsd) == kFALSE) continue;
323
324       fAcceptedEvents++;
325
326       for (int n = 0; n < fEsd->GetNumberOfTracks(); ++n) {
327
328            AliESDtrack *track = fEsd->GetTrack(n);
329            const Double_t sign = track->GetSign();
330            const Int_t prob = GetParticleType(track);
331
332            if (sign > 0)
333              fHistoposAllEvents->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
334
335            if (sign < 0)
336              fHistonegAllEvents->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
337
338            if (prob == 0)
339              fHistoElectronsAllEvents->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
340
341            if (prob == 1)
342              fHistoMuonsAllEvents->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
343
344            if (prob == 2)
345              fHistoPionsAllEvents->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
346
347            if (prob == 3)
348              fHistoKaonsAllEvents->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
349
350            if (prob == 4)
351              fHistoProtonsAllEvents->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
352         }
353    }
354
355    // Setting the current view to the first event
356    t->GetEntry(0);
357
358    // Usefull information,
359    // with this we can estimate the event efficiency
360    printf("Number of events loaded: %i, with AliPhysicsSelection: %i\n",fAcceptedEvents,fCollisionCandidatesOnly);
361
362    fDataAllEvents->DataChanged();
363
364    return fDataAllEvents;
365 }
366
367 //______________________________________________________________________________
368 TEveCaloDataHist* AliEveLego::FilterData()
369 {
370    // Tracks selection
371    if ( fTracksId == 2 )
372    {
373       fHistopos->Reset();
374       fHistoneg->Reset();
375       fHistoElectrons->Reset();
376       fHistoMuons->Reset();
377       fHistoPions->Reset();
378       fHistoKaons->Reset();
379       fHistoProtons->Reset();
380
381       const AliESDVertex *pv  = fEsd->GetPrimaryVertex();
382
383       for (Int_t n = 0; n < pv->GetNIndices(); n++ )
384       {
385          AliESDtrack *at = fEsd->GetTrack(pv->GetIndices()[n]);
386          const Double_t sign = at->GetSign();
387          const Int_t prob = GetParticleType(at);
388
389          if (sign > 0)
390            fHistopos->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
391
392          if (sign < 0)
393            fHistoneg->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
394
395          if (prob == 0)
396            fHistoElectrons->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
397
398          if (prob == 1)
399            fHistoMuons->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
400
401          if (prob == 2)
402            fHistoPions->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
403
404          if (prob == 3)
405            fHistoKaons->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
406
407          if (prob == 4)
408            fHistoProtons->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
409       }
410    }
411
412    fData->DataChanged();
413
414    // Max Pt threshold
415    if (GetPtMax() >= fMaxPt){
416       for (Int_t binx = 1; binx <= 100; binx++) {
417          for (Int_t biny = 1; biny <= 80; biny++) {
418
419            if (fHistopos->GetBinContent(binx, biny) >= fMaxPt)
420              fHistopos->SetBinContent(binx, biny, fMaxPt);
421
422            if (fHistoneg->GetBinContent(binx, biny) >= fMaxPt)
423              fHistoneg->SetBinContent(binx, biny, fMaxPt);
424
425            if (fHistoElectrons->GetBinContent(binx, biny) >= fMaxPt)
426              fHistoElectrons->SetBinContent(binx, biny, fMaxPt);
427
428            if (fHistoMuons->GetBinContent(binx, biny) >= fMaxPt)
429              fHistoMuons->SetBinContent(binx, biny, fMaxPt);
430
431            if (fHistoPions->GetBinContent(binx, biny) >= fMaxPt)
432              fHistoPions->SetBinContent(binx, biny, fMaxPt);
433
434            if (fHistoKaons->GetBinContent(binx, biny) >= fMaxPt)
435              fHistoKaons->SetBinContent(binx, biny, fMaxPt);
436
437            if (fHistoProtons->GetBinContent(binx, biny) >= fMaxPt)
438              fHistoProtons->SetBinContent(binx, biny, fMaxPt);
439          }
440       }
441    }
442
443    // Particle type filter
444    if ( fParticleTypeId[0] == kFALSE) fHistopos->Reset();
445    if ( fParticleTypeId[1] == kFALSE) fHistoneg->Reset();
446    if ( fParticleTypeId[2] == kFALSE) fHistoElectrons->Reset();
447    if ( fParticleTypeId[3] == kFALSE) fHistoMuons->Reset();
448    if ( fParticleTypeId[4] == kFALSE) fHistoPions->Reset();
449    if ( fParticleTypeId[5] == kFALSE) fHistoKaons->Reset();
450    if ( fParticleTypeId[6] == kFALSE) fHistoProtons->Reset();
451
452    fData->DataChanged();
453
454    return fData;
455 }
456
457 //______________________________________________________________________________
458 TEveCaloDataHist* AliEveLego::FilterAllData()
459 {
460    // Tracks selection
461    if ( fTracksIdAE == 2 )
462    {
463       fHistoposAllEvents->Reset();
464       fHistonegAllEvents->Reset();
465       fHistoElectronsAllEvents->Reset();
466       fHistoMuonsAllEvents->Reset();
467       fHistoPionsAllEvents->Reset();
468       fHistoKaonsAllEvents->Reset();
469       fHistoProtonsAllEvents->Reset();
470
471       TTree* t = AliEveEventManager::GetMaster()->GetESDTree();
472
473       // Getting current tracks for each event, filling histograms
474       Int_t fAcceptedEvents = 0;
475       for (int event = 0; event < t->GetEntries(); event++) {
476
477         t->GetEntry(event);
478         const AliESDVertex *pv  = fEsd->GetPrimaryVertex();
479
480         if (fCollisionCandidatesOnly == kTRUE)
481           if (fPhysicsSelection->IsCollisionCandidate(fEsd) == kFALSE) continue;
482
483         fAcceptedEvents++;
484
485         for (Int_t n = 0; n < pv->GetNIndices(); n++ )
486         {
487           AliESDtrack *at = fEsd->GetTrack(pv->GetIndices()[n]);
488           const Double_t sign = at->GetSign();
489           const Int_t prob = GetParticleType(at);
490
491           if (sign > 0)
492             fHistoposAllEvents->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
493
494           if (sign < 0)
495             fHistonegAllEvents->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
496
497           if (prob == 0)
498             fHistoElectronsAllEvents->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
499
500           if (prob == 1)
501             fHistoMuonsAllEvents->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
502
503           if (prob == 2)
504             fHistoPionsAllEvents->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
505
506           if (prob == 3)
507             fHistoKaonsAllEvents->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
508
509           if (prob == 4)
510             fHistoProtonsAllEvents->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
511         }
512       }
513       t->GetEntry(0);
514       printf("Number of events loaded: %i, with AliPhysicsSelection: %i\n",fAcceptedEvents,fCollisionCandidatesOnly);
515
516    } else {
517
518       LoadAllData();
519
520     }
521    
522    fDataAllEvents->DataChanged();
523
524    // Max Pt threshold
525    if (GetPtMaxAE() >= fMaxPtAE){
526       for (Int_t binx = 1; binx <= 100; binx++) {
527          for (Int_t biny = 1; biny <= 80; biny++) {
528
529             if (fHistoposAllEvents->GetBinContent(binx, biny) >= fMaxPtAE)            
530               fHistoposAllEvents->SetBinContent(binx, biny, fMaxPtAE);
531
532             if (fHistonegAllEvents->GetBinContent(binx, biny) >= fMaxPtAE)            
533               fHistonegAllEvents->SetBinContent(binx, biny, fMaxPtAE);
534
535             if (fHistoElectronsAllEvents->GetBinContent(binx, biny) >= fMaxPt)
536               fHistoElectronsAllEvents->SetBinContent(binx, biny, fMaxPt);
537
538             if (fHistoMuonsAllEvents->GetBinContent(binx, biny) >= fMaxPt)
539               fHistoMuonsAllEvents->SetBinContent(binx, biny, fMaxPt);
540
541             if (fHistoPionsAllEvents->GetBinContent(binx, biny) >= fMaxPt)
542               fHistoPionsAllEvents->SetBinContent(binx, biny, fMaxPt);
543
544             if (fHistoKaonsAllEvents->GetBinContent(binx, biny) >= fMaxPt)
545               fHistoKaonsAllEvents->SetBinContent(binx, biny, fMaxPt);
546
547             if (fHistoProtonsAllEvents->GetBinContent(binx, biny) >= fMaxPt)
548               fHistoProtonsAllEvents->SetBinContent(binx, biny, fMaxPt);
549          }
550       }
551    }
552
553    // Particles species and charges filter
554    if ( fParticleTypeIdAE[0] == kFALSE ) fHistoposAllEvents->Reset();
555    if ( fParticleTypeIdAE[1] == kFALSE ) fHistonegAllEvents->Reset();
556    if ( fParticleTypeIdAE[2] == kFALSE ) fHistoElectronsAllEvents->Reset();
557    if ( fParticleTypeIdAE[3] == kFALSE ) fHistoMuonsAllEvents->Reset();
558    if ( fParticleTypeIdAE[4] == kFALSE ) fHistoPionsAllEvents->Reset();
559    if ( fParticleTypeIdAE[5] == kFALSE ) fHistoKaonsAllEvents->Reset();
560    if ( fParticleTypeIdAE[6] == kFALSE ) fHistoProtonsAllEvents->Reset();
561
562    fDataAllEvents->DataChanged();
563
564    gEve->Redraw3D(kTRUE);
565
566    return fDataAllEvents;
567 }
568
569 //______________________________________________________________________________
570 void AliEveLego::Update()
571 {
572    // Load/Reload data
573    LoadData();
574
575    // Create new histo2d
576    CreateHistoLego();
577
578    // Create 3d view
579    Create3DView();
580
581    // Update the viewers
582    gEve->Redraw3D(kTRUE);
583 }
584
585 //______________________________________________________________________________
586 TEveCaloLego* AliEveLego::CreateHistoLego()
587 {
588    // Viewer initialization, tab creation
589    if (fHisto2dv == 0) {
590       TEveWindowSlot *fslot    = 0;
591       TEveBrowser    *fbrowser = gEve->GetBrowser();
592
593       fslot = TEveWindow::CreateWindowInTab(fbrowser->GetTabRight());
594       fslot->MakeCurrent();
595       fHisto2dv = gEve->SpawnNewViewer("2D Lego Histogram", "2D Lego Histogram");
596       fHisto2ds = gEve->SpawnNewScene("2D Lego Histogram", "2D Lego Histogram");
597       fHisto2dv->AddScene(fHisto2ds);
598       fHisto2dv->SetElementName("2D Lego Viewer");
599       fHisto2ds->SetElementName("2D Lego Scene");
600
601       fGlv = fHisto2dv->GetGLViewer();
602       fHisto2dLegoOverlay = new TEveCaloLegoOverlay();
603       fGlv->AddOverlayElement(fHisto2dLegoOverlay);
604       fGlv->SetCurrentCamera(TGLViewer::kCameraPerspXOY);
605
606       fHisto2ds->AddElement(fLego);
607
608       // move to real world coordinates
609       fLego->InitMainTrans();
610       Float_t sc = TMath::Min(fLego->GetEtaRng(), fLego->GetPhiRng());
611       fLego->RefMainTrans().SetScale(sc, sc, sc);
612
613       // set event handler to move from perspective to orthographic view.
614       fGlv->SetEventHandler(new TEveLegoEventHandler(fGlv->GetGLWidget(), fGlv, fLego));
615
616       fHisto2dLegoOverlay->SetCaloLego(fLego);
617    }
618
619    return fLego;
620 }
621
622 //______________________________________________________________________________
623 TEveCaloLego* AliEveLego::CreateHistoLego(TEveWindowSlot *slot)
624 {
625    // Viewer initialization, tab creation
626    if (fHisto2dAllEventsv0 == 0) {
627
628       slot->MakeCurrent();
629       fHisto2dAllEventsv0 = gEve->SpawnNewViewer("2D Lego Histogram", "2D Lego Histogram");
630       fHisto2dAllEventss0 = gEve->SpawnNewScene("2D Lego Histogram", "2D Lego Histogram");
631       fHisto2dAllEventsv0->AddScene(fHisto2dAllEventss0);
632       fHisto2dAllEventsv0->SetElementName("2D Lego Viewer");
633       fHisto2dAllEventss0->SetElementName("2D Lego Scene");
634
635       TGLViewer* glv = fHisto2dAllEventsv0->GetGLViewer();
636       fHisto2dAllEventsLegoOverlay = new TEveCaloLegoOverlay();
637       glv->AddOverlayElement(fHisto2dAllEventsLegoOverlay);
638       glv->SetCurrentCamera(TGLViewer::kCameraPerspXOY);
639
640       // Plotting histogram lego
641       fLegoAllEvents = new TEveCaloLego(fDataAllEvents);
642       fHisto2dAllEventss0->AddElement(fLegoAllEvents);
643
644       // Move to real world coordinates
645       fLegoAllEvents->InitMainTrans();
646       Float_t sc = TMath::Min(fLegoAllEvents->GetEtaRng(), fLegoAllEvents->GetPhiRng());
647       fLegoAllEvents->RefMainTrans().SetScale(sc, sc, sc);
648
649       // Set event handler to move from perspective to orthographic view.
650       glv->SetEventHandler(new TEveLegoEventHandler(glv->GetGLWidget(), glv, fLegoAllEvents));
651
652       fHisto2dAllEventsLegoOverlay->SetCaloLego(fLegoAllEvents);
653    }
654
655    return fLegoAllEvents;
656 }
657
658 //______________________________________________________________________________
659 TEveCalo3D* AliEveLego::Create3DView()
660 {
661    // Initialization
662    if (fHisto2ds2 == 0) {
663       fHisto2ds2 = gEve->SpawnNewScene("3D Histogram", "3D Histogram");
664       gEve->GetDefaultViewer()->AddScene(fHisto2ds2);
665       fHisto2ds2->SetElementName("3D Histogram Scene");
666       fHisto2ds2->AddElement(fCalo3d);
667    }
668
669    return fCalo3d;
670 }
671
672 //______________________________________________________________________________
673 TEveCalo3D* AliEveLego::Create3DView(TEveWindowSlot *slot)
674 {
675    // Creates a 3d view for the 3d histogram
676    if ( fHisto2dAllEventsv1 == 0 ) {
677
678       slot->MakeCurrent();
679       fHisto2dAllEventsv1 = gEve->SpawnNewViewer("3D Histogram", "3D Histogram");
680       fHisto2dAllEventss1 = gEve->SpawnNewScene("3D Histogram", "3D Histogram");
681       fHisto2dAllEventsv1->AddScene(fHisto2dAllEventss1);
682       fHisto2dAllEventsv1->SetElementName("3D Histogram Viewer");
683       fHisto2dAllEventss1->SetElementName("3D Histogram Scene");
684
685       fCalo3dAllEvents = new TEveCalo3D(fDataAllEvents);
686
687       fCalo3dAllEvents->SetBarrelRadius(550);
688       fCalo3dAllEvents->SetEndCapPos(550);
689       fHisto2dAllEventss1->AddElement(fCalo3dAllEvents);
690    }
691
692    return fCalo3dAllEvents;
693 }
694
695 //______________________________________________________________________________
696 void AliEveLego::CreateProjections(TEveWindowSlot* slot1, TEveWindowSlot* slot2)
697 {
698    // Create projections
699    if (fHisto2dAllEventsv2 == 0) {
700
701       slot1->MakeCurrent();
702       fHisto2dAllEventsv2 = gEve->SpawnNewViewer("RPhi projection", "RPhi projection");
703       fHisto2dAllEventss2 = gEve->SpawnNewScene("RPhi projection", "RPhi projection");
704       fHisto2dAllEventsv2->AddScene(fHisto2dAllEventss2);
705       fHisto2dAllEventsv2->SetElementName("RPhi Projection Viewer");
706       fHisto2dAllEventss2->SetElementName("RPhi Projection Scene");
707
708       TEveProjectionManager* mng1 = new TEveProjectionManager();
709       mng1->SetProjection(TEveProjection::kPT_RPhi);
710
711       TEveProjectionAxes* axeghisto2dAllEventss1 = new TEveProjectionAxes(mng1);
712       fHisto2dAllEventss2->AddElement(axeghisto2dAllEventss1);
713       TEveCalo2D* fcalo2d1 = (TEveCalo2D*) mng1->ImportElements(fCalo3dAllEvents);
714       fHisto2dAllEventss2->AddElement(fcalo2d1);
715
716       fHisto2dAllEventsv2->GetGLViewer()->SetCurrentCamera(TGLViewer::kCameraOrthoXOY);
717    }
718
719    if (fHisto2dAllEventsv3 == 0)
720    {
721       slot2->MakeCurrent();
722       fHisto2dAllEventsv3 = gEve->SpawnNewViewer("RhoZ projection", "RhoZ projection");
723       fHisto2dAllEventss3 = gEve->SpawnNewScene("RhoZ projection", "RhoZ projection");
724       fHisto2dAllEventsv3->AddScene(fHisto2dAllEventss3);
725       fHisto2dAllEventsv3->SetElementName("RhoZ Projection Viewer");
726       fHisto2dAllEventss3->SetElementName("RhoZ Projection Viewer");
727
728       TEveProjectionManager* mng2 = new TEveProjectionManager();
729       mng2->SetProjection(TEveProjection::kPT_RhoZ);
730
731       TEveProjectionAxes* axeghisto2dAllEventss2 = new TEveProjectionAxes(mng2);
732       fHisto2dAllEventss3->AddElement(axeghisto2dAllEventss2);
733       TEveCalo2D* fcalo2d2 = (TEveCalo2D*) mng2->ImportElements(fCalo3dAllEvents);
734       fHisto2dAllEventss3->AddElement(fcalo2d2);
735
736       fHisto2dAllEventsv3->GetGLViewer()->SetCurrentCamera(TGLViewer::kCameraOrthoXOY);
737    }
738    return;
739 }
740
741 //______________________________________________________________________________
742 TEveCaloDataHist* AliEveLego::LoadAllEvents()
743 {
744    // load all events data from ESD
745    if ( fHisto2dAllEventsSlot == 0 ) {
746
747       printf("Filling histogram...\n");
748       TStopwatch timer;
749       timer.Start();
750
751       // Creating 2D histograms
752       fHistoposAllEvents = new TH2F("fHistoposAllEvents","Histo 2d positive",
753                                  100,-1.5,1.5,80,-kPi,kPi);
754       fHistonegAllEvents = new TH2F("fHistonegAllEvents","Histo 2d negative",
755                                  100,-1.5,1.5,80,-kPi,kPi);
756       fHistoElectronsAllEvents = new TH2F("fHistoElectronsAllEvents","Histo 2d electrons",
757                                        100,-1.5,1.5,80,-kPi,kPi);
758       fHistoMuonsAllEvents = new TH2F("fHistoMuonsAllEvents","Histo 2d muons",
759                                        100,-1.5,1.5,80,-kPi,kPi);
760       fHistoPionsAllEvents = new TH2F("fHistoPionsAllEvents","Histo 2d pions",
761                                        100,-1.5,1.5,80,-kPi,kPi);
762       fHistoKaonsAllEvents = new TH2F("fHistoKaonsAllEvents","Histo 2d kaons",
763                                        100,-1.5,1.5,80,-kPi,kPi);
764       fHistoProtonsAllEvents = new TH2F("fHistoProtonsAllEvents","Histo 2d protons",
765                                        100,-1.5,1.5,80,-kPi,kPi);
766
767       fHistoposAllEvents->SetDirectory(0);
768       fHistonegAllEvents->SetDirectory(0);
769       fHistoElectronsAllEvents->SetDirectory(0);
770       fHistoMuonsAllEvents->SetDirectory(0);
771       fHistoPionsAllEvents->SetDirectory(0);
772       fHistoKaonsAllEvents->SetDirectory(0);
773       fHistoProtonsAllEvents->SetDirectory(0);
774
775       // colors from get_pdg_color() in /alice-macros/kine_tracks.C
776       static Color_t fElectro = 5;
777       static Color_t fMuon = 30;
778       static Color_t fPion = 3;
779       static Color_t fKaon = 38;
780       static Color_t fProton = 10;
781
782       fDataAllEvents = new TEveCaloDataHist();
783       fDataAllEvents->AddHistogram(fHistonegAllEvents);
784       fDataAllEvents->RefSliceInfo(0).Setup("NegCg:", 0, kBlue);
785       fDataAllEvents->AddHistogram(fHistoposAllEvents);
786       fDataAllEvents->RefSliceInfo(1).Setup("PosCg:", 0, kRed);
787       fDataAllEvents->AddHistogram(fHistoElectronsAllEvents);
788       fDataAllEvents->RefSliceInfo(2).Setup("Electrons:", 0, fElectro);
789       fDataAllEvents->AddHistogram(fHistoMuonsAllEvents);
790       fDataAllEvents->RefSliceInfo(3).Setup("Muons:", 0, fMuon);
791       fDataAllEvents->AddHistogram(fHistoPionsAllEvents);
792       fDataAllEvents->RefSliceInfo(4).Setup("Pions:", 0, fPion);
793       fDataAllEvents->AddHistogram(fHistoKaonsAllEvents);
794       fDataAllEvents->RefSliceInfo(5).Setup("Kaons:", 0, fKaon);
795       fDataAllEvents->AddHistogram(fHistoProtonsAllEvents);
796       fDataAllEvents->RefSliceInfo(6).Setup("Protons:", 0, fProton);
797
798       fDataAllEvents->GetEtaBins()->SetTitleFont(120);
799       fDataAllEvents->GetEtaBins()->SetTitle("h");
800       fDataAllEvents->GetPhiBins()->SetTitleFont(120);
801       fDataAllEvents->GetPhiBins()->SetTitle("f");
802       fDataAllEvents->IncDenyDestroy();
803
804       // Creating frames
805       fHisto2dAllEventsSlot = TEveWindow::CreateWindowInTab(gEve->GetBrowser()->GetTabRight());
806       TEveWindowPack* packH = fHisto2dAllEventsSlot->MakePack();
807       packH->SetElementName("Projections");
808       packH->SetHorizontal();
809       packH->SetShowTitleBar(kFALSE);
810
811       fHisto2dAllEventsSlot = packH->NewSlot();
812       TEveWindowPack* pack0 = fHisto2dAllEventsSlot->MakePack();
813       pack0->SetShowTitleBar(kFALSE);
814       TEveWindowSlot*  slotLeftTop   = pack0->NewSlot();
815       TEveWindowSlot* slotLeftBottom = pack0->NewSlot();
816
817       fHisto2dAllEventsSlot = packH->NewSlot();
818       TEveWindowPack* pack1 = fHisto2dAllEventsSlot->MakePack();
819       pack1->SetShowTitleBar(kFALSE);
820       TEveWindowSlot* slotRightTop    = pack1->NewSlot();
821       TEveWindowSlot* slotRightBottom = pack1->NewSlot();
822
823       // Creating viewers and scenes
824       Create3DView(slotLeftTop);
825       CreateHistoLego(slotLeftBottom);
826       CreateProjections(slotRightTop, slotRightBottom);
827
828       FilterAllData();
829
830       gEve->Redraw3D(kTRUE);
831
832       printf("Filling histogram... Finished\n");
833       timer.Stop();
834       timer.Print();
835    }
836    return fDataAllEvents;
837 }
838
839 //______________________________________________________________________________
840 void AliEveLego::ApplyParticleTypeSelectionAE()
841 {
842   // Reload all events applying particle type selection
843   FilterAllData();
844 }
845
846 //______________________________________________________________________________
847 Float_t AliEveLego::GetPtMax()
848 {
849    // Return pT maximum
850    return fData->GetMaxVal(fLego->GetPlotEt());
851 }
852
853 //______________________________________________________________________________
854 Float_t AliEveLego::GetPtMaxAE()
855 {
856    // Return pT max from all events
857    return fDataAllEvents->GetMaxVal(fLegoAllEvents->GetPlotEt());
858 }
859
860 //______________________________________________________________________________
861 Int_t AliEveLego::GetParticleType(AliESDtrack *track)
862 {
863   // Determine the particle type
864         Double_t prob[5]={0.};
865   track->GetESDpid(prob);
866
867   Double_t max = prob[0];
868   Int_t index = 0;
869
870   for (Int_t i = 1 ; i < 5; i++)
871   {
872     if (prob[i] > max){
873       max = prob[i];
874       index = i;
875     }
876   }
877   return index;
878 }
879
880 //______________________________________________________________________________
881 void AliEveLego::SetParticleType(Int_t id, Bool_t status)
882 {
883   // Activate/deactivate particles types
884   fParticleTypeId[id] = status;
885
886   Update();
887 }
888
889 //______________________________________________________________________________
890 void AliEveLego::SetParticleTypeAE(Int_t id, Bool_t status)
891 {
892   // Activate/deactivate particles types
893   fParticleTypeIdAE[id] = status;
894 }
895
896 //______________________________________________________________________________
897 void AliEveLego::SetMaxPt(Double_t val)
898 {
899    // Add new maximum
900    fMaxPt = val;   
901    Update();
902 }
903
904 //______________________________________________________________________________
905 void AliEveLego::SetMaxPtAE(Double_t val)
906 {
907    // Add new maximum
908    fMaxPtAE = val;
909    FilterAllData();
910 }
911
912 //______________________________________________________________________________
913 void AliEveLego::SetThreshold(Double_t val)
914 {  
915    // Setting up the new threshold for all histograms
916    fData->SetSliceThreshold(0,val);
917    fData->SetSliceThreshold(1,val);
918    fData->SetSliceThreshold(2,val);
919    fData->SetSliceThreshold(3,val);
920    fData->SetSliceThreshold(4,val);
921    fData->SetSliceThreshold(5,val);
922    fData->SetSliceThreshold(6,val);
923    fData->DataChanged();
924
925    gEve->Redraw3D(kTRUE);
926 }
927
928 //______________________________________________________________________________
929 void AliEveLego::SetThresholdAE(Double_t val)
930 {
931    // Setting up the new threshold for all histograms
932    fDataAllEvents->SetSliceThreshold(0,val);
933    fDataAllEvents->SetSliceThreshold(1,val);
934    fDataAllEvents->SetSliceThreshold(2,val);
935    fDataAllEvents->SetSliceThreshold(3,val);
936    fDataAllEvents->SetSliceThreshold(4,val);
937    fDataAllEvents->SetSliceThreshold(5,val);
938    fDataAllEvents->SetSliceThreshold(6,val);
939    fDataAllEvents->DataChanged();
940
941    gEve->Redraw3D(kTRUE);
942 }
943
944 //______________________________________________________________________________
945 void AliEveLego::SwitchDataType(Bool_t status)
946 {
947   // Activate/deactivate MC / real data type
948   fIsMC = status;
949
950   // Removing defaul physics selection
951   delete fPhysicsSelection;
952   fPhysicsSelection = NULL;
953
954   // Re-initialization of physics selection
955   fPhysicsSelection = new AliPhysicsSelection();
956   fPhysicsSelection->SetAnalyzeMC(fIsMC);
957   fPhysicsSelection->Initialize(fEsd);
958   FilterAllData();
959 }
960
961 //______________________________________________________________________________
962 void AliEveLego::SetCollisionCandidatesOnly()
963 {
964   // Activate/deactivate MC / real data type
965   if (fCollisionCandidatesOnly == 0)
966   {
967     fCollisionCandidatesOnly = 1;
968   } else {
969     fCollisionCandidatesOnly = 0;
970   }
971   FilterAllData();
972 }
973
974 /******************************************************************************/
975
976