]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliITSVertexerIons.cxx
AliHMPIDDigitN no longer needed
[u/mrichter/AliRoot.git] / ITS / AliITSVertexerIons.cxx
CommitLineData
7d62fb64 1#include "AliITSDetTypeRec.h"
b5d712e2 2#include "AliITSVertexerIons.h"
3#include "AliITSVertexerPPZ.h"
4#include "AliESDVertex.h"
c5f0f3c1 5#include "AliITSgeom.h"
88cb7938 6#include "AliITSLoader.h"
c5f0f3c1 7#include "AliITSRecPoint.h"
b5d712e2 8#include <TTree.h>
c5f0f3c1 9#include <TH1.h>
10#include <TF1.h>
b5d712e2 11
c5f0f3c1 12//////////////////////////////////////////////////////////////////////
13// AliITSVertexerIons is a class for full 3D primary vertex //
14// finding optimized for Ion-Ion interactions //
15// //
16// //
17// //
18// //
19// Written by Giuseppe Lo Re and Francesco Riggi //
20// Giuseppe.Lore@ct.infn.it //
c5f0f3c1 21// Franco.Riggi@ct.infn.it //
22// //
fe4280bf 23// Release date: Mar 2004 //
c5f0f3c1 24// //
25// //
26//////////////////////////////////////////////////////////////////////
27
c5f0f3c1 28
b5d712e2 29ClassImp(AliITSVertexerIons)
30
fe4280bf 31//______________________________________________________________________
8221b41b 32 AliITSVertexerIons::AliITSVertexerIons():AliITSVertexer(),
33fNpThreshold(0),
34fMaxDeltaPhi(0),
35fMaxDeltaZ(0){
c5f0f3c1 36 // Default Constructor
37
7d62fb64 38 //fITS = 0;
83028144 39 SetNpThreshold();
fe4280bf 40 SetMaxDeltaPhi();
41 SetMaxDeltaZ();
c5f0f3c1 42}
43
44//______________________________________________________________________
8221b41b 45AliITSVertexerIons::AliITSVertexerIons(TString fn):AliITSVertexer(fn),
46fNpThreshold(0),
47fMaxDeltaPhi(0),
48fMaxDeltaZ(0) {
c5f0f3c1 49 // Standard constructor
88cb7938 50
7d62fb64 51 //fITS = 0;
83028144 52 SetNpThreshold();
fe4280bf 53 SetMaxDeltaPhi();
b5d712e2 54 SetMaxDeltaZ();
55}
8221b41b 56/*
b5d712e2 57//______________________________________________________________________
58AliITSVertexerIons::AliITSVertexerIons(const AliITSVertexerIons &source):AliITSVertexer(source) {
59 // Copy constructor
60 // Copies are not allowed. The method is protected to avoid misuse.
61 Error("AliITSVertexerIons","Copy constructor not allowed\n");
62}
8221b41b 63*/
b5d712e2 64//_________________________________________________________________________
8221b41b 65//AliITSVertexerIons& AliITSVertexerIons::operator=(const AliITSVertexerIons &/*source*/) {
b5d712e2 66 // Assignment operator
67 // Assignment is not allowed. The method is protected to avoid misuse.
8221b41b 68 //Error("= operator","Assignment operator not allowed\n");
69 //return *this;
70//}
c5f0f3c1 71
72
73//______________________________________________________________________
74AliITSVertexerIons::~AliITSVertexerIons() {
75 // Default Destructor
7d62fb64 76 //fITS = 0;
c5f0f3c1 77}
78
79//______________________________________________________________________
b5d712e2 80AliESDVertex* AliITSVertexerIons::FindVertexForCurrentEvent(Int_t evnumber){
81// Defines the AliESDVertex for the current event
fe4280bf 82
c5f0f3c1 83 fCurrentVertex = 0;
fe4280bf 84
88cb7938 85 AliRunLoader *rl = AliRunLoader::GetRunLoader();
7d62fb64 86 AliITSLoader* itsloader = (AliITSLoader*)rl->GetLoader("ITSLoader");
87 TDirectory * olddir = gDirectory;
88 rl->CdGAFile();
89 AliITSgeom* g2 = (AliITSgeom*)gDirectory->Get("AliITSgeom");
90 olddir->cd();
91
fe4280bf 92 TTree *tr = itsloader->TreeR();
7d62fb64 93 AliITSDetTypeRec detTypeRec;
94
95 detTypeRec.SetTreeAddressR(tr);
96
97 TClonesArray *recpoints = detTypeRec.RecPoints();
98 AliITSRecPoint *pnt;
99
88cb7938 100
fe4280bf 101 Int_t npoints=0;
102 Int_t nopoints1=40000;
103 Int_t nopoints2=40000;
88cb7938 104 Float_t l[3], p[3];
c5f0f3c1 105
83028144 106 Double_t *z1, *z2, *y1, *y2, *x1, *x2, *phi1, *phi2, *r1, *r2;
107 z1=new Double_t[nopoints1];
108 z2=new Double_t[nopoints2];
109 y1=new Double_t[nopoints1];
110 y2=new Double_t[nopoints2];
111 x1=new Double_t[nopoints1];
112 x2=new Double_t[nopoints2];
113 phi1=new Double_t[nopoints1];
114 phi2=new Double_t[nopoints2];
115 r1=new Double_t[nopoints1];
116 r2=new Double_t[nopoints2];
c5f0f3c1 117
fe4280bf 118 Double_t mxpiu = 0;
119 Double_t mxmeno = 0;
120 Double_t mypiu = 0;
121 Double_t mymeno = 0;
122 Double_t r=0;
123
124 Int_t np1=0, np2=0;
125 for(Int_t i=g2->GetStartSPD();i<=g2->GetLastSPD();i++) {
7d62fb64 126 detTypeRec.ResetRecPoints();
fe4280bf 127 tr->GetEvent(i);
128 npoints = recpoints->GetEntries();
129 for (Int_t ipoint=0;ipoint<npoints;ipoint++) {
130 pnt = (AliITSRecPoint*)recpoints->UncheckedAt(ipoint);
00a7cc50 131 l[0]=pnt->GetDetLocalX();
fe4280bf 132 l[1]=0;
00a7cc50 133 l[2]=pnt->GetDetLocalZ();
fe4280bf 134 g2->LtoG(i, l, p);
135 r=TMath::Sqrt(TMath::Power(p[0],2)+TMath::Power(p[1],2));
136
137 if(i<80 && TMath::Abs(p[2])<14.35) {
138 y1[np1]=p[1];
139 x1[np1]=p[0];
140 z1[np1]=p[2];
141 if(p[0]>0) mxpiu++;
142 if(p[0]<0) mxmeno++;
143 if(p[1]>0) mypiu++;
144 if(p[1]<0) mymeno++;
145 np1++;
2257f27e 146 }
fe4280bf 147 if(i>=80 && TMath::Abs(p[2])<14.35) {
148 y2[np2]=p[1];
149 x2[np2]=p[0];
150 z2[np2]=p[2];
151 np2++;
88cb7938 152 }
153 }
fe4280bf 154 }
c5f0f3c1 155
83028144 156 if(np1<fNpThreshold) {
83028144 157 Warning("FindVertexForCurrentEvent","AliITSVertexerIons finder is not reliable for low multiplicity events. Switching to AliITSVertexerPPZ with default parameters...\n");
158 Warning("FindVertexForCurrentEvent","N rec points = %d - Threshold is %d",np1,fNpThreshold);
fe4280bf 159 AliITSVertexerPPZ *dovert = new AliITSVertexerPPZ("default");
83028144 160 fCurrentVertex =dovert->FindVertexForCurrentEvent(rl->GetEventNumber());
161 delete dovert;
162 return fCurrentVertex;
163 }
164
fe4280bf 165 if(!np1 || !np2) {
166 Error("FindVertexForCurrentEvent","No points in the pixer layers");
83028144 167 return fCurrentVertex;
fe4280bf 168 }
169 if(fDebug>0) cout << "N. points layer 1 and 2 : " << np1 << " " << np2 << endl;
170
171 Double_t asparx = (mxpiu-mxmeno)/(mxpiu+mxmeno);
172 Double_t aspary = (mypiu-mymeno)/(mypiu+mymeno);
173
174 Double_t x0 = 2.86*asparx;
175 Double_t y0 = 2.86*aspary;
176
177 if(fDebug>0) cout << "Rough xy vertex = " << x0 << " " << y0 << endl;
178
179 for(Int_t i=0;i<np1;i++) {
180 x1[i]-=x0;
181 y1[i]-=y0;
182 PhiFunc(x1[i],y1[i],phi1[i]);
183 r1[i]=TMath::Sqrt(x1[i]*x1[i]+y1[i]*y1[i]);
88cb7938 184 }
fe4280bf 185 for(Int_t i=0;i<np2;i++) {
186 x2[i]-=x0;
187 y2[i]-=y0;
188 PhiFunc(x2[i],y2[i],phi2[i]);
189 r2[i]=TMath::Sqrt(x2[i]*x2[i]+y2[i]*y2[i]);
83028144 190 }
fe4280bf 191
192 Int_t nbinxy=400;
193 Int_t nbinz=400;
194 TH1F *hxv=new TH1F("hxv","",nbinxy,-5,5);
195 TH1F *hyv=new TH1F("hyv","",nbinxy,-5,5);
196 TH1F *hzv=new TH1F("hzv","",nbinz,-15,15);
197
198 Double_t dphi;
199 for(Int_t j=0; j<np1; j++) {
200 for(Int_t k=0; k<np2; k++) {
201 dphi=TMath::Abs(phi2[k]-phi1[j]);
88cb7938 202 if(dphi>180) dphi = 360-dphi;
fe4280bf 203 if(dphi>fMaxDeltaPhi) continue;
204 hzv->Fill(z1[j]-(z2[k]-z1[j])/((r2[k]/r1[j])-1));
205 }
88cb7938 206 }
fe4280bf 207
208 // Fitting ...
209 Double_t max;
7d62fb64 210 Int_t binMax;
b5d712e2 211 Double_t maxcenter;
83028144 212
fe4280bf 213 max = hzv->GetMaximum();
7d62fb64 214 binMax=hzv->GetMaximumBin();
215 maxcenter=hzv->GetBinCenter(binMax);
fe4280bf 216 Double_t dxy=1.5;
83028144 217
b5d712e2 218 TF1 *fz = new TF1 ("fz","([0]*exp(-0.5*((x-[1])/[2])*((x-[1])/[2])))+[3]",maxcenter-dxy,maxcenter+dxy);
fe4280bf 219 fz->SetParameter(0,max);
b5d712e2 220 fz->SetParameter(1,maxcenter);
fe4280bf 221 fz->SetParameter(2,0.1);
222 fz->SetLineColor(kRed);
223 hzv->Fit("fz","RQ0");
224
225 if(TMath::Abs(x0)>0.4 || TMath::Abs(y0)>0.4) {
226 Warning("FindVertexForCurrentEvent","AliITSVertexerIonsGeneral finder is not reliable for events with x and y vertex coordinates too far from the origin of the reference system. Their values will be set to 0.\n");
227 Double_t position[3]={0,0,fz->GetParameter(1)};
228 Double_t resolution[3]={0,0,0};
229 Double_t snr[3]={0,0,0};
230 Char_t name[30];
231 if(fDebug>0)Info("FindVertexForCurrentEvent","Vertex found for event %d",evnumber);
232 sprintf(name,"Vertex");
233 fCurrentVertex = new AliESDVertex(position,resolution,snr,name);
83028144 234 return fCurrentVertex;
88cb7938 235 }
83028144 236
fe4280bf 237 for(Int_t j=0; j<np1; j++) {
238 for(Int_t k=0; k<np2; k++) {
239 if(TMath::Abs((z1[j]-(z2[k]-z1[j])/((r2[k]/r1[j])-1))-fz->GetParameter(1))>fMaxDeltaZ) continue;
240 if(y2[k]==y1[j]) continue;
241 hxv->Fill(x0+(x2[k]-((x2[k]-x1[j])/(y2[k]-y1[j]))*y2[k]));
242 if(x2[k]==x1[j]) continue;
243 hyv->Fill(y0+(y2[k]-((y2[k]-y1[j])/(x2[k]-x1[j]))*x2[k]));
88cb7938 244 }
245 }
fe4280bf 246
247 TF1 *fx = new TF1 ("fx","([0]*exp(-0.5*((x-[1])/[2])*((x-[1])/[2])))+[3]+[4]*x+[5]*x*x",x0-dxy,x0+dxy);
248 fx->SetParameter(0,100);
249 Double_t dist=0.3;
b5d712e2 250 Double_t xapprox=FindMaxAround(x0,hxv,dist);
251 if(fDebug>0) cout << "xapprox = " << xapprox << endl;
252 fx->SetParameter(1,xapprox);
253 Double_t difcentroid=0.07;
254 fx->SetParLimits(1,xapprox-difcentroid,xapprox+difcentroid);
fe4280bf 255 fx->SetParameter(2,0.1);
256 fx->SetLineColor(kRed);
257 hxv->Fit("fx","RQW0");
83028144 258
fe4280bf 259 TF1 *fy = new TF1 ("fy","([0]*exp(-0.5*((x-[1])/[2])*((x-[1])/[2])))+[3]+[4]*x+[5]*x*x",y0-dxy,y0+dxy);
260 fy->SetParameter(0,100);
b5d712e2 261 Double_t yapprox=FindMaxAround(y0,hyv,dist);
262 if(fDebug>0) cout << "yapprox = " << yapprox << endl;
263 fy->SetParameter(1,yapprox);
264 fy->SetParLimits(1,yapprox-difcentroid,yapprox+difcentroid);
fe4280bf 265 fy->SetParameter(2,0.1);
266 fy->SetLineColor(kRed);
267 hyv->Fit("fy","RQW0");
83028144 268
83028144 269 delete [] z1;
270 delete [] z2;
271 delete [] y1;
272 delete [] y2;
273 delete [] x1;
274 delete [] x2;
88cb7938 275 delete [] r1;
276 delete [] r2;
277 delete [] phi1;
278 delete [] phi2;
fe4280bf 279 delete hxv;
280 delete hyv;
281 delete hzv;
282
283 Double_t position[3]={fx->GetParameter(1),fy->GetParameter(1),fz->GetParameter(1)};
284 Double_t resolution[3]={fx->GetParameter(2),fy->GetParameter(2),fz->GetParameter(2)};
285 Double_t snr[3]={0,0,0};
286
88cb7938 287 Char_t name[30];
088e0b8d 288 if(fDebug>0)Info("FindVertexForCurrentEvent","Vertex found for event %d",evnumber);
88cb7938 289 sprintf(name,"Vertex");
d681bb2d 290 fCurrentVertex = new AliESDVertex(position,resolution,snr,name);
7d62fb64 291
88cb7938 292 return fCurrentVertex;
c5f0f3c1 293}
294
295//______________________________________________________________________
fe4280bf 296void AliITSVertexerIons::PhiFunc(Double_t &x,Double_t &y,Double_t &phi) {
b5d712e2 297 // Method for the computation of the phi angle given x and y. The result is in degrees.
fe4280bf 298 if(y>0 && x>0) phi=(TMath::ATan((Double_t)(y/x))*57.29578);
299 if(y>0 && x<0) phi=(TMath::ATan((Double_t)(y/x))*57.29578)+180;
300 if(y<0 && x<0) phi=(TMath::ATan((Double_t)(y/x))*57.29578)+180;
301 if(y<0 && x>0) phi=(TMath::ATan((Double_t)(y/x))*57.29578)+360;;
c5f0f3c1 302}
303
304//______________________________________________________________________
305void AliITSVertexerIons::FindVertices(){
306 // computes the vertices of the events in the range FirstEvent - LastEvent
88cb7938 307 AliRunLoader *rl = AliRunLoader::GetRunLoader();
83028144 308 AliITSLoader* itsloader = (AliITSLoader*) rl->GetLoader("ITSLoader");
309 itsloader->LoadRecPoints("read");
ab6a6f40 310 for(Int_t i=fFirstEvent;i<=fLastEvent;i++){
88cb7938 311 rl->GetEvent(i);
c5f0f3c1 312 FindVertexForCurrentEvent(i);
313 if(fCurrentVertex){
314 WriteCurrentVertex();
315 }
316 else {
317 if(fDebug>0){
318 cout<<"Vertex not found for event "<<i<<endl;
319 }
320 }
321 }
c5f0f3c1 322}
323
324//________________________________________________________
325void AliITSVertexerIons::PrintStatus() const {
326 // Print current status
88cb7938 327 cout <<"=======================================================\n";
328 cout <<" Debug flag: "<<fDebug<<endl;
c5f0f3c1 329 cout<<"First event to be processed "<<fFirstEvent;
330 cout<<"\n Last event to be processed "<<fLastEvent<<endl;
331 if(fCurrentVertex)fCurrentVertex->PrintStatus();
332}
fe4280bf 333
334//________________________________________________________
335Double_t AliITSVertexerIons::FindMaxAround(Double_t point, TH1F *h, Double_t distance) {
b5d712e2 336 // It finds the maximum of h within point-distance and distance+point.
337 Int_t maxcontent=0;
338 Int_t maxbin=0;
fe4280bf 339 for(Int_t i=0;i<h->GetNbinsX();i++) {
340 Int_t content=(Int_t)h->GetBinContent(i);
341 Double_t center=(Double_t)h->GetBinCenter(i);
1360bc6f 342 if(TMath::Abs(center-point)>distance) continue;
b5d712e2 343 if(content>maxcontent) {maxcontent=content;maxbin=i;}
fe4280bf 344 }
b5d712e2 345 Double_t max=h->GetBinCenter(maxbin);
fe4280bf 346 return max;
347}
348
349