]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSVertexer3D.cxx
vertexerZ called when vertexer3D fails (F. Prino)
[u/mrichter/AliRoot.git] / ITS / AliITSVertexer3D.cxx
1 /**************************************************************************
2  * Copyright(c) 2006-2008, 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 #include <TH3F.h>
16 #include <TTree.h>
17 #include "AliESDVertex.h"
18 #include "AliLog.h"
19 #include "AliStrLine.h"
20 #include "AliTracker.h"
21 #include "AliITSDetTypeRec.h"
22 #include "AliITSRecPoint.h"
23 #include "AliITSgeomTGeo.h"
24 #include "AliVertexerTracks.h"
25 #include "AliITSVertexer3D.h"
26 #include "AliITSVertexerZ.h"
27 /////////////////////////////////////////////////////////////////
28 // this class implements a method to determine
29 // the 3 coordinates of the primary vertex
30 // for p-p collisions 
31 // It can be used successfully with Pb-Pb collisions
32 ////////////////////////////////////////////////////////////////
33
34 ClassImp(AliITSVertexer3D)
35
36 /* $Id$ */
37
38 //______________________________________________________________________
39 AliITSVertexer3D::AliITSVertexer3D():AliITSVertexer(),
40 fLines("AliStrLine",1000),
41 fVert3D(),
42 fCoarseDiffPhiCut(0.),
43 fCoarseMaxRCut(0.),
44 fMaxRCut(0.),
45 fZCutDiamond(0.),
46 fMaxZCut(0.),
47 fDCAcut(0.),
48 fDiffPhiMax(0.),
49 fMeanPSelTrk(0.),
50 fMeanPtSelTrk(0.)
51 {
52   // Default constructor
53   SetCoarseDiffPhiCut();
54   SetCoarseMaxRCut();
55   SetMaxRCut();
56   SetZCutDiamond();
57   SetMaxZCut();
58   SetDCAcut();
59   SetDiffPhiMax();
60   SetMeanPSelTracks();
61   SetMeanPtSelTracks();
62 }
63
64 //______________________________________________________________________
65 AliITSVertexer3D::~AliITSVertexer3D() {
66   // Destructor
67   fLines.Clear("C");
68 }
69
70 //______________________________________________________________________
71 void AliITSVertexer3D::ResetVert3D(){
72   //
73   fVert3D.SetXv(0.);
74   fVert3D.SetYv(0.);
75   fVert3D.SetZv(0.);
76   fVert3D.SetDispersion(0.);
77   fVert3D.SetNContributors(0);
78 }
79 //______________________________________________________________________
80 AliESDVertex* AliITSVertexer3D::FindVertexForCurrentEvent(TTree *itsClusterTree){
81   // Defines the AliESDVertex for the current event
82   ResetVert3D();
83   AliDebug(1,"FindVertexForCurrentEvent - 3D - PROCESSING NEXT EVENT");
84   fLines.Clear();
85
86   Int_t nolines = FindTracklets(itsClusterTree,0);
87   fCurrentVertex = 0;
88   if(nolines>=2){
89     Int_t rc=Prepare3DVertex(0);
90     if(rc==0) fVert3D=AliVertexerTracks::TrackletVertexFinder(&fLines,0);
91   /*  uncomment to debug
92       printf("Vertex found in first iteration:\n");
93       fVert3D.Print();
94       printf("Start second iteration\n");
95   end of debug lines  */
96     if(fVert3D.GetNContributors()>0){
97       fLines.Clear("C");
98       nolines = FindTracklets(itsClusterTree,1);
99       if(nolines>=2){
100         rc=Prepare3DVertex(1);
101         if(rc==0) fVert3D=AliVertexerTracks::TrackletVertexFinder(&fLines,0);
102       }
103     }
104   /*  uncomment to debug 
105       printf("Vertex found in second iteration:\n");
106       fVert3D.Print();
107       end of debug lines  */ 
108  
109     Float_t vRadius=TMath::Sqrt(fVert3D.GetXv()*fVert3D.GetXv()+fVert3D.GetYv()*fVert3D.GetYv());
110     if(vRadius<GetPipeRadius() && fVert3D.GetNContributors()>0){
111       Double_t position[3]={fVert3D.GetXv(),fVert3D.GetYv(),fVert3D.GetZv()};
112       Double_t covmatrix[6];
113       fVert3D.GetCovMatrix(covmatrix);
114       Double_t chi2=99999.;
115       Int_t    nContr=fVert3D.GetNContributors();
116       fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);    
117       fCurrentVertex->SetTitle("vertexer: 3D");
118       fCurrentVertex->SetName("SPDVertex3D");
119       fCurrentVertex->SetDispersion(fVert3D.GetDispersion());
120     }
121   }
122   if(!fCurrentVertex){
123     AliITSVertexerZ vertz(GetNominalPos()[0],GetNominalPos()[1]);
124     AliDebug(1,"Call Vertexer Z\n");
125     AliESDVertex* vtxz = vertz.FindVertexForCurrentEvent(itsClusterTree);
126     if(vtxz){
127       Double_t position[3]={GetNominalPos()[0],GetNominalPos()[1],vtxz->GetZv()};
128       Double_t covmatrix[6];
129       vtxz->GetCovMatrix(covmatrix);
130       Double_t chi2=99999.;
131       Int_t    nContr=vtxz->GetNContributors();
132       fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);    
133       fCurrentVertex->SetTitle("vertexer: Z");
134       fCurrentVertex->SetName("SPDVertexZ");
135       delete vtxz;
136       //      printf("Vertexer Z success\n");
137     }
138
139   }
140   FindMultiplicity(itsClusterTree);
141   return fCurrentVertex;
142 }  
143
144 //______________________________________________________________________
145 Int_t AliITSVertexer3D::FindTracklets(TTree *itsClusterTree, Int_t optCuts){
146   // All the possible combinations between recpoints on layer 1and 2 are
147   // considered. Straight lines (=tracklets)are formed. 
148   // The tracklets are processed in Prepare3DVertex
149   AliITSDetTypeRec detTypeRec;
150
151   TTree *tR = itsClusterTree;
152   detTypeRec.SetTreeAddressR(tR);
153   TClonesArray *itsRec  = 0;
154   // lc1 and gc1 are local and global coordinates for layer 1
155   //  Float_t lc1[3]={0.,0.,0.};
156   Float_t gc1[3]={0.,0.,0.};
157   // lc2 and gc2 are local and global coordinates for layer 2
158   //  Float_t lc2[3]={0.,0.,0.};
159   Float_t gc2[3]={0.,0.,0.};
160
161   itsRec = detTypeRec.RecPoints();
162   TBranch *branch;
163   branch = tR->GetBranch("ITSRecPoints");
164
165   // Set values for cuts
166   Float_t xbeam=0., ybeam=0.;
167   Float_t zvert=0.;
168   Float_t deltaPhi=fCoarseDiffPhiCut;
169   Float_t deltaR=fCoarseMaxRCut;
170   Float_t dZmax=fZCutDiamond;
171   if(optCuts){
172     xbeam=fVert3D.GetXv();
173     ybeam=fVert3D.GetYv();
174     zvert=fVert3D.GetZv();
175     deltaPhi = fDiffPhiMax; 
176     deltaR=fMaxRCut;
177     dZmax=fMaxZCut;
178   }
179   Int_t nrpL1 = 0;    // number of rec points on layer 1
180   Int_t nrpL2 = 0;    // number of rec points on layer 2
181
182   // By default irstL1=0 and lastL1=79
183   Int_t firstL1 = AliITSgeomTGeo::GetModuleIndex(1,1,1);
184   Int_t lastL1 = AliITSgeomTGeo::GetModuleIndex(2,1,1)-1;
185   for(Int_t module= firstL1; module<=lastL1;module++){  // count number of recopints on layer 1
186     branch->GetEvent(module);
187     nrpL1+= itsRec->GetEntries();
188     detTypeRec.ResetRecPoints();
189   }
190   //By default firstL2=80 and lastL2=239
191   Int_t firstL2 = AliITSgeomTGeo::GetModuleIndex(2,1,1);
192   Int_t lastL2 = AliITSgeomTGeo::GetModuleIndex(3,1,1)-1;
193   for(Int_t module= firstL2; module<=lastL2;module++){  // count number of recopints on layer 2
194     branch->GetEvent(module);
195     nrpL2+= itsRec->GetEntries();
196     detTypeRec.ResetRecPoints();
197   }
198   if(nrpL1 == 0 || nrpL2 == 0){
199     return -1;
200   }
201   AliDebug(1,Form("RecPoints on Layer 1,2 = %d, %d\n",nrpL1,nrpL2));
202
203   Double_t a[3]={xbeam,ybeam,0.}; 
204   Double_t b[3]={xbeam,ybeam,10.};
205   AliStrLine zeta(a,b,kTRUE);
206   Float_t bField=AliTracker::GetBz()/10.; //T
207   SetMeanPPtSelTracks(bField);
208
209   Int_t nolines = 0;
210   // Loop on modules of layer 1
211   for(Int_t modul1= firstL1; modul1<=lastL1;modul1++){   // Loop on modules of layer 1
212     UShort_t ladder=int(modul1/4)+1; // ladders are numbered starting from 1
213     branch->GetEvent(modul1);
214     Int_t nrecp1 = itsRec->GetEntries();
215     static TClonesArray prpl1("AliITSRecPoint",nrecp1);
216     prpl1.SetOwner();
217     for(Int_t j=0;j<nrecp1;j++){
218       AliITSRecPoint *recp = (AliITSRecPoint*)itsRec->At(j);
219       new(prpl1[j])AliITSRecPoint(*recp);
220     }
221     detTypeRec.ResetRecPoints();
222     for(Int_t j=0;j<nrecp1;j++){
223       AliITSRecPoint *recp1 = (AliITSRecPoint*)prpl1.At(j);
224       // Local coordinates of this recpoint
225       /*
226       lc[0]=recp1->GetDetLocalX();
227       lc[2]=recp1->GetDetLocalZ();
228       */
229       recp1->GetGlobalXYZ(gc1);
230       Double_t phi1 = TMath::ATan2(gc1[1]-ybeam,gc1[0]-xbeam);
231       if(phi1<0)phi1=2*TMath::Pi()+phi1;
232       for(Int_t ladl2=0 ; ladl2<fLadOnLay2*2+1;ladl2++){
233         for(Int_t k=0;k<4;k++){
234           Int_t ladmod=fLadders[ladder-1]+ladl2;
235           if(ladmod>AliITSgeomTGeo::GetNLadders(2)) ladmod=ladmod-AliITSgeomTGeo::GetNLadders(2);
236           Int_t modul2=AliITSgeomTGeo::GetModuleIndex(2,ladmod,k+1);
237           branch->GetEvent(modul2);
238           Int_t nrecp2 = itsRec->GetEntries();
239           for(Int_t j2=0;j2<nrecp2;j2++){
240             AliITSRecPoint *recp2 = (AliITSRecPoint*)itsRec->At(j2);
241             /*
242             lc2[0]=recp2->GetDetLocalX();
243             lc2[2]=recp2->GetDetLocalZ();
244             */
245             recp2->GetGlobalXYZ(gc2);
246             Double_t phi2 = TMath::ATan2(gc2[1]-ybeam,gc2[0]-xbeam);
247             if(phi2<0)phi2=2*TMath::Pi()+phi2;
248             Double_t diff = TMath::Abs(phi2-phi1); 
249             if(diff>TMath::Pi())diff=2.*TMath::Pi()-diff; 
250             if(diff>deltaPhi)continue;
251             AliStrLine line(gc1,gc2,kTRUE);
252             Double_t cp[3];
253             Int_t retcode = line.Cross(&zeta,cp);
254             if(retcode<0)continue;
255             Double_t dca = line.GetDCA(&zeta);
256             if(dca<0.) continue;
257             if(dca>deltaR)continue;
258             Double_t deltaZ=cp[2]-zvert;
259             if(TMath::Abs(deltaZ)>dZmax)continue;
260
261             if(nolines == 0){
262               if(fLines.GetEntriesFast()>0)fLines.Clear("C");
263             }
264             Float_t cov[6];
265             recp2->GetGlobalCov(cov);
266
267             
268             Float_t rad1=TMath::Sqrt(gc1[0]*gc1[0]+gc1[1]*gc1[1]);
269             Float_t rad2=TMath::Sqrt(gc2[0]*gc2[0]+gc2[1]*gc2[1]);
270             Float_t factor=(rad1+rad2)/(rad2-rad1); //factor to account for error on tracklet direction 
271
272             Float_t curvErr=0;
273             if(bField>0.00001){
274               Float_t curvRadius=fMeanPtSelTrk/(0.3*bField)*100; //cm 
275               Float_t dRad=TMath::Sqrt(TMath::Power((gc1[0]-gc2[0]),2)+TMath::Power((gc1[1]-gc2[1]),2));
276               Float_t aux=dRad/2.+rad1;
277               curvErr=TMath::Sqrt(curvRadius*curvRadius-dRad*dRad/4.)-TMath::Sqrt(curvRadius*curvRadius-aux*aux); //cm
278             }
279
280             Float_t sigmasq[3];
281             sigmasq[0]=(cov[0]+curvErr*curvErr/2.)*factor*factor;
282             sigmasq[1]=(cov[3]+curvErr*curvErr/2.)*factor*factor;
283             sigmasq[2]=cov[5]*factor*factor;
284
285             // Multiple scattering
286             Float_t beta=1.;
287             Float_t beta2=beta*beta;
288             Float_t p2=fMeanPSelTrk*fMeanPSelTrk;
289             Float_t rBP=GetPipeRadius();
290             Float_t dBP=0.08/35.3; // 800 um of Be
291             Float_t dL1=0.01; //approx. 1% of radiation length  
292             Float_t theta2BP=14.1*14.1/(beta2*p2*1e6)*TMath::Abs(dBP);
293             Float_t theta2L1=14.1*14.1/(beta2*p2*1e6)*TMath::Abs(dL1);
294             Float_t thetaBP=TMath::Sqrt(theta2BP);
295             Float_t thetaL1=TMath::Sqrt(theta2L1);
296 //          Float_t geomfac[3];
297 //          geomfac[0]=sin(phi1)*sin(phi1);
298 //          geomfac[1]=cos(phi1)*cos(phi1);
299 //          Float_t tgth=(gc2[2]-gc1[2])/(rad2-rad1);       
300 //          geomfac[2]=1+tgth*tgth;
301             for(Int_t ico=0; ico<3;ico++){    
302 //            printf("Error on coord. %d due to cov matrix+curvErr=%f\n",ico,sigmasq[ico]);
303 //            //              sigmasq[ico]+=rad1*rad1*geomfac[ico]*theta2L1/2; // multiple scattering in layer 1
304 //            //  sigmasq[ico]+=rBP*rBP*geomfac[ico]*theta2BP/2; // multiple scattering in beam pipe
305               sigmasq[ico]+=TMath::Power(rad1*TMath::Tan(thetaL1),2)/3.;
306               sigmasq[ico]+=TMath::Power(rBP*TMath::Tan(thetaBP),2)/3.;
307
308 //            printf("Multipl. scatt. contr %d = %f (LAY1), %f (BP)\n",ico,rad1*rad1*geomfac[ico]*theta2L1/2,rBP*rBP*geomfac[ico]*theta2BP/2);
309 //            printf("Total error on coord %d = %f\n",ico,sigmasq[ico]);
310             }
311             Float_t wmat[9]={1.,0.,0.,0.,1.,0.,0.,0.,1.};
312             if(sigmasq[0]!=0.) wmat[0]=1./sigmasq[0];
313             if(sigmasq[1]!=0.) wmat[4]=1./sigmasq[1];
314             if(sigmasq[2]!=0.) wmat[8]=1./sigmasq[2];
315             new(fLines[nolines++])AliStrLine(gc1,sigmasq,wmat,gc2,kTRUE);
316
317           }
318           detTypeRec.ResetRecPoints();
319         }
320       }
321     }
322     prpl1.Clear();
323   }
324   if(nolines == 0)return -2;
325   return nolines;
326 }
327
328 //______________________________________________________________________
329 Int_t  AliITSVertexer3D::Prepare3DVertex(Int_t optCuts){
330   // Finds the 3D vertex information using tracklets
331   Int_t retcode = -1;
332
333   Float_t xbeam=0.;
334   Float_t ybeam=0.;
335   Float_t zvert=0.;
336   Float_t deltaR=fCoarseMaxRCut;
337   Float_t dZmax=fZCutDiamond;
338   if(optCuts){
339     xbeam=fVert3D.GetXv();
340     ybeam=fVert3D.GetYv();
341     zvert=fVert3D.GetZv();
342     deltaR=fMaxRCut;
343     dZmax=fMaxZCut;
344   }
345
346   Int_t nbr=50;
347   Float_t rl=-fCoarseMaxRCut;
348   Float_t rh=fCoarseMaxRCut;
349   Int_t nbz=100;
350   Float_t zl=-fZCutDiamond;
351   Float_t zh=fZCutDiamond;
352   Float_t binsizer=(rh-rl)/nbr;
353   Float_t binsizez=(zh-zl)/nbz;
354   TH3F *h3d = new TH3F("h3d","xyz distribution",nbr,rl,rh,nbr,rl,rh,nbz,zl,zh);
355
356   // cleanup of the TCLonesArray of tracklets (i.e. fakes are removed)
357   Int_t *validate = new Int_t [fLines.GetEntriesFast()];
358   for(Int_t i=0; i<fLines.GetEntriesFast();i++)validate[i]=0;
359   for(Int_t i=0; i<fLines.GetEntriesFast()-1;i++){
360     if(validate[i]==1)continue;
361     AliStrLine *l1 = (AliStrLine*)fLines.At(i);
362     for(Int_t j=i+1;j<fLines.GetEntriesFast();j++){
363       AliStrLine *l2 = (AliStrLine*)fLines.At(j);
364       Float_t dca=l1->GetDCA(l2);
365       if(dca > fDCAcut || dca<0.00001) continue;
366       Double_t point[3];
367       Int_t retc = l1->Cross(l2,point);
368       if(retc<0)continue;
369       Double_t deltaZ=point[2]-zvert;
370      if(TMath::Abs(deltaZ)>dZmax)continue;
371       Double_t rad=TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
372       if(rad>fCoarseMaxRCut)continue;
373       Double_t deltaX=point[0]-xbeam;
374       Double_t deltaY=point[1]-ybeam;
375       Double_t raddist=TMath::Sqrt(deltaX*deltaX+deltaY*deltaY);
376       if(raddist>deltaR)continue;
377       validate[i]=1;
378       validate[j]=1;
379       h3d->Fill(point[0],point[1],point[2]);
380     }
381   }
382
383
384
385   Int_t numbtracklets=0;
386   for(Int_t i=0; i<fLines.GetEntriesFast();i++)if(validate[i]>=1)numbtracklets++;
387   if(numbtracklets<2){delete [] validate; delete h3d; return retcode; }
388
389   for(Int_t i=0; i<fLines.GetEntriesFast();i++){
390     if(validate[i]<1)fLines.RemoveAt(i);
391   }
392   fLines.Compress();
393   AliDebug(1,Form("Number of tracklets (after compress)%d ",fLines.GetEntriesFast()));
394   delete [] validate;
395
396
397   // finds peak in histo
398   TAxis *xax = h3d->GetXaxis();  
399   TAxis *yax = h3d->GetYaxis();
400   TAxis *zax = h3d->GetZaxis();
401   Double_t peak[3]={0.,0.,0.};
402   Float_t contref = 0.;
403   for(Int_t i=xax->GetFirst();i<=xax->GetLast();i++){
404     Float_t xval = xax->GetBinCenter(i);
405     for(Int_t j=yax->GetFirst();j<=yax->GetLast();j++){
406       Float_t yval = yax->GetBinCenter(j);
407       for(Int_t k=zax->GetFirst();k<=zax->GetLast();k++){
408         Float_t bc = h3d->GetBinContent(i,j,k);
409         Float_t zval = zax->GetBinCenter(k);
410         if(bc>contref){
411           contref = bc;
412           peak[2] = zval;
413           peak[1] = yval;
414           peak[0] = xval;
415         }
416       }
417     }
418   }
419   delete h3d;
420
421   //         Second selection loop
422   Float_t bs=(binsizer+binsizez)/2.;
423   for(Int_t i=0; i<fLines.GetEntriesFast();i++){
424     AliStrLine *l1 = (AliStrLine*)fLines.At(i);
425     if(l1->GetDistFromPoint(peak)>2.5*bs)fLines.RemoveAt(i);
426   }
427   fLines.Compress();
428   AliDebug(1,Form("Number of tracklets (after 2nd compression) %d",fLines.GetEntriesFast()));
429
430   if(fLines.GetEntriesFast()>1){
431     //  find a first candidate for the primary vertex
432     fVert3D=AliVertexerTracks::TrackletVertexFinder(&fLines,0); 
433     // make a further selection on tracklets based on this first candidate
434     fVert3D.GetXYZ(peak);
435     AliDebug(1,Form("FIRST V candidate: %f ; %f ; %f",peak[0],peak[1],peak[2]));
436     for(Int_t i=0; i<fLines.GetEntriesFast();i++){
437       AliStrLine *l1 = (AliStrLine*)fLines.At(i);
438       if(l1->GetDistFromPoint(peak)> fDCAcut)fLines.RemoveAt(i);
439     }
440     fLines.Compress();
441     AliDebug(1,Form("Number of tracklets (after 3rd compression) %d",fLines.GetEntriesFast()));
442     if(fLines.GetEntriesFast()>1) retcode=0; // this new tracklet selection is used
443     else retcode =1; // the previous tracklet selection will be used
444   }
445   else {
446     retcode = 0;
447   }
448   return retcode;  
449 }
450
451 //________________________________________________________
452 void AliITSVertexer3D::SetMeanPPtSelTracks(Float_t fieldTesla){
453   // Sets mean values of P and Pt based on the field
454   if(TMath::Abs(fieldTesla-0.5)<0.01){
455     SetMeanPSelTracks(0.885);
456     SetMeanPtSelTracks(0.630);
457   }else if(TMath::Abs(fieldTesla-0.4)<0.01){
458     SetMeanPSelTracks(0.805);
459     SetMeanPtSelTracks(0.580);
460   }else if(TMath::Abs(fieldTesla-0.2)<0.01){
461     SetMeanPSelTracks(0.740);
462     SetMeanPtSelTracks(0.530);
463   }else if(fieldTesla<0.00001){
464     SetMeanPSelTracks(0.730);
465     SetMeanPtSelTracks(0.510);
466   }else{
467     SetMeanPSelTracks();
468     SetMeanPtSelTracks();
469   }
470 }
471
472
473 //________________________________________________________
474 void AliITSVertexer3D::PrintStatus() const {
475   // Print current status
476   printf("=======================================================\n");
477   printf("Loose cut on Delta Phi %f\n",fCoarseDiffPhiCut);
478   printf("Cut on tracklet DCA to Z axis %f\n",fCoarseMaxRCut);
479   printf("Cut on tracklet DCA to beam axis %f\n",fMaxRCut);
480   printf("Cut on diamond (Z) %f\n",fZCutDiamond);
481   printf("Cut on DCA - tracklet to tracklet and to vertex %f\n",fDCAcut);
482   printf(" Max Phi difference: %f\n",fDiffPhiMax);
483   printf("=======================================================\n");
484 }