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