]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSVertexer3D.cxx
19cd1a65b1fc63d988268ccf1d3b29502d315460
[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     AliDebug(1,"Call Vertexer Z\n");
124     AliESDVertex* vtxz = vertz.FindVertexForCurrentEvent(itsClusterTree);
125     if(vtxz){
126       Double_t position[3]={GetNominalPos()[0],GetNominalPos()[1],vtxz->GetZv()};
127       Double_t covmatrix[6];
128       vtxz->GetCovMatrix(covmatrix);
129       Double_t chi2=99999.;
130       Int_t    nContr=vtxz->GetNContributors();
131       fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);    
132       fCurrentVertex->SetTitle("vertexer: Z");
133       fCurrentVertex->SetName("SPDVertexZ");
134       delete vtxz;
135       //      printf("Vertexer Z success\n");
136     }
137
138   }
139   FindMultiplicity(itsClusterTree);
140   return fCurrentVertex;
141 }  
142
143 //______________________________________________________________________
144 Int_t AliITSVertexer3D::FindTracklets(TTree *itsClusterTree, Int_t optCuts){
145   // All the possible combinations between recpoints on layer 1and 2 are
146   // considered. Straight lines (=tracklets)are formed. 
147   // The tracklets are processed in Prepare3DVertex
148   AliITSDetTypeRec detTypeRec;
149
150   TTree *tR = itsClusterTree;
151   detTypeRec.SetTreeAddressR(tR);
152   TClonesArray *itsRec  = 0;
153   // lc1 and gc1 are local and global coordinates for layer 1
154   //  Float_t lc1[3]={0.,0.,0.};
155   Float_t gc1[3]={0.,0.,0.};
156   // lc2 and gc2 are local and global coordinates for layer 2
157   //  Float_t lc2[3]={0.,0.,0.};
158   Float_t gc2[3]={0.,0.,0.};
159
160   itsRec = detTypeRec.RecPoints();
161   TBranch *branch;
162   branch = tR->GetBranch("ITSRecPoints");
163
164   // Set values for cuts
165   Float_t xbeam=GetNominalPos()[0]; 
166   Float_t ybeam=GetNominalPos()[1];
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     if(!fUseModule[modul1]) continue;
213     UShort_t ladder=int(modul1/4)+1; // ladders are numbered starting from 1
214     branch->GetEvent(modul1);
215     Int_t nrecp1 = itsRec->GetEntries();
216     static TClonesArray prpl1("AliITSRecPoint",nrecp1);
217     prpl1.SetOwner();
218     for(Int_t j=0;j<nrecp1;j++){
219       AliITSRecPoint *recp = (AliITSRecPoint*)itsRec->At(j);
220       new(prpl1[j])AliITSRecPoint(*recp);
221     }
222     detTypeRec.ResetRecPoints();
223     for(Int_t j=0;j<nrecp1;j++){
224       AliITSRecPoint *recp1 = (AliITSRecPoint*)prpl1.At(j);
225       // Local coordinates of this recpoint
226       /*
227       lc[0]=recp1->GetDetLocalX();
228       lc[2]=recp1->GetDetLocalZ();
229       */
230       recp1->GetGlobalXYZ(gc1);
231       Double_t phi1 = TMath::ATan2(gc1[1]-ybeam,gc1[0]-xbeam);
232       if(phi1<0)phi1=2*TMath::Pi()+phi1;
233       for(Int_t ladl2=0 ; ladl2<fLadOnLay2*2+1;ladl2++){
234         for(Int_t k=0;k<4;k++){
235           Int_t ladmod=fLadders[ladder-1]+ladl2;
236           if(ladmod>AliITSgeomTGeo::GetNLadders(2)) ladmod=ladmod-AliITSgeomTGeo::GetNLadders(2);
237           Int_t modul2=AliITSgeomTGeo::GetModuleIndex(2,ladmod,k+1);
238           if(!fUseModule[modul2]) continue;
239           branch->GetEvent(modul2);
240           Int_t nrecp2 = itsRec->GetEntries();
241           for(Int_t j2=0;j2<nrecp2;j2++){
242             AliITSRecPoint *recp2 = (AliITSRecPoint*)itsRec->At(j2);
243             /*
244             lc2[0]=recp2->GetDetLocalX();
245             lc2[2]=recp2->GetDetLocalZ();
246             */
247             recp2->GetGlobalXYZ(gc2);
248             Double_t phi2 = TMath::ATan2(gc2[1]-ybeam,gc2[0]-xbeam);
249             if(phi2<0)phi2=2*TMath::Pi()+phi2;
250             Double_t diff = TMath::Abs(phi2-phi1); 
251             if(diff>TMath::Pi())diff=2.*TMath::Pi()-diff; 
252             if(diff>deltaPhi)continue;
253             AliStrLine line(gc1,gc2,kTRUE);
254             Double_t cp[3];
255             Int_t retcode = line.Cross(&zeta,cp);
256             if(retcode<0)continue;
257             Double_t dca = line.GetDCA(&zeta);
258             if(dca<0.) continue;
259             if(dca>deltaR)continue;
260             Double_t deltaZ=cp[2]-zvert;
261             if(TMath::Abs(deltaZ)>dZmax)continue;
262
263             if(nolines == 0){
264               if(fLines.GetEntriesFast()>0)fLines.Clear("C");
265             }
266             Float_t cov[6];
267             recp2->GetGlobalCov(cov);
268
269             
270             Float_t rad1=TMath::Sqrt(gc1[0]*gc1[0]+gc1[1]*gc1[1]);
271             Float_t rad2=TMath::Sqrt(gc2[0]*gc2[0]+gc2[1]*gc2[1]);
272             Float_t factor=(rad1+rad2)/(rad2-rad1); //factor to account for error on tracklet direction 
273
274             Float_t curvErr=0;
275             if(bField>0.00001){
276               Float_t curvRadius=fMeanPtSelTrk/(0.3*bField)*100; //cm 
277               Float_t dRad=TMath::Sqrt(TMath::Power((gc1[0]-gc2[0]),2)+TMath::Power((gc1[1]-gc2[1]),2));
278               Float_t aux=dRad/2.+rad1;
279               curvErr=TMath::Sqrt(curvRadius*curvRadius-dRad*dRad/4.)-TMath::Sqrt(curvRadius*curvRadius-aux*aux); //cm
280             }
281
282             Float_t sigmasq[3];
283             sigmasq[0]=(cov[0]+curvErr*curvErr/2.)*factor*factor;
284             sigmasq[1]=(cov[3]+curvErr*curvErr/2.)*factor*factor;
285             sigmasq[2]=cov[5]*factor*factor;
286
287             // Multiple scattering
288             Float_t beta=1.;
289             Float_t beta2=beta*beta;
290             Float_t p2=fMeanPSelTrk*fMeanPSelTrk;
291             Float_t rBP=GetPipeRadius();
292             Float_t dBP=0.08/35.3; // 800 um of Be
293             Float_t dL1=0.01; //approx. 1% of radiation length  
294             Float_t theta2BP=14.1*14.1/(beta2*p2*1e6)*TMath::Abs(dBP);
295             Float_t theta2L1=14.1*14.1/(beta2*p2*1e6)*TMath::Abs(dL1);
296             Float_t thetaBP=TMath::Sqrt(theta2BP);
297             Float_t thetaL1=TMath::Sqrt(theta2L1);
298 //          Float_t geomfac[3];
299 //          geomfac[0]=sin(phi1)*sin(phi1);
300 //          geomfac[1]=cos(phi1)*cos(phi1);
301 //          Float_t tgth=(gc2[2]-gc1[2])/(rad2-rad1);       
302 //          geomfac[2]=1+tgth*tgth;
303             for(Int_t ico=0; ico<3;ico++){    
304 //            printf("Error on coord. %d due to cov matrix+curvErr=%f\n",ico,sigmasq[ico]);
305 //            //              sigmasq[ico]+=rad1*rad1*geomfac[ico]*theta2L1/2; // multiple scattering in layer 1
306 //            //  sigmasq[ico]+=rBP*rBP*geomfac[ico]*theta2BP/2; // multiple scattering in beam pipe
307               sigmasq[ico]+=TMath::Power(rad1*TMath::Tan(thetaL1),2)/3.;
308               sigmasq[ico]+=TMath::Power(rBP*TMath::Tan(thetaBP),2)/3.;
309
310 //            printf("Multipl. scatt. contr %d = %f (LAY1), %f (BP)\n",ico,rad1*rad1*geomfac[ico]*theta2L1/2,rBP*rBP*geomfac[ico]*theta2BP/2);
311 //            printf("Total error on coord %d = %f\n",ico,sigmasq[ico]);
312             }
313             Float_t wmat[9]={1.,0.,0.,0.,1.,0.,0.,0.,1.};
314             if(sigmasq[0]!=0.) wmat[0]=1./sigmasq[0];
315             if(sigmasq[1]!=0.) wmat[4]=1./sigmasq[1];
316             if(sigmasq[2]!=0.) wmat[8]=1./sigmasq[2];
317             new(fLines[nolines++])AliStrLine(gc1,sigmasq,wmat,gc2,kTRUE);
318
319           }
320           detTypeRec.ResetRecPoints();
321         }
322       }
323     }
324     prpl1.Clear();
325   }
326   if(nolines == 0)return -2;
327   return nolines;
328 }
329
330 //______________________________________________________________________
331 Int_t  AliITSVertexer3D::Prepare3DVertex(Int_t optCuts){
332   // Finds the 3D vertex information using tracklets
333   Int_t retcode = -1;
334
335   Float_t xbeam=GetNominalPos()[0];
336   Float_t ybeam=GetNominalPos()[1];
337   Float_t zvert=0.;
338   Float_t deltaR=fCoarseMaxRCut;
339   Float_t dZmax=fZCutDiamond;
340   if(optCuts){
341     xbeam=fVert3D.GetXv();
342     ybeam=fVert3D.GetYv();
343     zvert=fVert3D.GetZv();
344     deltaR=fMaxRCut;
345     dZmax=fMaxZCut;
346   }
347
348   Int_t nbr=50;
349   Float_t rl=-fCoarseMaxRCut;
350   Float_t rh=fCoarseMaxRCut;
351   Int_t nbz=100;
352   Float_t zl=-fZCutDiamond;
353   Float_t zh=fZCutDiamond;
354   Float_t binsizer=(rh-rl)/nbr;
355   Float_t binsizez=(zh-zl)/nbz;
356   TH3F *h3d = new TH3F("h3d","xyz distribution",nbr,rl,rh,nbr,rl,rh,nbz,zl,zh);
357   Int_t nbrcs=25;
358   Int_t nbzcs=50;
359   TH3F *h3dcs = new TH3F("h3dcs","xyz distribution",nbrcs,rl,rh,nbrcs,rl,rh,nbzcs,zl,zh);
360
361   // cleanup of the TCLonesArray of tracklets (i.e. fakes are removed)
362   Int_t *validate = new Int_t [fLines.GetEntriesFast()];
363   for(Int_t i=0; i<fLines.GetEntriesFast();i++)validate[i]=0;
364   for(Int_t i=0; i<fLines.GetEntriesFast()-1;i++){
365     AliStrLine *l1 = (AliStrLine*)fLines.At(i);
366     for(Int_t j=i+1;j<fLines.GetEntriesFast();j++){
367       AliStrLine *l2 = (AliStrLine*)fLines.At(j);
368       Float_t dca=l1->GetDCA(l2);
369       if(dca > fDCAcut || dca<0.00001) continue;
370       Double_t point[3];
371       Int_t retc = l1->Cross(l2,point);
372       if(retc<0)continue;
373       Double_t deltaZ=point[2]-zvert;
374      if(TMath::Abs(deltaZ)>dZmax)continue;
375       Double_t rad=TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
376       if(rad>fCoarseMaxRCut)continue;
377       Double_t deltaX=point[0]-xbeam;
378       Double_t deltaY=point[1]-ybeam;
379       Double_t raddist=TMath::Sqrt(deltaX*deltaX+deltaY*deltaY);
380       if(raddist>deltaR)continue;
381       validate[i]=1;
382       validate[j]=1;
383       h3d->Fill(point[0],point[1],point[2]);
384       h3dcs->Fill(point[0],point[1],point[2]);
385     }
386   }
387
388
389
390   Int_t numbtracklets=0;
391   for(Int_t i=0; i<fLines.GetEntriesFast();i++)if(validate[i]>=1)numbtracklets++;
392   if(numbtracklets<2){delete [] validate; delete h3d; return retcode; }
393
394   for(Int_t i=0; i<fLines.GetEntriesFast();i++){
395     if(validate[i]<1)fLines.RemoveAt(i);
396   }
397   fLines.Compress();
398   AliDebug(1,Form("Number of tracklets (after compress)%d ",fLines.GetEntriesFast()));
399   delete [] validate;
400
401   //        Find peaks in histos
402
403   Double_t peak[3]={0.,0.,0.};
404   Int_t ntrkl,ntimes;
405   FindPeaks(h3d,peak,ntrkl,ntimes);  
406   delete h3d;
407
408   if(optCuts==0 && ntrkl<=2){
409     ntrkl=0;
410     ntimes=0;
411     FindPeaks(h3dcs,peak,ntrkl,ntimes);  
412     binsizer=(rh-rl)/nbrcs;
413     binsizez=(zh-zl)/nbzcs;
414     if(ntrkl==1 || ntimes>1){delete h3dcs; return retcode;}
415   }
416   delete h3dcs;
417
418   //         Second selection loop
419
420   Float_t bs=(binsizer+binsizez)/2.;
421   for(Int_t i=0; i<fLines.GetEntriesFast();i++){
422     AliStrLine *l1 = (AliStrLine*)fLines.At(i);
423     if(l1->GetDistFromPoint(peak)>2.5*bs)fLines.RemoveAt(i);
424   }
425   fLines.Compress();
426   AliDebug(1,Form("Number of tracklets (after 2nd compression) %d",fLines.GetEntriesFast()));
427
428   if(fLines.GetEntriesFast()>1){
429     //  find a first candidate for the primary vertex
430     fVert3D=AliVertexerTracks::TrackletVertexFinder(&fLines,0); 
431     // make a further selection on tracklets based on this first candidate
432     fVert3D.GetXYZ(peak);
433     AliDebug(1,Form("FIRST V candidate: %f ; %f ; %f",peak[0],peak[1],peak[2]));
434     for(Int_t i=0; i<fLines.GetEntriesFast();i++){
435       AliStrLine *l1 = (AliStrLine*)fLines.At(i);
436       if(l1->GetDistFromPoint(peak)> fDCAcut)fLines.RemoveAt(i);
437     }
438     fLines.Compress();
439     AliDebug(1,Form("Number of tracklets (after 3rd compression) %d",fLines.GetEntriesFast()));
440     if(fLines.GetEntriesFast()>1) retcode=0; // this new tracklet selection is used
441     else retcode =1; // the previous tracklet selection will be used
442   }
443   else {
444     retcode = 0;
445   }
446   return retcode;  
447 }
448
449 //________________________________________________________
450 void AliITSVertexer3D::SetMeanPPtSelTracks(Float_t fieldTesla){
451   // Sets mean values of P and Pt based on the field
452   if(TMath::Abs(fieldTesla-0.5)<0.01){
453     SetMeanPSelTracks(0.885);
454     SetMeanPtSelTracks(0.630);
455   }else if(TMath::Abs(fieldTesla-0.4)<0.01){
456     SetMeanPSelTracks(0.805);
457     SetMeanPtSelTracks(0.580);
458   }else if(TMath::Abs(fieldTesla-0.2)<0.01){
459     SetMeanPSelTracks(0.740);
460     SetMeanPtSelTracks(0.530);
461   }else if(fieldTesla<0.00001){
462     SetMeanPSelTracks(0.730);
463     SetMeanPtSelTracks(0.510);
464   }else{
465     SetMeanPSelTracks();
466     SetMeanPtSelTracks();
467   }
468 }
469
470 //________________________________________________________
471 void AliITSVertexer3D::FindPeaks(TH3F* histo, Double_t *peak, Int_t &nOfTracklets, Int_t &nOfTimes){
472   // Finds bin with max contents in 3D histo of tracket intersections
473   TAxis *xax = histo->GetXaxis();  
474   TAxis *yax = histo->GetYaxis();
475   TAxis *zax = histo->GetZaxis();
476   peak[0]=0.;
477   peak[1]=0.;
478   peak[2]=0.;
479   nOfTracklets = 0;
480   nOfTimes=0;
481   for(Int_t i=xax->GetFirst();i<=xax->GetLast();i++){
482     Float_t xval = xax->GetBinCenter(i);
483     for(Int_t j=yax->GetFirst();j<=yax->GetLast();j++){
484       Float_t yval = yax->GetBinCenter(j);
485       for(Int_t k=zax->GetFirst();k<=zax->GetLast();k++){
486         Float_t zval = zax->GetBinCenter(k);
487         Int_t bc =(Int_t)histo->GetBinContent(i,j,k);
488         if(bc>nOfTracklets){
489           nOfTracklets = bc;
490           peak[2] = zval;
491           peak[1] = yval;
492           peak[0] = xval;
493           nOfTimes = 1;
494         }
495         if(bc==nOfTracklets){
496           nOfTimes++;
497         }
498       }
499     }
500   }
501   
502 }
503
504 //________________________________________________________
505 void AliITSVertexer3D::PrintStatus() const {
506   // Print current status
507   printf("=======================================================\n");
508   printf("Loose cut on Delta Phi %f\n",fCoarseDiffPhiCut);
509   printf("Cut on tracklet DCA to Z axis %f\n",fCoarseMaxRCut);
510   printf("Cut on tracklet DCA to beam axis %f\n",fMaxRCut);
511   printf("Cut on diamond (Z) %f\n",fZCutDiamond);
512   printf("Cut on DCA - tracklet to tracklet and to vertex %f\n",fDCAcut);
513   printf(" Max Phi difference: %f\n",fDiffPhiMax);
514   printf("=======================================================\n");
515 }