]> git.uio.no Git - u/mrichter/AliRoot.git/blame - DISPLAY/AliDisplayClusters.cxx
Fix for memory leak in fClusters
[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>
024a7e64 23#include <TFile.h>
7bb7ac14 24#include <TPolyMarker3D.h>
25
024a7e64 26#include "AliClusters.h"
27#include "AliDisplay2.h"
7bb7ac14 28#include "AliDisplayClusters.h"
7bb7ac14 29#include "AliITS.h"
30#include "AliITSLoader.h"
7bb7ac14 31#include "AliITSclusterV2.h"
024a7e64 32#include "AliITSgeom.h"
33#include "AliModuleInfo.h"
34#include "AliRun.h"
35#include "AliRunLoader.h"
7bb7ac14 36#include "AliTPCLoader.h"
024a7e64 37#include "AliTPCParam.h"
7bb7ac14 38#include "AliTPCcluster.h"
39
7bb7ac14 40ClassImp(AliDisplayClusters);
41
42//_____________________________________________________________
43AliDisplayClusters::AliDisplayClusters()
44{
45 //Default constructor
46 fPoints = new TPolyMarker3D[gAliDisplay2->GetNbModules()];
47 fName = new char*[gAliDisplay2->GetNbModules()];
48 fNb=0;
49 for(Int_t i=0;i<gAliDisplay2->GetNbModules();i++){
50 fPoints[i].SetMarkerSize(0.2);
51 fPoints[i].SetMarkerColor(2);
52 fPoints[i].SetMarkerStyle(1);
53 }
54}
55
56//_____________________________________________________________
57AliDisplayClusters::~AliDisplayClusters()
58{
59 // Destructor
60 delete [] fPoints;
61}
62
63//_____________________________________________________________
64Int_t AliDisplayClusters::GetNbClusters()
65{
66 // Returns the number of clusters
67 Int_t r=0;
68 for(Int_t i=0;i<fNb;i++){
69 if(gAliDisplay2->GetModuleInfo()->IsEnabled(fName[i])) r+=fPoints[i].GetN();
70 }
71 return r;
72}
73
74//_____________________________________________________________
75void AliDisplayClusters::LoadClusters(const char *name,Int_t nevent)
76{
77 // Loads ITS and TPC clusters
78 if(strstr(name,"ITS")) LoadITSClusters(nevent);
79 if(strstr(name,"TPC")) LoadTPCClusters(nevent);
80}
81
82//_____________________________________________________________
83void AliDisplayClusters::LoadITSClusters(Int_t nevent)
84{
85 // Loads ITS clusters
86 fName[fNb]=new char[strlen("ITS")];
87 strcpy(fName[fNb],"ITS");
88 AliRunLoader *rl = AliRunLoader::Open("galice.root");
89 if(!rl) {
90 cerr<<"Can't open galice.root";
91 return;
92 }
93 AliITSLoader *itsl = (AliITSLoader*)rl->GetLoader("ITSLoader");
94 AliITS *its = (AliITS*)gAlice->GetModule("ITS");
95
96 rl->GetEvent(nevent);
97 itsl->LoadRecPoints();
98 TTree *cTree=itsl->TreeR();
99 if(!cTree)
100 {
101 cerr<<"Error occured during ITS clusters load";
102 return;
103 }
104
105 AliITSgeom *geom=its->GetITSgeom();
106 Int_t count = 0;
107
108 TClonesArray *clusters=new TClonesArray("AliITSclusterV2",10000);
109 TBranch *branch=cTree->GetBranch("Clusters");
110 branch->SetAddress(&clusters);
111 Int_t nentr=(Int_t)cTree->GetEntries();
112 for (Int_t i=0; i<nentr; i++) {
113 if (!cTree->GetEvent(i)) continue;
114
115 Double_t rot[9];
116 geom->GetRotMatrix(i,rot);
117 Int_t lay,lad,det;
118 geom->GetModuleId(i,lay,lad,det);
119 Float_t tx,ty,tz;
120 geom->GetTrans(lay,lad,det,tx,ty,tz);
121
122 Double_t r=-tx*rot[1]+ty*rot[0];
123 if (lay==1) r=-r;
124 Double_t phi=TMath::ATan2(rot[1],rot[0]);
125 if (lay==1) phi-=3.1415927;
126 Double_t cp=TMath::Cos(phi), sp=TMath::Sin(phi);
127
128 Int_t ncl=clusters->GetEntriesFast();
129 while (ncl--) {
130 AliITSclusterV2 *c=(AliITSclusterV2*)clusters->UncheckedAt(ncl);
131 Double_t g[3];
132 g[0]= r*cp + c->GetY()*sp; //
133 g[1]=-r*sp + c->GetY()*cp; //
134 g[2]=c->GetZ();
135 fPoints[fNb].SetPoint(count,g[0],g[1],g[2]);
136 count++;
137 }
138 }
139 fNb++;
140}
141
142//_____________________________________________________________
143void AliDisplayClusters::LoadTPCClusters(Int_t nevent)
144{
145 // Loads TPC clusters
146 fName[fNb]=new char[strlen("TPC")];
147 strcpy(fName[fNb],"TPC");
148 TFile *file = TFile::Open("galice.root");
149 AliTPCParam *dig=(AliTPCParam *)file->Get("75x40_100x60_150x60");
150 if (!dig) {cerr<<"TPC parameters have not been found !\n";}
151 file->Close();
152
153 AliRunLoader *rl = AliRunLoader::Open("galice.root");
154 if(!rl) {
155 cerr<<"Can't open galice.root";
156 return;
157 }
158 AliTPCLoader *itsl = (AliTPCLoader*)rl->GetLoader("TPCLoader");
159 if(!itsl){
160 cerr<<"Can't find Loader";
161 return;
162 }
163
164 rl->GetEvent(nevent);
165 itsl->LoadRecPoints();
166 TTree *cTree=itsl->TreeR();
167 if(!cTree)
168 {
169 cerr<<"Error during TPC clusters load";
170 return;
171 }
172
173 Int_t count = 0;
174 Float_t noiseth = 10;
175
176 AliClusters *clusters=new AliClusters();
177 clusters->SetClass("AliTPCcluster");
178
179 cTree->SetBranchAddress("Segment",&clusters);
180
181 Int_t nrows=Int_t(cTree->GetEntries());
182 for (Int_t n=0; n<nrows; n++) {
183 cTree->GetEvent(n);
184 Int_t sec,row;
185 dig->AdjustSectorRow(clusters->GetID(),sec,row);
186 TClonesArray &clrow=*clusters->GetArray();
187 Int_t ncl=clrow.GetEntriesFast();
188 while (ncl--) {
189 AliTPCcluster *cl=(AliTPCcluster*)clrow[ncl];
190 Double_t x=dig->GetPadRowRadii(sec,row), y=cl->GetY(), z=cl->GetZ();
191 if (cl->GetQ()<noiseth) continue;
192 Float_t cs, sn, tmp;
193 dig->AdjustCosSin(sec,cs,sn);
194 tmp = x*cs-y*sn; y= x*sn+y*cs; x=tmp;
195 fPoints[fNb].SetPoint(count,x,y,z);
196 count++;
197 }
198 clrow.Clear();
199
200 }
201 delete cTree;
202 delete dig;
203 fNb++;
204
205}
206
207//_____________________________________________________________
208void AliDisplayClusters::Draw()
209{
210 // Draws clusters
211 for(Int_t i=0;i<fNb;i++){
212 if(gAliDisplay2->GetModuleInfo()->IsEnabled(fName[i])) fPoints[i].Draw();
213 }
214}
215