]> git.uio.no Git - u/mrichter/AliRoot.git/blame - DISPLAY/AliDisplayClusters.cxx
OpengAliceFile method removed, loaders data members set in the constructor, loading...
[u/mrichter/AliRoot.git] / DISPLAY / AliDisplayClusters.cxx
CommitLineData
7bb7ac14 1/**************************************************************************
2 * Copyright(c) 1998-2003, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15/////////////////////////////////////////////////////////////////////////
16// ALICE DISPLAY CLUSTERS CLASS //
17// Author: Mayeul ROUSSELET //
18// e-mail: Mayeul.Rousselet@cern.ch //
19// Last update:26/08/2003 //
20/////////////////////////////////////////////////////////////////////////
21
22#include <Riostream.h>
23#include <TPolyMarker3D.h>
24
25#include "AliDisplayClusters.h"
26#include "AliModuleInfo.h"
27
28#include "AliRun.h"
29#include "AliRunLoader.h"
30#include "AliITS.h"
31#include "AliITSLoader.h"
32#include "AliITSgeom.h"
33#include "AliITSclusterV2.h"
34#include "AliTPCParam.h"
35#include "AliTPCLoader.h"
36#include "AliClusters.h"
37#include "AliTPCcluster.h"
38
39#include "AliDisplay2.h"
40
41ClassImp(AliDisplayClusters);
42
43//_____________________________________________________________
44AliDisplayClusters::AliDisplayClusters()
45{
46 //Default constructor
47 fPoints = new TPolyMarker3D[gAliDisplay2->GetNbModules()];
48 fName = new char*[gAliDisplay2->GetNbModules()];
49 fNb=0;
50 for(Int_t i=0;i<gAliDisplay2->GetNbModules();i++){
51 fPoints[i].SetMarkerSize(0.2);
52 fPoints[i].SetMarkerColor(2);
53 fPoints[i].SetMarkerStyle(1);
54 }
55}
56
57//_____________________________________________________________
58AliDisplayClusters::~AliDisplayClusters()
59{
60 // Destructor
61 delete [] fPoints;
62}
63
64//_____________________________________________________________
65Int_t AliDisplayClusters::GetNbClusters()
66{
67 // Returns the number of clusters
68 Int_t r=0;
69 for(Int_t i=0;i<fNb;i++){
70 if(gAliDisplay2->GetModuleInfo()->IsEnabled(fName[i])) r+=fPoints[i].GetN();
71 }
72 return r;
73}
74
75//_____________________________________________________________
76void AliDisplayClusters::LoadClusters(const char *name,Int_t nevent)
77{
78 // Loads ITS and TPC clusters
79 if(strstr(name,"ITS")) LoadITSClusters(nevent);
80 if(strstr(name,"TPC")) LoadTPCClusters(nevent);
81}
82
83//_____________________________________________________________
84void AliDisplayClusters::LoadITSClusters(Int_t nevent)
85{
86 // Loads ITS clusters
87 fName[fNb]=new char[strlen("ITS")];
88 strcpy(fName[fNb],"ITS");
89 AliRunLoader *rl = AliRunLoader::Open("galice.root");
90 if(!rl) {
91 cerr<<"Can't open galice.root";
92 return;
93 }
94 AliITSLoader *itsl = (AliITSLoader*)rl->GetLoader("ITSLoader");
95 AliITS *its = (AliITS*)gAlice->GetModule("ITS");
96
97 rl->GetEvent(nevent);
98 itsl->LoadRecPoints();
99 TTree *cTree=itsl->TreeR();
100 if(!cTree)
101 {
102 cerr<<"Error occured during ITS clusters load";
103 return;
104 }
105
106 AliITSgeom *geom=its->GetITSgeom();
107 Int_t count = 0;
108
109 TClonesArray *clusters=new TClonesArray("AliITSclusterV2",10000);
110 TBranch *branch=cTree->GetBranch("Clusters");
111 branch->SetAddress(&clusters);
112 Int_t nentr=(Int_t)cTree->GetEntries();
113 for (Int_t i=0; i<nentr; i++) {
114 if (!cTree->GetEvent(i)) continue;
115
116 Double_t rot[9];
117 geom->GetRotMatrix(i,rot);
118 Int_t lay,lad,det;
119 geom->GetModuleId(i,lay,lad,det);
120 Float_t tx,ty,tz;
121 geom->GetTrans(lay,lad,det,tx,ty,tz);
122
123 Double_t r=-tx*rot[1]+ty*rot[0];
124 if (lay==1) r=-r;
125 Double_t phi=TMath::ATan2(rot[1],rot[0]);
126 if (lay==1) phi-=3.1415927;
127 Double_t cp=TMath::Cos(phi), sp=TMath::Sin(phi);
128
129 Int_t ncl=clusters->GetEntriesFast();
130 while (ncl--) {
131 AliITSclusterV2 *c=(AliITSclusterV2*)clusters->UncheckedAt(ncl);
132 Double_t g[3];
133 g[0]= r*cp + c->GetY()*sp; //
134 g[1]=-r*sp + c->GetY()*cp; //
135 g[2]=c->GetZ();
136 fPoints[fNb].SetPoint(count,g[0],g[1],g[2]);
137 count++;
138 }
139 }
140 fNb++;
141}
142
143//_____________________________________________________________
144void AliDisplayClusters::LoadTPCClusters(Int_t nevent)
145{
146 // Loads TPC clusters
147 fName[fNb]=new char[strlen("TPC")];
148 strcpy(fName[fNb],"TPC");
149 TFile *file = TFile::Open("galice.root");
150 AliTPCParam *dig=(AliTPCParam *)file->Get("75x40_100x60_150x60");
151 if (!dig) {cerr<<"TPC parameters have not been found !\n";}
152 file->Close();
153
154 AliRunLoader *rl = AliRunLoader::Open("galice.root");
155 if(!rl) {
156 cerr<<"Can't open galice.root";
157 return;
158 }
159 AliTPCLoader *itsl = (AliTPCLoader*)rl->GetLoader("TPCLoader");
160 if(!itsl){
161 cerr<<"Can't find Loader";
162 return;
163 }
164
165 rl->GetEvent(nevent);
166 itsl->LoadRecPoints();
167 TTree *cTree=itsl->TreeR();
168 if(!cTree)
169 {
170 cerr<<"Error during TPC clusters load";
171 return;
172 }
173
174 Int_t count = 0;
175 Float_t noiseth = 10;
176
177 AliClusters *clusters=new AliClusters();
178 clusters->SetClass("AliTPCcluster");
179
180 cTree->SetBranchAddress("Segment",&clusters);
181
182 Int_t nrows=Int_t(cTree->GetEntries());
183 for (Int_t n=0; n<nrows; n++) {
184 cTree->GetEvent(n);
185 Int_t sec,row;
186 dig->AdjustSectorRow(clusters->GetID(),sec,row);
187 TClonesArray &clrow=*clusters->GetArray();
188 Int_t ncl=clrow.GetEntriesFast();
189 while (ncl--) {
190 AliTPCcluster *cl=(AliTPCcluster*)clrow[ncl];
191 Double_t x=dig->GetPadRowRadii(sec,row), y=cl->GetY(), z=cl->GetZ();
192 if (cl->GetQ()<noiseth) continue;
193 Float_t cs, sn, tmp;
194 dig->AdjustCosSin(sec,cs,sn);
195 tmp = x*cs-y*sn; y= x*sn+y*cs; x=tmp;
196 fPoints[fNb].SetPoint(count,x,y,z);
197 count++;
198 }
199 clrow.Clear();
200
201 }
202 delete cTree;
203 delete dig;
204 fNb++;
205
206}
207
208//_____________________________________________________________
209void AliDisplayClusters::Draw()
210{
211 // Draws clusters
212 for(Int_t i=0;i<fNb;i++){
213 if(gAliDisplay2->GetModuleInfo()->IsEnabled(fName[i])) fPoints[i].Draw();
214 }
215}
216