]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVE/alice-macros/histo2d.C
removing obsolete 2007-beamtest macros
[u/mrichter/AliRoot.git] / EVE / alice-macros / histo2d.C
CommitLineData
e516f157 1/**************************************************************************
2 * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
3 * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for *
4 * full copyright notice. *
5 **************************************************************************/
6
7double pi = TMath::Pi();
bfd5d6c1 8TEveViewer *g_histo2d_v = 0;
9TEveScene *g_histo2d_s = 0;
10TEveScene *g_histo2d_s2 = 0;
11TEveCaloLegoOverlay* g_histo2d_lego_overlay = 0;
e516f157 12
13
14TEveCaloDataHist* histo2d()
bfd5d6c1 15{
e516f157 16
17 // Access to esdTree
18 AliESDEvent* esd = AliEveEventManager::AssertESD();
19
20 // Creating 2D histograms
21 TH2F *histopos = new TH2F("histopos","Histo 2d positive",100,-1.5,1.5,80,-pi,pi);
22 TH2F *histoneg = new TH2F("histoneg","Histo 2d negative",100,-1.5,1.5,80,-pi,pi);
23
24 cout<<"Event: "<<AliEveEventManager::GetMaster()->GetEventId()
25 <<", Number of tracks: "<<esd->GetNumberOfTracks()<<endl;
26
27 // Getting current tracks, filling histograms
28 for (int n = 0; n < esd->GetNumberOfTracks(); ++n) {
29
30 if (esd->GetTrack(n)->GetSign() > 0) {
31 histopos->Fill(esd->GetTrack(n)->Eta(),
32 getphi(esd->GetTrack(n)->Phi()),
33 fabs(esd->GetTrack(n)->Pt()));
34 } else {
35 histoneg->Fill(esd->GetTrack(n)->Eta(),
36 getphi(esd->GetTrack(n)->Phi()),
37 fabs(esd->GetTrack(n)->Pt()));
38 }
39 }
40
41 TEveCaloDataHist* data = new TEveCaloDataHist();
42 AliEveEventManager::RegisterTransient(data);
43
44 data->AddHistogram(histoneg);
45 data->RefSliceInfo(0).Setup("NegCg:", 0, kBlue);
46 data->AddHistogram(histopos);
47 data->RefSliceInfo(1).Setup("PosCg:", 0, kRed);
48 data->GetEtaBins()->SetTitleFont(120);
49 data->GetEtaBins()->SetTitle("h");
50 data->GetPhiBins()->SetTitleFont(120);
51 data->GetPhiBins()->SetTitle("f");
52 data->IncDenyDestroy();
53
54 // Plotting the lego histogram in a new tab
55 Create_histo_lego(data);
56
57 // Plotting the 3D histogram
58 TEveCalo3D *calo3d = Create_3D_view(data);
59
60 // Plotting projections RPhi and RhoZ
61 Create_projections(data, calo3d);
62
63 gEve->Redraw3D(kTRUE);
64
65 return data;
66}
67
68//______________________________________________________________________________
69Double_t getphi(Double_t phi)
70{
71 if (phi > pi) {
72 phi -= 2*pi;
73 }
74 return phi;
75}
76
77//______________________________________________________________________________
78TEveCaloLego* Create_histo_lego(TEveCaloData* data){
79
80 TGLViewer* glv;
81
82 // Viewer initialization, tab creation
bfd5d6c1 83 if (g_histo2d_v == 0) {
e516f157 84 TEveWindowSlot *slot = 0;
85 TEveBrowser *browser = gEve->GetBrowser();
86
87 slot = TEveWindow::CreateWindowInTab(browser->GetTabRight());
88 slot->MakeCurrent();
bfd5d6c1 89 g_histo2d_v = gEve->SpawnNewViewer("2D Lego Histogram", "2D Lego Histogram");
90 g_histo2d_s = gEve->SpawnNewScene("2D Lego Histogram", "2D Lego Histogram");
91 g_histo2d_v->AddScene(g_histo2d_s);
92 g_histo2d_v->SetElementName("2D Lego Viewer");
93 g_histo2d_s->SetElementName("2D Lego Scene");
94
95 glv = g_histo2d_v->GetGLViewer();
96 g_histo2d_lego_overlay = new TEveCaloLegoOverlay();
97 glv->AddOverlayElement(g_histo2d_lego_overlay);
e516f157 98 glv->SetCurrentCamera(TGLViewer::kCameraPerspXOY);
99 } else {
bfd5d6c1 100 glv = g_histo2d_v->GetGLViewer();
e516f157 101 }
102
103 //plotting histo
104 TEveCaloLego* lego = new TEveCaloLego(data);
bfd5d6c1 105 g_histo2d_s->AddElement(lego);
e516f157 106 AliEveEventManager::RegisterTransient(lego);
107
108 // move to real world coordinates
109 lego->InitMainTrans();
110 Float_t sc = TMath::Min(lego->GetEtaRng(), lego->GetPhiRng());
111 lego->RefMainTrans().SetScale(sc, sc, sc);
112
113 // set event handler to move from perspective to orthographic view.
114 glv->SetEventHandler(new TEveLegoEventHandler(glv->GetGLWidget(), glv, lego));
115
bfd5d6c1 116 g_histo2d_lego_overlay->SetCaloLego(lego);
e516f157 117
118 return lego;
119}
120
121//______________________________________________________________________________
122TEveCalo3D* Create_3D_view(TEveCaloData* data){
123
124 //initialization
bfd5d6c1 125 if (g_histo2d_s2 == 0) {
126 g_histo2d_s2 = gEve->SpawnNewScene("3D Histogram", "3D Histogram");
127 gEve->GetDefaultViewer()->AddScene(g_histo2d_s2);
128 g_histo2d_s2->SetElementName("3D Histogram Scene");
e516f157 129 }
130
131 TEveCalo3D* calo3d = new TEveCalo3D(data);
132 AliEveEventManager::RegisterTransient(calo3d);
133
134 calo3d->SetBarrelRadius(550);
135 calo3d->SetEndCapPos(550);
bfd5d6c1 136 g_histo2d_s2->AddElement(calo3d);
e516f157 137
138 return calo3d;
139}
140
141//______________________________________________________________________________
142AliEveMultiView* Create_projections(TEveCaloData* data, TEveCalo3D *calo3d){
143
144 AliEveMultiView *al = AliEveMultiView::Instance();
145 al->ImportEventRPhi(calo3d);
146 al->ImportEventRhoZ(calo3d);
147
148 return al;
149}