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