2 // Author: Stefano Carrazza 2010
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 **************************************************************************/
10 #include "AliESDEvent.h"
11 #include "AliEveLego.h"
12 #include "AliEveEventManager.h"
13 #include "AliEveMultiView.h"
14 #include "AliPhysicsSelection.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"
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"
34 //______________________________________________________________________________
35 // Allow histograms visualization in 2D and 3D.
39 Double_t kPi = TMath::Pi();
41 //______________________________________________________________________________
42 AliEveLego::AliEveLego(const char* name) :
43 TEveElementList(name),
45 fCollisionCandidatesOnly(kFALSE),
55 fHistoposAllEvents(0),
57 fHistonegAllEvents(0),
59 fHistoElectronsAllEvents(0),
61 fHistoMuonsAllEvents(0),
63 fHistoPionsAllEvents(0),
65 fHistoKaonsAllEvents(0),
67 fHistoProtonsAllEvents(0),
78 fHisto2dAllEventsv0(0),
79 fHisto2dAllEventsv1(0),
80 fHisto2dAllEventsv2(0),
81 fHisto2dAllEventsv3(0),
82 fHisto2dAllEventss0(0),
83 fHisto2dAllEventss1(0),
84 fHisto2dAllEventss2(0),
85 fHisto2dAllEventss3(0),
87 fHisto2dLegoOverlay(0),
88 fHisto2dAllEventsLegoOverlay(0),
89 fHisto2dAllEventsSlot(0)
92 gEve->AddToListTree(this,0);
94 // Get Current ESD event
95 fEsd = AliEveEventManager::AssertESD();
97 // Particles types per default
98 fParticleTypeId = new Bool_t[7];
99 fParticleTypeIdAE = new Bool_t[7];
101 for (Int_t s = 0; s < 7; s++)
104 fParticleTypeId[s] = kFALSE;
105 fParticleTypeIdAE[s] = kFALSE;
107 fParticleTypeId[s] = kTRUE;
108 fParticleTypeIdAE[s] = kTRUE;
112 fPhysicsSelection = new AliPhysicsSelection();
113 fPhysicsSelection->Initialize(fEsd->GetRunNumber());
115 fHistopos = new TH2F("histopos","Histo 2d positive", 100, -1.5, 1.5, 80, -kPi, kPi);
116 fHistoneg = new TH2F("histoneg","Histo 2d negative", 100, -1.5, 1.5, 80, -kPi, kPi);
117 fHistoElectrons = new TH2F("histoele","Histo 2d electron", 100, -1.5, 1.5, 80, -kPi, kPi);
118 fHistoMuons = new TH2F("histomuo","Histo 2d muons ", 100, -1.5, 1.5, 80, -kPi, kPi);
119 fHistoPions = new TH2F("histopio","Histo 2d pions ", 100, -1.5, 1.5, 80, -kPi, kPi);
120 fHistoKaons = new TH2F("histokao","Histo 2d kaons ", 100, -1.5, 1.5, 80, -kPi, kPi);
121 fHistoProtons = new TH2F("histopro","Histo 2d protons ", 100, -1.5, 1.5, 80, -kPi, kPi);
123 fHistopos->SetDirectory(0);
124 fHistoneg->SetDirectory(0);
125 fHistoElectrons->SetDirectory(0);
126 fHistoMuons->SetDirectory(0);
127 fHistoPions->SetDirectory(0);
128 fHistoKaons->SetDirectory(0);
129 fHistoProtons->SetDirectory(0);
131 // colors from get_pdg_color() in /alice-macros/kine_tracks.C
132 static Color_t fElectro = 5;
133 static Color_t fMuon = 30;
134 static Color_t fPion = 3;
135 static Color_t fKaon = 38;
136 static Color_t fProton = 10;
138 fData = new TEveCaloDataHist();
140 fData->AddHistogram(fHistoneg);
141 fData->RefSliceInfo(0).Setup("NegCg:", 0, kBlue);
142 fData->AddHistogram(fHistopos);
143 fData->RefSliceInfo(1).Setup("PosCg:", 0, kRed);
144 fData->AddHistogram(fHistoElectrons);
145 fData->RefSliceInfo(2).Setup("Electrons:", 0, fElectro);
146 fData->AddHistogram(fHistoMuons);
147 fData->RefSliceInfo(3).Setup("Muons:", 0, fMuon);
148 fData->AddHistogram(fHistoPions);
149 fData->RefSliceInfo(4).Setup("Pions:", 0, fPion);
150 fData->AddHistogram(fHistoKaons);
151 fData->RefSliceInfo(5).Setup("Kaons:", 0, fKaon);
152 fData->AddHistogram(fHistoProtons);
153 fData->RefSliceInfo(6).Setup("Protons:", 0, fProton);
155 fData->GetEtaBins()->SetTitleFont(120);
156 fData->GetEtaBins()->SetTitle("h");
157 fData->GetPhiBins()->SetTitleFont(120);
158 fData->GetPhiBins()->SetTitle("f");
159 fData->IncDenyDestroy();
161 fCalo3d = new TEveCalo3D(fData);
162 fCalo3d->SetBarrelRadius(550);
163 fCalo3d->SetEndCapPos(550);
166 fLego = new TEveCaloLego(fData);
169 fAl = AliEveMultiView::Instance();
170 fAl->ImportEventRPhi(fCalo3d);
171 fAl->ImportEventRhoZ(fCalo3d);
177 //______________________________________________________________________________
178 AliEveLego::~AliEveLego()
180 // deleting variables
182 delete fPhysicsSelection;
184 delete fHistoposAllEvents;
186 delete fHistonegAllEvents;
187 delete fHistoElectrons;
188 delete fHistoElectronsAllEvents;
190 delete fHistoMuonsAllEvents;
192 delete fHistoPionsAllEvents;
194 delete fHistoKaonsAllEvents;
195 delete fHistoProtons;
196 delete fHistoProtonsAllEvents;
198 delete[] fParticleTypeId;
199 delete[] fParticleTypeIdAE;
202 delete fDataAllEvents;
204 delete fLegoAllEvents;
206 delete fCalo3dAllEvents;
212 delete fHisto2dAllEventsv0;
213 delete fHisto2dAllEventsv1;
214 delete fHisto2dAllEventsv2;
215 delete fHisto2dAllEventsv3;
216 delete fHisto2dAllEventss0;
217 delete fHisto2dAllEventss1;
218 delete fHisto2dAllEventss2;
219 delete fHisto2dAllEventss3;
222 delete fHisto2dLegoOverlay;
223 delete fHisto2dAllEventsLegoOverlay;
224 delete fHisto2dAllEventsSlot;
230 //____________________________________________________________________________
231 Double_t getphi(Double_t phi)
233 // phi correction for alice
235 if (phi > TMath::Pi()) {
236 phi -= TMath::TwoPi();
242 //______________________________________________________________________________
243 TEveCaloDataHist* AliEveLego::LoadData()
245 // Load data from ESD tree
248 fHistoElectrons->Reset();
249 fHistoMuons->Reset();
250 fHistoPions->Reset();
251 fHistoKaons->Reset();
252 fHistoProtons->Reset();
254 // Getting current tracks, filling histograms
255 for (int n = 0; n < fEsd->GetNumberOfTracks(); ++n) {
257 AliESDtrack *track = fEsd->GetTrack(n);
258 const Double_t sign = fEsd->GetTrack(n)->GetSign();
259 const Int_t prob = GetParticleType(track);
262 fHistopos->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
265 fHistoneg->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
268 fHistoElectrons->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
271 fHistoMuons->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
274 fHistoPions->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
277 fHistoKaons->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
280 fHistoProtons->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
283 fData->DataChanged();
290 //______________________________________________________________________________
291 TEveCaloDataHist* AliEveLego::LoadAllData()
293 // load data from all events ESD
294 fHistoposAllEvents->Reset();
295 fHistonegAllEvents->Reset();
296 fHistoElectronsAllEvents->Reset();
297 fHistoMuonsAllEvents->Reset();
298 fHistoPionsAllEvents->Reset();
299 fHistoKaonsAllEvents->Reset();
300 fHistoProtonsAllEvents->Reset();
302 TTree* t = AliEveEventManager::GetMaster()->GetESDTree();
304 // Getting current tracks for each event, filling histograms
305 Int_t fAcceptedEvents = 0;
306 for (int event = 0; event < t->GetEntries(); event++) {
309 if (fCollisionCandidatesOnly == kTRUE)
310 if (fPhysicsSelection->IsCollisionCandidate(fEsd) == kFALSE) continue;
314 for (int n = 0; n < fEsd->GetNumberOfTracks(); ++n) {
316 AliESDtrack *track = fEsd->GetTrack(n);
317 const Double_t sign = track->GetSign();
318 const Int_t prob = GetParticleType(track);
321 fHistoposAllEvents->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
324 fHistonegAllEvents->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
327 fHistoElectronsAllEvents->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
330 fHistoMuonsAllEvents->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
333 fHistoPionsAllEvents->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
336 fHistoKaonsAllEvents->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
339 fHistoProtonsAllEvents->Fill(track->Eta(), getphi(track->Phi()), fabs(track->Pt()));
343 printf("Number of events loaded: %i, with AliPhysicsSelection: %i\n",fAcceptedEvents,fCollisionCandidatesOnly);
345 fDataAllEvents->DataChanged();
347 return fDataAllEvents;
350 //______________________________________________________________________________
351 TEveCaloDataHist* AliEveLego::FilterData()
354 if ( fTracksId == 2 )
358 fHistoElectrons->Reset();
359 fHistoMuons->Reset();
360 fHistoPions->Reset();
361 fHistoKaons->Reset();
362 fHistoProtons->Reset();
364 const AliESDVertex *pv = fEsd->GetPrimaryVertex();
366 for (Int_t n = 0; n < pv->GetNIndices(); n++ )
368 AliESDtrack *at = fEsd->GetTrack(pv->GetIndices()[n]);
369 const Double_t sign = at->GetSign();
370 const Int_t prob = GetParticleType(at);
373 fHistopos->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
376 fHistoneg->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
379 fHistoElectrons->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
382 fHistoMuons->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
385 fHistoPions->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
388 fHistoKaons->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
391 fHistoProtons->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
395 fData->DataChanged();
398 if (GetPtMax() >= fMaxPt){
399 for (Int_t binx = 1; binx <= 100; binx++) {
400 for (Int_t biny = 1; biny <= 80; biny++) {
402 if (fHistopos->GetBinContent(binx, biny) >= fMaxPt)
403 fHistopos->SetBinContent(binx, biny, fMaxPt);
405 if (fHistoneg->GetBinContent(binx, biny) >= fMaxPt)
406 fHistoneg->SetBinContent(binx, biny, fMaxPt);
408 if (fHistoElectrons->GetBinContent(binx, biny) >= fMaxPt)
409 fHistoElectrons->SetBinContent(binx, biny, fMaxPt);
411 if (fHistoMuons->GetBinContent(binx, biny) >= fMaxPt)
412 fHistoMuons->SetBinContent(binx, biny, fMaxPt);
414 if (fHistoPions->GetBinContent(binx, biny) >= fMaxPt)
415 fHistoPions->SetBinContent(binx, biny, fMaxPt);
417 if (fHistoKaons->GetBinContent(binx, biny) >= fMaxPt)
418 fHistoKaons->SetBinContent(binx, biny, fMaxPt);
420 if (fHistoProtons->GetBinContent(binx, biny) >= fMaxPt)
421 fHistoProtons->SetBinContent(binx, biny, fMaxPt);
426 // Particle type filter
427 if ( fParticleTypeId[0] == kFALSE) fHistopos->Reset();
428 if ( fParticleTypeId[1] == kFALSE) fHistoneg->Reset();
429 if ( fParticleTypeId[2] == kFALSE) fHistoElectrons->Reset();
430 if ( fParticleTypeId[3] == kFALSE) fHistoMuons->Reset();
431 if ( fParticleTypeId[4] == kFALSE) fHistoPions->Reset();
432 if ( fParticleTypeId[5] == kFALSE) fHistoKaons->Reset();
433 if ( fParticleTypeId[6] == kFALSE) fHistoProtons->Reset();
435 fData->DataChanged();
440 //______________________________________________________________________________
441 TEveCaloDataHist* AliEveLego::FilterAllData()
444 if ( fTracksIdAE == 2 )
446 fHistoposAllEvents->Reset();
447 fHistonegAllEvents->Reset();
448 fHistoElectronsAllEvents->Reset();
449 fHistoMuonsAllEvents->Reset();
450 fHistoPionsAllEvents->Reset();
451 fHistoKaonsAllEvents->Reset();
452 fHistoProtonsAllEvents->Reset();
454 TTree* t = AliEveEventManager::GetMaster()->GetESDTree();
456 // Getting current tracks for each event, filling histograms
457 Int_t fAcceptedEvents = 0;
458 for (int event = 0; event < t->GetEntries(); event++) {
461 const AliESDVertex *pv = fEsd->GetPrimaryVertex();
463 if (fCollisionCandidatesOnly == kTRUE)
464 if (fPhysicsSelection->IsCollisionCandidate(fEsd) == kFALSE) continue;
468 for (Int_t n = 0; n < pv->GetNIndices(); n++ )
470 AliESDtrack *at = fEsd->GetTrack(pv->GetIndices()[n]);
471 const Double_t sign = at->GetSign();
472 const Int_t prob = GetParticleType(at);
475 fHistoposAllEvents->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
478 fHistonegAllEvents->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
481 fHistoElectronsAllEvents->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
484 fHistoMuonsAllEvents->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
487 fHistoPionsAllEvents->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
490 fHistoKaonsAllEvents->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
493 fHistoProtonsAllEvents->Fill(at->Eta(), getphi(at->Phi()), fabs(at->Pt()));
497 printf("Number of events loaded: %i, with AliPhysicsSelection: %i\n",fAcceptedEvents,fCollisionCandidatesOnly);
505 fDataAllEvents->DataChanged();
508 if (GetPtMaxAE() >= fMaxPtAE){
509 for (Int_t binx = 1; binx <= 100; binx++) {
510 for (Int_t biny = 1; biny <= 80; biny++) {
512 if (fHistoposAllEvents->GetBinContent(binx, biny) >= fMaxPtAE)
513 fHistoposAllEvents->SetBinContent(binx, biny, fMaxPtAE);
515 if (fHistonegAllEvents->GetBinContent(binx, biny) >= fMaxPtAE)
516 fHistonegAllEvents->SetBinContent(binx, biny, fMaxPtAE);
518 if (fHistoElectronsAllEvents->GetBinContent(binx, biny) >= fMaxPt)
519 fHistoElectronsAllEvents->SetBinContent(binx, biny, fMaxPt);
521 if (fHistoMuonsAllEvents->GetBinContent(binx, biny) >= fMaxPt)
522 fHistoMuonsAllEvents->SetBinContent(binx, biny, fMaxPt);
524 if (fHistoPionsAllEvents->GetBinContent(binx, biny) >= fMaxPt)
525 fHistoPionsAllEvents->SetBinContent(binx, biny, fMaxPt);
527 if (fHistoKaonsAllEvents->GetBinContent(binx, biny) >= fMaxPt)
528 fHistoKaonsAllEvents->SetBinContent(binx, biny, fMaxPt);
530 if (fHistoProtonsAllEvents->GetBinContent(binx, biny) >= fMaxPt)
531 fHistoProtonsAllEvents->SetBinContent(binx, biny, fMaxPt);
537 if ( fParticleTypeIdAE[0] == kFALSE ) fHistoposAllEvents->Reset();
538 if ( fParticleTypeIdAE[1] == kFALSE ) fHistonegAllEvents->Reset();
539 if ( fParticleTypeIdAE[2] == kFALSE ) fHistoElectronsAllEvents->Reset();
540 if ( fParticleTypeIdAE[3] == kFALSE ) fHistoMuonsAllEvents->Reset();
541 if ( fParticleTypeIdAE[4] == kFALSE ) fHistoPionsAllEvents->Reset();
542 if ( fParticleTypeIdAE[5] == kFALSE ) fHistoKaonsAllEvents->Reset();
543 if ( fParticleTypeIdAE[6] == kFALSE ) fHistoProtonsAllEvents->Reset();
545 fDataAllEvents->DataChanged();
547 gEve->Redraw3D(kTRUE);
549 return fDataAllEvents;
552 //______________________________________________________________________________
553 void AliEveLego::Update()
558 // Create new histo2d
564 // Update the viewers
565 gEve->Redraw3D(kTRUE);
568 //______________________________________________________________________________
569 TEveCaloLego* AliEveLego::CreateHistoLego()
571 // Viewer initialization, tab creation
572 if (fHisto2dv == 0) {
573 TEveWindowSlot *fslot = 0;
574 TEveBrowser *fbrowser = gEve->GetBrowser();
576 fslot = TEveWindow::CreateWindowInTab(fbrowser->GetTabRight());
577 fslot->MakeCurrent();
578 fHisto2dv = gEve->SpawnNewViewer("2D Lego Histogram", "2D Lego Histogram");
579 fHisto2ds = gEve->SpawnNewScene("2D Lego Histogram", "2D Lego Histogram");
580 fHisto2dv->AddScene(fHisto2ds);
581 fHisto2dv->SetElementName("2D Lego Viewer");
582 fHisto2ds->SetElementName("2D Lego Scene");
584 fGlv = fHisto2dv->GetGLViewer();
585 fHisto2dLegoOverlay = new TEveCaloLegoOverlay();
586 fGlv->AddOverlayElement(fHisto2dLegoOverlay);
587 fGlv->SetCurrentCamera(TGLViewer::kCameraPerspXOY);
589 fHisto2ds->AddElement(fLego);
591 // move to real world coordinates
592 fLego->InitMainTrans();
593 Float_t sc = TMath::Min(fLego->GetEtaRng(), fLego->GetPhiRng());
594 fLego->RefMainTrans().SetScale(sc, sc, sc);
596 // set event handler to move from perspective to orthographic view.
597 fGlv->SetEventHandler(new TEveLegoEventHandler(fGlv->GetGLWidget(), fGlv, fLego));
599 fHisto2dLegoOverlay->SetCaloLego(fLego);
605 //______________________________________________________________________________
606 TEveCaloLego* AliEveLego::CreateHistoLego(TEveWindowSlot *slot)
608 // Viewer initialization, tab creation
609 if (fHisto2dAllEventsv0 == 0) {
612 fHisto2dAllEventsv0 = gEve->SpawnNewViewer("2D Lego Histogram", "2D Lego Histogram");
613 fHisto2dAllEventss0 = gEve->SpawnNewScene("2D Lego Histogram", "2D Lego Histogram");
614 fHisto2dAllEventsv0->AddScene(fHisto2dAllEventss0);
615 fHisto2dAllEventsv0->SetElementName("2D Lego Viewer");
616 fHisto2dAllEventss0->SetElementName("2D Lego Scene");
618 TGLViewer* glv = fHisto2dAllEventsv0->GetGLViewer();
619 fHisto2dAllEventsLegoOverlay = new TEveCaloLegoOverlay();
620 glv->AddOverlayElement(fHisto2dAllEventsLegoOverlay);
621 glv->SetCurrentCamera(TGLViewer::kCameraPerspXOY);
623 // Plotting histogram lego
624 fLegoAllEvents = new TEveCaloLego(fDataAllEvents);
625 fHisto2dAllEventss0->AddElement(fLegoAllEvents);
627 // Move to real world coordinates
628 fLegoAllEvents->InitMainTrans();
629 Float_t sc = TMath::Min(fLegoAllEvents->GetEtaRng(), fLegoAllEvents->GetPhiRng());
630 fLegoAllEvents->RefMainTrans().SetScale(sc, sc, sc);
632 // Set event handler to move from perspective to orthographic view.
633 glv->SetEventHandler(new TEveLegoEventHandler(glv->GetGLWidget(), glv, fLegoAllEvents));
635 fHisto2dAllEventsLegoOverlay->SetCaloLego(fLegoAllEvents);
638 return fLegoAllEvents;
641 //______________________________________________________________________________
642 TEveCalo3D* AliEveLego::Create3DView()
645 if (fHisto2ds2 == 0) {
646 fHisto2ds2 = gEve->SpawnNewScene("3D Histogram", "3D Histogram");
647 gEve->GetDefaultViewer()->AddScene(fHisto2ds2);
648 fHisto2ds2->SetElementName("3D Histogram Scene");
649 fHisto2ds2->AddElement(fCalo3d);
655 //______________________________________________________________________________
656 TEveCalo3D* AliEveLego::Create3DView(TEveWindowSlot *slot)
658 // creates a 3d view for the 3d histogram
659 if ( fHisto2dAllEventsv1 == 0 ) {
662 fHisto2dAllEventsv1 = gEve->SpawnNewViewer("3D Histogram", "3D Histogram");
663 fHisto2dAllEventss1 = gEve->SpawnNewScene("3D Histogram", "3D Histogram");
664 fHisto2dAllEventsv1->AddScene(fHisto2dAllEventss1);
665 fHisto2dAllEventsv1->SetElementName("3D Histogram Viewer");
666 fHisto2dAllEventss1->SetElementName("3D Histogram Scene");
668 fCalo3dAllEvents = new TEveCalo3D(fDataAllEvents);
670 fCalo3dAllEvents->SetBarrelRadius(550);
671 fCalo3dAllEvents->SetEndCapPos(550);
672 fHisto2dAllEventss1->AddElement(fCalo3dAllEvents);
675 return fCalo3dAllEvents;
678 //______________________________________________________________________________
679 void AliEveLego::CreateProjections(TEveWindowSlot* slot1, TEveWindowSlot* slot2)
681 // create projections
682 if (fHisto2dAllEventsv2 == 0) {
684 slot1->MakeCurrent();
685 fHisto2dAllEventsv2 = gEve->SpawnNewViewer("RPhi projection", "RPhi projection");
686 fHisto2dAllEventss2 = gEve->SpawnNewScene("RPhi projection", "RPhi projection");
687 fHisto2dAllEventsv2->AddScene(fHisto2dAllEventss2);
688 fHisto2dAllEventsv2->SetElementName("RPhi Projection Viewer");
689 fHisto2dAllEventss2->SetElementName("RPhi Projection Scene");
691 TEveProjectionManager* mng1 = new TEveProjectionManager();
692 mng1->SetProjection(TEveProjection::kPT_RPhi);
694 TEveProjectionAxes* axeghisto2dAllEventss1 = new TEveProjectionAxes(mng1);
695 fHisto2dAllEventss2->AddElement(axeghisto2dAllEventss1);
696 TEveCalo2D* fcalo2d1 = (TEveCalo2D*) mng1->ImportElements(fCalo3dAllEvents);
697 fHisto2dAllEventss2->AddElement(fcalo2d1);
699 fHisto2dAllEventsv2->GetGLViewer()->SetCurrentCamera(TGLViewer::kCameraOrthoXOY);
702 if (fHisto2dAllEventsv3 == 0) {
704 slot2->MakeCurrent();
705 fHisto2dAllEventsv3 = gEve->SpawnNewViewer("RhoZ projection", "RhoZ projection");
706 fHisto2dAllEventss3 = gEve->SpawnNewScene("RhoZ projection", "RhoZ projection");
707 fHisto2dAllEventsv3->AddScene(fHisto2dAllEventss3);
708 fHisto2dAllEventsv3->SetElementName("RhoZ Projection Viewer");
709 fHisto2dAllEventss3->SetElementName("RhoZ Projection Viewer");
711 TEveProjectionManager* mng2 = new TEveProjectionManager();
712 mng2->SetProjection(TEveProjection::kPT_RhoZ);
714 TEveProjectionAxes* axeghisto2dAllEventss2 = new TEveProjectionAxes(mng2);
715 fHisto2dAllEventss3->AddElement(axeghisto2dAllEventss2);
716 TEveCalo2D* fcalo2d2 = (TEveCalo2D*) mng2->ImportElements(fCalo3dAllEvents);
717 fHisto2dAllEventss3->AddElement(fcalo2d2);
719 fHisto2dAllEventsv3->GetGLViewer()->SetCurrentCamera(TGLViewer::kCameraOrthoXOY);
725 //______________________________________________________________________________
726 TEveCaloDataHist* AliEveLego::LoadAllEvents()
728 // load all events data from ESD
729 if ( fHisto2dAllEventsSlot == 0 ) {
731 printf("Filling histogram...\n");
735 // Creating 2D histograms
736 fHistoposAllEvents = new TH2F("fHistoposAllEvents","Histo 2d positive",
737 100,-1.5,1.5,80,-kPi,kPi);
738 fHistonegAllEvents = new TH2F("fHistonegAllEvents","Histo 2d negative",
739 100,-1.5,1.5,80,-kPi,kPi);
740 fHistoElectronsAllEvents = new TH2F("fHistoElectronsAllEvents","Histo 2d electrons",
741 100,-1.5,1.5,80,-kPi,kPi);
742 fHistoMuonsAllEvents = new TH2F("fHistoMuonsAllEvents","Histo 2d muons",
743 100,-1.5,1.5,80,-kPi,kPi);
744 fHistoPionsAllEvents = new TH2F("fHistoPionsAllEvents","Histo 2d pions",
745 100,-1.5,1.5,80,-kPi,kPi);
746 fHistoKaonsAllEvents = new TH2F("fHistoKaonsAllEvents","Histo 2d kaons",
747 100,-1.5,1.5,80,-kPi,kPi);
748 fHistoProtonsAllEvents = new TH2F("fHistoProtonsAllEvents","Histo 2d protons",
749 100,-1.5,1.5,80,-kPi,kPi);
751 fHistoposAllEvents->SetDirectory(0);
752 fHistonegAllEvents->SetDirectory(0);
753 fHistoElectronsAllEvents->SetDirectory(0);
754 fHistoMuonsAllEvents->SetDirectory(0);
755 fHistoPionsAllEvents->SetDirectory(0);
756 fHistoKaonsAllEvents->SetDirectory(0);
757 fHistoProtonsAllEvents->SetDirectory(0);
759 // colors from get_pdg_color() in /alice-macros/kine_tracks.C
760 static Color_t fElectro = 5;
761 static Color_t fMuon = 30;
762 static Color_t fPion = 3;
763 static Color_t fKaon = 38;
764 static Color_t fProton = 10;
766 fDataAllEvents = new TEveCaloDataHist();
767 fDataAllEvents->AddHistogram(fHistonegAllEvents);
768 fDataAllEvents->RefSliceInfo(0).Setup("NegCg:", 0, kBlue);
769 fDataAllEvents->AddHistogram(fHistoposAllEvents);
770 fDataAllEvents->RefSliceInfo(1).Setup("PosCg:", 0, kRed);
771 fDataAllEvents->AddHistogram(fHistoElectronsAllEvents);
772 fDataAllEvents->RefSliceInfo(2).Setup("Electrons:", 0, fElectro);
773 fDataAllEvents->AddHistogram(fHistoMuonsAllEvents);
774 fDataAllEvents->RefSliceInfo(3).Setup("Muons:", 0, fMuon);
775 fDataAllEvents->AddHistogram(fHistoPionsAllEvents);
776 fDataAllEvents->RefSliceInfo(4).Setup("Pions:", 0, fPion);
777 fDataAllEvents->AddHistogram(fHistoKaonsAllEvents);
778 fDataAllEvents->RefSliceInfo(5).Setup("Kaons:", 0, fKaon);
779 fDataAllEvents->AddHistogram(fHistoProtonsAllEvents);
780 fDataAllEvents->RefSliceInfo(6).Setup("Protons:", 0, fProton);
782 fDataAllEvents->GetEtaBins()->SetTitleFont(120);
783 fDataAllEvents->GetEtaBins()->SetTitle("h");
784 fDataAllEvents->GetPhiBins()->SetTitleFont(120);
785 fDataAllEvents->GetPhiBins()->SetTitle("f");
786 fDataAllEvents->IncDenyDestroy();
789 fHisto2dAllEventsSlot = TEveWindow::CreateWindowInTab(gEve->GetBrowser()->GetTabRight());
790 TEveWindowPack* packH = fHisto2dAllEventsSlot->MakePack();
791 packH->SetElementName("Projections");
792 packH->SetHorizontal();
793 packH->SetShowTitleBar(kFALSE);
795 fHisto2dAllEventsSlot = packH->NewSlot();
796 TEveWindowPack* pack0 = fHisto2dAllEventsSlot->MakePack();
797 pack0->SetShowTitleBar(kFALSE);
798 TEveWindowSlot* slotLeftTop = pack0->NewSlot();
799 TEveWindowSlot* slotLeftBottom = pack0->NewSlot();
801 fHisto2dAllEventsSlot = packH->NewSlot();
802 TEveWindowPack* pack1 = fHisto2dAllEventsSlot->MakePack();
803 pack1->SetShowTitleBar(kFALSE);
804 TEveWindowSlot* slotRightTop = pack1->NewSlot();
805 TEveWindowSlot* slotRightBottom = pack1->NewSlot();
807 // Creating viewers and scenes
808 Create3DView(slotLeftTop);
809 CreateHistoLego(slotLeftBottom);
810 CreateProjections(slotRightTop, slotRightBottom);
814 gEve->Redraw3D(kTRUE);
816 printf("Filling histogram... Finished\n");
822 return fDataAllEvents;
825 //______________________________________________________________________________
826 void AliEveLego::ApplyParticleTypeSelectionAE()
828 // reload all events applying particle type selection
832 //______________________________________________________________________________
833 Float_t AliEveLego::GetPtMax()
836 return fData->GetMaxVal(fLego->GetPlotEt());
839 //______________________________________________________________________________
840 Float_t AliEveLego::GetPtMaxAE()
842 // return pT max from all events
843 return fDataAllEvents->GetMaxVal(fLegoAllEvents->GetPlotEt());
846 //______________________________________________________________________________
847 Int_t AliEveLego::GetParticleType(AliESDtrack *track)
849 // determine the particle type
850 Double_t *prob = new Double_t[5];
851 track->GetESDpid(prob);
853 Double_t max = prob[0];
856 for (Int_t i = 1 ; i < 5; i++)
866 //______________________________________________________________________________
867 void AliEveLego::SetParticleType(Int_t id)
869 // activate/deactivate particles types
870 if (fParticleTypeId[id] == 0)
872 fParticleTypeId[id] = 1;
874 fParticleTypeId[id] = 0;
880 //______________________________________________________________________________
881 void AliEveLego::SetParticleTypeAE(Int_t id)
883 // activate/deactivate particles types
884 if (fParticleTypeIdAE[id] == 0)
886 fParticleTypeIdAE[id] = 1;
888 fParticleTypeIdAE[id] = 0;
892 //______________________________________________________________________________
893 void AliEveLego::SetMaxPt(Double_t val)
900 //______________________________________________________________________________
901 void AliEveLego::SetMaxPtAE(Double_t val)
908 //______________________________________________________________________________
909 void AliEveLego::SetThreshold(Double_t val)
911 // Setting up the new threshold for all histograms
912 fData->SetSliceThreshold(0,val);
913 fData->SetSliceThreshold(1,val);
914 fData->SetSliceThreshold(2,val);
915 fData->SetSliceThreshold(3,val);
916 fData->SetSliceThreshold(4,val);
917 fData->SetSliceThreshold(5,val);
918 fData->SetSliceThreshold(6,val);
919 fData->DataChanged();
921 gEve->Redraw3D(kTRUE);
924 //______________________________________________________________________________
925 void AliEveLego::SetThresholdAE(Double_t val)
927 // Setting up the new threshold for all histograms
928 fDataAllEvents->SetSliceThreshold(0,val);
929 fDataAllEvents->SetSliceThreshold(1,val);
930 fDataAllEvents->SetSliceThreshold(2,val);
931 fDataAllEvents->SetSliceThreshold(3,val);
932 fDataAllEvents->SetSliceThreshold(4,val);
933 fDataAllEvents->SetSliceThreshold(5,val);
934 fDataAllEvents->SetSliceThreshold(6,val);
935 fDataAllEvents->DataChanged();
937 gEve->Redraw3D(kTRUE);
940 //______________________________________________________________________________
941 void AliEveLego::SwitchDataType()
943 // activate/deactivate MC / real data type
950 //Removing defaul physics selection
952 delete fPhysicsSelection;
953 fPhysicsSelection = NULL;
955 fPhysicsSelection = new AliPhysicsSelection();
956 fPhysicsSelection->SetAnalyzeMC(fIsMC);
957 fPhysicsSelection->Initialize(fEsd->GetRunNumber());
961 //______________________________________________________________________________
962 void AliEveLego::SetCollisionCandidatesOnly()
964 // activate/deactivate MC / real data type
965 if (fCollisionCandidatesOnly == 0)
967 fCollisionCandidatesOnly = 1;
969 fCollisionCandidatesOnly = 0;
974 /******************************************************************************/