1 // Author: Stefano Carrazza 2010
3 /**************************************************************************
4 * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
5 * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for *
6 * full copyright notice. *
7 **************************************************************************/
9 #if !defined(__CINT__) || defined(__MAKECINT__)
10 #include <TGLViewer.h>
11 #include <TGLWidget.h>
15 #include <TEveBrowser.h>
17 #include <TEveCaloData.h>
18 #include <TEveCaloLegoOverlay.h>
19 #include <TEveLegoEventHandler.h>
20 #include <TEveManager.h>
21 #include <TEveProjectionManager.h>
22 #include <TEveProjectionAxes.h>
23 #include <TEveScene.h>
24 #include <TEveTrans.h>
25 #include <TEveViewer.h>
26 #include <TEveWindow.h>
28 #include <AliESDEvent.h>
29 #include <AliEveEventManager.h>
32 double pi = TMath::Pi();
33 TEveViewer *g_histo2d_all_events_v0 = 0;
34 TEveViewer *g_histo2d_all_events_v1 = 0;
35 TEveViewer *g_histo2d_all_events_v2 = 0;
36 TEveViewer *g_histo2d_all_events_v3 = 0;
37 TEveScene *g_histo2d_all_events_s0 = 0;
38 TEveScene *g_histo2d_all_events_s1 = 0;
39 TEveScene *g_histo2d_all_events_s2 = 0;
40 TEveScene *g_histo2d_all_events_s3 = 0;
41 TEveCaloLegoOverlay* g_histo2d_all_events_lego_overlay = 0;
42 TEveWindowSlot* g_histo2d_all_events_slot = 0;
44 Double_t GetPhi(Double_t phi);
45 TEveCaloLego* CreateHistoLego(TEveCaloData* data, TEveWindowSlot* slot);
46 TEveCalo3D* Create3DView(TEveCaloData* data, TEveWindowSlot* slot);
47 void CreateProjections(TEveCaloData* data, TEveCalo3D *calo3d, TEveWindowSlot* slot1, TEveWindowSlot* slot2);
49 TEveCaloDataHist* histo2d_all_events()
52 TEveCaloDataHist* data_t;
54 if ( g_histo2d_all_events_slot == 0 ) {
55 Info("histo2d_all_events", "Filling histogram...");
58 AliESDEvent* esd = AliEveEventManager::AssertESD();
59 TTree* t = AliEveEventManager::GetMaster()->GetESDTree();
61 // Creating 2D histograms
62 TH2F *histopos_t = new TH2F("histopos_t","Histo 2d positive",
63 100,-1.5,1.5,80,-pi,pi);
64 TH2F *histoneg_t = new TH2F("histoneg_t","Histo 2d negative",
65 100,-1.5,1.5,80,-pi,pi);
67 // Getting current tracks for each event, filling histograms
68 for ( int event = 0; event < t->GetEntries(); event++ ) {
70 for ( int n = 0; n < esd->GetNumberOfTracks(); ++n ) {
72 if ( esd->GetTrack(n)->GetSign() > 0 ) {
73 histopos_t->Fill(esd->GetTrack(n)->Eta(),
74 GetPhi(esd->GetTrack(n)->Phi()),
75 fabs(esd->GetTrack(n)->Pt()));
77 histoneg_t->Fill(esd->GetTrack(n)->Eta(),
78 GetPhi(esd->GetTrack(n)->Phi()),
79 fabs(esd->GetTrack(n)->Pt()));
84 data_t = new TEveCaloDataHist();
85 data_t->AddHistogram(histoneg_t);
86 data_t->RefSliceInfo(0).Setup("NegCg:", 0, kBlue);
87 data_t->AddHistogram(histopos_t);
88 data_t->RefSliceInfo(1).Setup("PosCg:", 0, kRed);
89 data_t->GetEtaBins()->SetTitleFont(120);
90 data_t->GetEtaBins()->SetTitle("h");
91 data_t->GetPhiBins()->SetTitleFont(120);
92 data_t->GetPhiBins()->SetTitle("f");
93 data_t->IncDenyDestroy();
96 g_histo2d_all_events_slot = TEveWindow::CreateWindowInTab(gEve->GetBrowser()->GetTabRight());
97 TEveWindowPack* packH = g_histo2d_all_events_slot->MakePack();
98 packH->SetElementName("Projections");
99 packH->SetHorizontal();
100 packH->SetShowTitleBar(kFALSE);
102 g_histo2d_all_events_slot = packH->NewSlot();
103 TEveWindowPack* pack0 = g_histo2d_all_events_slot->MakePack();
104 pack0->SetShowTitleBar(kFALSE);
105 TEveWindowSlot* slotLeftTop = pack0->NewSlot();
106 TEveWindowSlot* slotLeftBottom = pack0->NewSlot();
108 g_histo2d_all_events_slot = packH->NewSlot();
109 TEveWindowPack* pack1 = g_histo2d_all_events_slot->MakePack();
110 pack1->SetShowTitleBar(kFALSE);
111 TEveWindowSlot* slotRightTop = pack1->NewSlot();
112 TEveWindowSlot* slotRightBottom = pack1->NewSlot();
114 // Creating viewers and scenes
115 TEveCalo3D* calo3d = Create3DView(data_t, slotLeftTop);
116 CreateHistoLego(data_t, slotLeftBottom);
117 CreateProjections(data_t, calo3d, slotRightTop, slotRightBottom);
119 gEve->Redraw3D(kTRUE);
121 Info("histo2d_all_events", "...Finished");
127 //______________________________________________________________________________
128 Double_t GetPhi(Double_t phi)
136 //______________________________________________________________________________
137 TEveCaloLego* CreateHistoLego(TEveCaloData* data, TEveWindowSlot* slot){
141 // Viewer initialization, tab creation
142 if ( g_histo2d_all_events_v0 == 0 ) {
144 TEveBrowser *browser = gEve->GetBrowser();
146 g_histo2d_all_events_v0 = gEve->SpawnNewViewer("2D Lego Histogram", "2D Lego Histogram");
147 g_histo2d_all_events_s0 = gEve->SpawnNewScene("2D Lego Histogram", "2D Lego Histogram");
148 g_histo2d_all_events_v0->AddScene(g_histo2d_all_events_s0);
149 g_histo2d_all_events_v0->SetElementName("2D Lego Viewer");
150 g_histo2d_all_events_s0->SetElementName("2D Lego Scene");
152 TGLViewer* glv = g_histo2d_all_events_v0->GetGLViewer();
153 g_histo2d_all_events_lego_overlay = new TEveCaloLegoOverlay();
154 glv->AddOverlayElement(g_histo2d_all_events_lego_overlay);
155 glv->SetCurrentCamera(TGLViewer::kCameraPerspXOY);
157 // Plotting histogram lego
158 lego = new TEveCaloLego(data);
159 g_histo2d_all_events_s0->AddElement(lego);
161 // Move to real world coordinates
162 lego->InitMainTrans();
163 Float_t sc = TMath::Min(lego->GetEtaRng(), lego->GetPhiRng());
164 lego->RefMainTrans().SetScale(sc, sc, sc);
166 // Set event handler to move from perspective to orthographic view.
167 glv->SetEventHandler(new TEveLegoEventHandler(glv->GetGLWidget(), glv, lego));
169 g_histo2d_all_events_lego_overlay->SetCaloLego(lego);
175 //______________________________________________________________________________
176 TEveCalo3D* Create3DView(TEveCaloData* data, TEveWindowSlot* slot){
180 if ( g_histo2d_all_events_v1 == 0 ) {
182 TEveBrowser *browser = gEve->GetBrowser();
184 g_histo2d_all_events_v1 = gEve->SpawnNewViewer("3D Histogram", "3D Histogram");
185 g_histo2d_all_events_s1 = gEve->SpawnNewScene("3D Histogram", "3D Histogram");
186 g_histo2d_all_events_v1->AddScene(g_histo2d_all_events_s1);
187 g_histo2d_all_events_v1->SetElementName("3D Histogram Viewer");
188 g_histo2d_all_events_s1->SetElementName("3D Histogram Scene");
190 calo3d = new TEveCalo3D(data);
192 calo3d->SetBarrelRadius(550);
193 calo3d->SetEndCapPos(550);
194 g_histo2d_all_events_s1->AddElement(calo3d);
200 //______________________________________________________________________________
201 void CreateProjections(TEveCaloData* data, TEveCalo3D *calo3d, TEveWindowSlot* slot1, TEveWindowSlot* slot2){
203 if ( g_histo2d_all_events_v2 == 0 ) {
205 TEveBrowser *browser = gEve->GetBrowser();
206 slot1->MakeCurrent();
207 g_histo2d_all_events_v2 = gEve->SpawnNewViewer("RPhi projection", "RPhi projection");
208 g_histo2d_all_events_s2 = gEve->SpawnNewScene("RPhi projection", "RPhi projection");
209 g_histo2d_all_events_v2->AddScene(g_histo2d_all_events_s2);
210 g_histo2d_all_events_v2->SetElementName("RPhi Projection Viewer");
211 g_histo2d_all_events_s2->SetElementName("RPhi Projection Scene");
213 TEveProjectionManager* mng1 = new TEveProjectionManager();
214 mng1->SetProjection(TEveProjection::kPT_RPhi);
216 TEveProjectionAxes* axeg_histo2d_all_events_s1 = new TEveProjectionAxes(mng1);
217 g_histo2d_all_events_s2->AddElement(axeg_histo2d_all_events_s1);
218 TEveCalo2D* calo2d1 = (TEveCalo2D*) mng1->ImportElements(calo3d);
219 g_histo2d_all_events_s2->AddElement(calo2d1);
221 g_histo2d_all_events_v2->GetGLViewer()->SetCurrentCamera(TGLViewer::kCameraOrthoXOY);
224 if ( g_histo2d_all_events_v3 == 0 ) {
226 TEveBrowser *browser = gEve->GetBrowser();
227 slot2->MakeCurrent();
228 g_histo2d_all_events_v3 = gEve->SpawnNewViewer("RhoZ projection", "RhoZ projection");
229 g_histo2d_all_events_s3 = gEve->SpawnNewScene("RhoZ projection", "RhoZ projection");
230 g_histo2d_all_events_v3->AddScene(g_histo2d_all_events_s3);
231 g_histo2d_all_events_v3->SetElementName("RhoZ Projection Viewer");
232 g_histo2d_all_events_s3->SetElementName("RhoZ Projection Viewer");
234 TEveProjectionManager* mng2 = new TEveProjectionManager();
235 mng2->SetProjection(TEveProjection::kPT_RhoZ);
237 TEveProjectionAxes* axeg_histo2d_all_events_s2 = new TEveProjectionAxes(mng2);
238 g_histo2d_all_events_s3->AddElement(axeg_histo2d_all_events_s2);
239 TEveCalo2D* calo2d2 = (TEveCalo2D*) mng2->ImportElements(calo3d);
240 g_histo2d_all_events_s3->AddElement(calo2d2);
242 g_histo2d_all_events_v3->GetGLViewer()->SetCurrentCamera(TGLViewer::kCameraOrthoXOY);