]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVE/alice-macros/vplot_tpc.C
cosmetics
[u/mrichter/AliRoot.git] / EVE / alice-macros / vplot_tpc.C
CommitLineData
08b0f222 1// $Id$
2// Main authors: Adam Jacholkowski & Matevz Tadel: 2009
3
4/**************************************************************************
5 * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
6 * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for *
7 * full copyright notice. *
8 **************************************************************************/
9
176fc00c 10#if !defined(__CINT__) || defined(__MAKECINT__)
08b0f222 11#include <TGLViewer.h>
12#include <TGLCameraOverlay.h>
08b0f222 13#include <TEveManager.h>
14#include <TEveBrowser.h>
15#include <TEveViewer.h>
16#include <TEveScene.h>
17#include <TEvePointSet.h>
18
6c49a8e1 19#include <AliCluster.h>
20#include <AliESDEvent.h>
21#include <AliESDVertex.h>
22#include <AliRunLoader.h>
23#include <AliTPCClustersRow.h>
24#include <AliEveEventManager.h>
08b0f222 25#endif
26
27TEveViewer *gVPTPCView = 0;
28TEveScene *gVPTPCScene = 0;
29
30TEvePointSet* vplot_tpc(TEveElement* cont=0, Float_t maxR=270)
31{
32 if (gVPTPCView == 0)
33 {
34 TEveWindowSlot *slot = 0;
35 TEveBrowser *browser = gEve->GetBrowser();
36
37 slot = TEveWindow::CreateWindowInTab(browser->GetTabRight());
38 slot->MakeCurrent();
39 gVPTPCView = gEve->SpawnNewViewer("V-Plot", "");
40 gVPTPCScene = gEve->SpawnNewScene("V-Plot", "Scene holding elements for the V-Plot TPC.");
41 gVPTPCView->AddScene(gVPTPCScene);
42
43 TGLViewer *glv = gVPTPCView->GetGLViewer();
44 glv->SetCurrentCamera(TGLViewer::kCameraOrthoXOY);
45 glv->ResetCamerasAfterNextUpdate();
46
47 TGLCameraOverlay* co = glv->GetCameraOverlay();
48 co->SetShowOrthographic(true); //(false);
49 co->SetOrthographicMode(TGLCameraOverlay::kAxis); // ::kPlaneIntersect or ::kBar
50 }
51
52 const Int_t kMaxCl=100*160;
53
54 AliEveEventManager::AssertGeometry();
55
b1916ff9 56 Double_t pvert[3] = { 0, 0, 0 };
57 if (AliEveEventManager::HasESD())
58 {
59 AliESDEvent *esd = AliEveEventManager::AssertESD();
60 const AliESDVertex *tpcv = esd->GetPrimaryVertexTPC();
61 if (tpcv->GetStatus())
62 tpcv->GetXYZ(pvert);
63 else
64 Info("vplot_tpc", "Primary vertex TPC not available, using 0.");
65 }
66 else
67 {
68 Info("vplot_tpc", "ESD not available, using 0 for primary vertex.");
69 }
70
08b0f222 71 AliRunLoader* rl = AliEveEventManager::AssertRunLoader();
72 rl->LoadRecPoints("TPC");
73
74 TTree *cTree = rl->GetTreeR("TPC", false);
75 if (cTree == 0)
76 return 0;
77
78 AliTPCClustersRow *clrow = new AliTPCClustersRow();
79 clrow->SetClass("AliTPCclusterMI");
80 clrow->SetArray(kMaxCl);
81 cTree->SetBranchAddress("Segment", &clrow);
82
83 TEvePointSet* vplot = new TEvePointSet(kMaxCl);
84 vplot->SetOwnIds(kTRUE);
85
86 const Float_t phimin = -3.15;
87 const Float_t phimax = 3.15;
88 const Float_t etamin = -1.2;
89 const Float_t etamax = 1.2;
90 const Float_t vconst = 0.0003;
91
92 Float_t rhomax = 246.6;
93 Float_t rholim, rhoup;
94 Float_t zmax = 250.0;
95 Float_t rho, eta, phi, theta;
96 Float_t r3d, r3dmax, r3d1, r3d2;
97 //
98 Float_t maxRsqr = maxR*maxR;
99
100 Int_t nentr = (Int_t) cTree->GetEntries();
101 for (Int_t i = 0; i < nentr; ++i)
102 {
103 if (!cTree->GetEvent(i)) continue;
104
105 TClonesArray *cl = clrow->GetArray();
106 Int_t ncl = cl->GetEntriesFast();
107
108 while (ncl--)
109 {
110 AliCluster *c = (AliCluster*) cl->UncheckedAt(ncl);
111 Float_t g[3]; //global coordinates
112 c->GetGlobalXYZ(g);
b1916ff9 113 g[0] -= pvert[0];
114 g[1] -= pvert[1];
115 g[2] -= pvert[2];
08b0f222 116 if (g[0]*g[0] + g[1]*g[1] < maxRsqr)
117 {
118 phi = TMath::ATan2(g[1], g[0]);
119 rho = TMath::Sqrt(g[0]*g[0] + g[1]*g[1]);
120 theta = TMath::ATan2(rho, g[2]);
121 eta = -1.0*TMath::Log(TMath::Tan(0.5*theta));
122 rhoup = zmax*rho / TMath::Abs(g[2]);
123 rholim = TMath::Min(rhoup,rhomax);
124 // a version using rather r3d
125 r3d = TMath::Sqrt(rho*rho + g[2]*g[2]);
126 r3d1 = rhomax / TMath::Sin(theta);
127 r3d2 = TMath::Abs(zmax/TMath::Cos(theta));
128 r3dmax = TMath::Min(r3d1, r3d2);
129
130 if (eta>etamin && eta<etamax && phi>phimin && phi<phimax)
131 {
132 Float_t deta = vconst*(r3dmax - r3d);
133
134 vplot->SetNextPoint(eta + deta, phi, 0);
135 vplot->SetNextPoint(eta - deta, phi, 0);
136 }
137 }
138 }
139 cl->Clear();
140 }
141
142 delete clrow;
143
144 rl->UnloadRecPoints("TPC");
145
146 if (vplot->Size() == 0 && gEve->GetKeepEmptyCont() == kFALSE)
147 {
148 Warning("vplot_tpc.C", "No TPC clusters were found.");
149 delete vplot;
150 return 0;
151 }
152
153 vplot->SetName("V=Plot TPC");
154 vplot->SetTitle(Form("N=%d", vplot->Size() / 2));
155
156 vplot->SetMainColor (kOrange);
157 vplot->SetMarkerSize(0.2);
158 vplot->SetMarkerStyle(1);
159 // vplot->ApplyVizTag("V-Plot TPC", "V-Plot");
160
161 if (cont)
162 {
163 cont->AddElement(vplot);
164 }
165 else
166 {
167 gVPTPCScene->AddElement(vplot);
168 }
169 AliEveEventManager::RegisterTransient(vplot);
170
171 gEve->Redraw3D();
172
173 return vplot;
174}