1 #include "AliITSDetTypeRec.h"
2 #include "AliITSVertexerIons.h"
3 #include "AliITSVertexerZ.h"
4 #include "AliESDVertex.h"
5 #include "AliITSgeom.h"
6 #include "AliITSLoader.h"
7 #include "AliITSRecPoint.h"
13 //////////////////////////////////////////////////////////////////////
14 // AliITSVertexerIons is a class for full 3D primary vertex //
15 // finding optimized for Ion-Ion interactions //
20 // Written by Giuseppe Lo Re and Francesco Riggi //
21 // Giuseppe.Lore@ct.infn.it //
22 // Franco.Riggi@ct.infn.it //
24 // Release date: Mar 2004 //
27 //////////////////////////////////////////////////////////////////////
30 ClassImp(AliITSVertexerIons)
32 //______________________________________________________________________
33 AliITSVertexerIons::AliITSVertexerIons():AliITSVertexer(),
37 // Default Constructor
45 //______________________________________________________________________
46 AliITSVertexerIons::AliITSVertexerIons(TString fn):AliITSVertexer(fn),
50 // Standard constructor
58 //______________________________________________________________________
59 AliITSVertexerIons::AliITSVertexerIons(const AliITSVertexerIons &source):AliITSVertexer(source) {
61 // Copies are not allowed. The method is protected to avoid misuse.
62 Error("AliITSVertexerIons","Copy constructor not allowed\n");
65 //_________________________________________________________________________
66 //AliITSVertexerIons& AliITSVertexerIons::operator=(const AliITSVertexerIons &/*source*/) {
67 // Assignment operator
68 // Assignment is not allowed. The method is protected to avoid misuse.
69 //Error("= operator","Assignment operator not allowed\n");
74 //______________________________________________________________________
75 AliITSVertexerIons::~AliITSVertexerIons() {
80 //______________________________________________________________________
81 AliESDVertex* AliITSVertexerIons::FindVertexForCurrentEvent(Int_t evnumber){
82 // Defines the AliESDVertex for the current event
86 AliRunLoader *rl = AliRunLoader::GetRunLoader();
87 AliITSLoader* itsloader = (AliITSLoader*)rl->GetLoader("ITSLoader");
88 TDirectory * olddir = gDirectory;
90 AliITSgeom* g2 = (AliITSgeom*)gDirectory->Get("AliITSgeom");
93 TTree *tr = itsloader->TreeR();
94 AliITSDetTypeRec detTypeRec;
96 detTypeRec.SetTreeAddressR(tr);
98 TClonesArray *recpoints = detTypeRec.RecPoints();
103 Int_t nopoints1=40000;
104 Int_t nopoints2=40000;
107 Double_t *z1, *z2, *y1, *y2, *x1, *x2, *phi1, *phi2, *r1, *r2;
108 z1=new Double_t[nopoints1];
109 z2=new Double_t[nopoints2];
110 y1=new Double_t[nopoints1];
111 y2=new Double_t[nopoints2];
112 x1=new Double_t[nopoints1];
113 x2=new Double_t[nopoints2];
114 phi1=new Double_t[nopoints1];
115 phi2=new Double_t[nopoints2];
116 r1=new Double_t[nopoints1];
117 r2=new Double_t[nopoints2];
126 for(Int_t i=g2->GetStartSPD();i<=g2->GetLastSPD();i++) {
127 detTypeRec.ResetRecPoints();
129 npoints = recpoints->GetEntries();
130 for (Int_t ipoint=0;ipoint<npoints;ipoint++) {
131 pnt = (AliITSRecPoint*)recpoints->UncheckedAt(ipoint);
132 l[0]=pnt->GetDetLocalX();
134 l[2]=pnt->GetDetLocalZ();
136 r=TMath::Sqrt(TMath::Power(p[0],2)+TMath::Power(p[1],2));
138 if(i<80 && TMath::Abs(p[2])<14.35) {
148 if(i>=80 && TMath::Abs(p[2])<14.35) {
157 if(np1<fNpThreshold) {
158 Warning("FindVertexForCurrentEvent","AliITSVertexerIons finder is not reliable for low multiplicity events. Switching to AliITSVertexerZ with default parameters...\n");
159 Warning("FindVertexForCurrentEvent","N rec points = %d - Threshold is %d",np1,fNpThreshold);
160 AliITSVertexerZ *dovert = new AliITSVertexerZ("default");
161 fCurrentVertex =dovert->FindVertexForCurrentEvent(rl->GetEventNumber());
163 return fCurrentVertex;
167 Error("FindVertexForCurrentEvent","No points in the pixer layers");
168 return fCurrentVertex;
170 AliDebug(1,Form("N. points layer 1 and 2 : %d %d",np1,np2));
172 Double_t asparx = (mxpiu-mxmeno)/(mxpiu+mxmeno);
173 Double_t aspary = (mypiu-mymeno)/(mypiu+mymeno);
175 Double_t x0 = 2.86*asparx;
176 Double_t y0 = 2.86*aspary;
178 AliDebug(1,Form("Rough xy vertex = %f %f",x0,y0));
180 for(Int_t i=0;i<np1;i++) {
183 PhiFunc(x1[i],y1[i],phi1[i]);
184 r1[i]=TMath::Sqrt(x1[i]*x1[i]+y1[i]*y1[i]);
186 for(Int_t i=0;i<np2;i++) {
189 PhiFunc(x2[i],y2[i],phi2[i]);
190 r2[i]=TMath::Sqrt(x2[i]*x2[i]+y2[i]*y2[i]);
195 TH1F *hxv=new TH1F("hxv","",nbinxy,-5,5);
196 TH1F *hyv=new TH1F("hyv","",nbinxy,-5,5);
197 TH1F *hzv=new TH1F("hzv","",nbinz,-15,15);
200 for(Int_t j=0; j<np1; j++) {
201 for(Int_t k=0; k<np2; k++) {
202 dphi=TMath::Abs(phi2[k]-phi1[j]);
203 if(dphi>180) dphi = 360-dphi;
204 if(dphi>fMaxDeltaPhi) continue;
205 hzv->Fill(z1[j]-(z2[k]-z1[j])/((r2[k]/r1[j])-1));
214 max = hzv->GetMaximum();
215 binMax=hzv->GetMaximumBin();
216 maxcenter=hzv->GetBinCenter(binMax);
219 TF1 *fz = new TF1 ("fz","([0]*exp(-0.5*((x-[1])/[2])*((x-[1])/[2])))+[3]",maxcenter-dxy,maxcenter+dxy);
220 fz->SetParameter(0,max);
221 fz->SetParameter(1,maxcenter);
222 fz->SetParameter(2,0.1);
223 fz->SetLineColor(kRed);
224 hzv->Fit("fz","RQ0");
226 if(TMath::Abs(x0)>0.4 || TMath::Abs(y0)>0.4) {
227 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");
228 Double_t position[3]={0,0,fz->GetParameter(1)};
229 Double_t resolution[3]={0,0,0};
230 Double_t snr[3]={0,0,0};
232 AliDebug(1,Form("Vertex found for event %d",evnumber));
233 sprintf(name,"Vertex");
234 fCurrentVertex = new AliESDVertex(position,resolution,snr,name);
235 return fCurrentVertex;
238 for(Int_t j=0; j<np1; j++) {
239 for(Int_t k=0; k<np2; k++) {
240 if(TMath::Abs((z1[j]-(z2[k]-z1[j])/((r2[k]/r1[j])-1))-fz->GetParameter(1))>fMaxDeltaZ) continue;
241 if(y2[k]==y1[j]) continue;
242 hxv->Fill(x0+(x2[k]-((x2[k]-x1[j])/(y2[k]-y1[j]))*y2[k]));
243 if(x2[k]==x1[j]) continue;
244 hyv->Fill(y0+(y2[k]-((y2[k]-y1[j])/(x2[k]-x1[j]))*x2[k]));
248 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);
249 fx->SetParameter(0,100);
251 Double_t xapprox=FindMaxAround(x0,hxv,dist);
252 AliDebug(1,Form("xapprox = %f",xapprox));
253 fx->SetParameter(1,xapprox);
254 Double_t difcentroid=0.07;
255 fx->SetParLimits(1,xapprox-difcentroid,xapprox+difcentroid);
256 fx->SetParameter(2,0.1);
257 fx->SetLineColor(kRed);
258 hxv->Fit("fx","RQW0");
260 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);
261 fy->SetParameter(0,100);
262 Double_t yapprox=FindMaxAround(y0,hyv,dist);
263 AliDebug(1,Form("yapprox = %f",yapprox));
264 fy->SetParameter(1,yapprox);
265 fy->SetParLimits(1,yapprox-difcentroid,yapprox+difcentroid);
266 fy->SetParameter(2,0.1);
267 fy->SetLineColor(kRed);
268 hyv->Fit("fy","RQW0");
284 Double_t position[3]={fx->GetParameter(1),fy->GetParameter(1),fz->GetParameter(1)};
285 Double_t resolution[3]={fx->GetParameter(2),fy->GetParameter(2),fz->GetParameter(2)};
286 Double_t snr[3]={0,0,0};
289 AliDebug(1,Form("Vertex found for event %d",evnumber));
290 sprintf(name,"Vertex");
291 fCurrentVertex = new AliESDVertex(position,resolution,snr,name);
293 return fCurrentVertex;
296 //______________________________________________________________________
297 void AliITSVertexerIons::PhiFunc(Double_t &x,Double_t &y,Double_t &phi) {
298 // Method for the computation of the phi angle given x and y. The result is in degrees.
299 if(y>0 && x>0) phi=(TMath::ATan((Double_t)(y/x))*57.29578);
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)+180;
302 if(y<0 && x>0) phi=(TMath::ATan((Double_t)(y/x))*57.29578)+360;;
305 //______________________________________________________________________
306 void AliITSVertexerIons::FindVertices(){
307 // computes the vertices of the events in the range FirstEvent - LastEvent
308 AliRunLoader *rl = AliRunLoader::GetRunLoader();
309 AliITSLoader* itsloader = (AliITSLoader*) rl->GetLoader("ITSLoader");
310 itsloader->LoadRecPoints("read");
311 for(Int_t i=fFirstEvent;i<=fLastEvent;i++){
313 FindVertexForCurrentEvent(i);
315 WriteCurrentVertex();
318 AliDebug(1,Form("Vertex not found for event %d",i));
323 //________________________________________________________
324 void AliITSVertexerIons::PrintStatus() const {
325 // Print current status
326 cout <<"=======================================================\n";
327 cout<<"First event to be processed "<<fFirstEvent;
328 cout<<"\n Last event to be processed "<<fLastEvent<<endl;
329 if(fCurrentVertex)fCurrentVertex->PrintStatus();
332 //________________________________________________________
333 Double_t AliITSVertexerIons::FindMaxAround(Double_t point, TH1F *h, Double_t distance) {
334 // It finds the maximum of h within point-distance and distance+point.
337 for(Int_t i=0;i<h->GetNbinsX();i++) {
338 Int_t content=(Int_t)h->GetBinContent(i);
339 Double_t center=(Double_t)h->GetBinCenter(i);
340 if(TMath::Abs(center-point)>distance) continue;
341 if(content>maxcontent) {maxcontent=content;maxbin=i;}
343 Double_t max=h->GetBinCenter(maxbin);