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