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