]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - ITS/AliITSVertexer3D.cxx
Memory leak fixed.
[u/mrichter/AliRoot.git] / ITS / AliITSVertexer3D.cxx
... / ...
CommitLineData
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
33ClassImp(AliITSVertexer3D)
34
35/* $Id$ */
36
37//______________________________________________________________________
38AliITSVertexer3D::AliITSVertexer3D():AliITSVertexer(),
39fLines("AliStrLine",1000),
40fVert3D(),
41fCoarseDiffPhiCut(0.),
42fCoarseMaxRCut(0.),
43fMaxRCut(0.),
44fZCutDiamond(0.),
45fMaxZCut(0.),
46fDCAcut(0.),
47fDiffPhiMax(0.),
48fMeanPSelTrk(0.),
49fMeanPtSelTrk(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//______________________________________________________________________
64AliITSVertexer3D::~AliITSVertexer3D() {
65 // Destructor
66 fLines.Clear("C");
67}
68
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}
78//______________________________________________________________________
79AliESDVertex* 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 vertz.SetLowLimit(-fZCutDiamond);
125 vertz.SetHighLimit(fZCutDiamond);
126 AliESDVertex* vtxz = vertz.FindVertexForCurrentEvent(itsClusterTree);
127 if(vtxz){
128 Double_t position[3]={GetNominalPos()[0],GetNominalPos()[1],vtxz->GetZv()};
129 Double_t covmatrix[6];
130 vtxz->GetCovMatrix(covmatrix);
131 Double_t chi2=99999.;
132 Int_t nContr=vtxz->GetNContributors();
133 fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);
134 fCurrentVertex->SetTitle("vertexer: Z");
135 fCurrentVertex->SetName("SPDVertexZ");
136 delete vtxz;
137 // printf("Vertexer Z success\n");
138 }
139
140 }
141 FindMultiplicity(itsClusterTree);
142 return fCurrentVertex;
143}
144
145//______________________________________________________________________
146Int_t AliITSVertexer3D::FindTracklets(TTree *itsClusterTree, Int_t optCuts){
147 // All the possible combinations between recpoints on layer 1and 2 are
148 // considered. Straight lines (=tracklets)are formed.
149 // The tracklets are processed in Prepare3DVertex
150 AliITSDetTypeRec detTypeRec;
151
152 TTree *tR = itsClusterTree;
153 detTypeRec.SetTreeAddressR(tR);
154 TClonesArray *itsRec = 0;
155 // lc1 and gc1 are local and global coordinates for layer 1
156 // Float_t lc1[3]={0.,0.,0.};
157 Float_t gc1[3]={0.,0.,0.};
158 // lc2 and gc2 are local and global coordinates for layer 2
159 // Float_t lc2[3]={0.,0.,0.};
160 Float_t gc2[3]={0.,0.,0.};
161
162 itsRec = detTypeRec.RecPoints();
163 TBranch *branch;
164 branch = tR->GetBranch("ITSRecPoints");
165
166 // Set values for cuts
167 Float_t xbeam=GetNominalPos()[0];
168 Float_t ybeam=GetNominalPos()[1];
169 Float_t zvert=0.;
170 Float_t deltaPhi=fCoarseDiffPhiCut;
171 Float_t deltaR=fCoarseMaxRCut;
172 Float_t dZmax=fZCutDiamond;
173 if(optCuts){
174 xbeam=fVert3D.GetXv();
175 ybeam=fVert3D.GetYv();
176 zvert=fVert3D.GetZv();
177 deltaPhi = fDiffPhiMax;
178 deltaR=fMaxRCut;
179 dZmax=fMaxZCut;
180 }
181 Int_t nrpL1 = 0; // number of rec points on layer 1
182 Int_t nrpL2 = 0; // number of rec points on layer 2
183
184 // By default irstL1=0 and lastL1=79
185 Int_t firstL1 = AliITSgeomTGeo::GetModuleIndex(1,1,1);
186 Int_t lastL1 = AliITSgeomTGeo::GetModuleIndex(2,1,1)-1;
187 for(Int_t module= firstL1; module<=lastL1;module++){ // count number of recopints on layer 1
188 branch->GetEvent(module);
189 nrpL1+= itsRec->GetEntries();
190 detTypeRec.ResetRecPoints();
191 }
192 //By default firstL2=80 and lastL2=239
193 Int_t firstL2 = AliITSgeomTGeo::GetModuleIndex(2,1,1);
194 Int_t lastL2 = AliITSgeomTGeo::GetModuleIndex(3,1,1)-1;
195 for(Int_t module= firstL2; module<=lastL2;module++){ // count number of recopints on layer 2
196 branch->GetEvent(module);
197 nrpL2+= itsRec->GetEntries();
198 detTypeRec.ResetRecPoints();
199 }
200 if(nrpL1 == 0 || nrpL2 == 0){
201 return -1;
202 }
203 AliDebug(1,Form("RecPoints on Layer 1,2 = %d, %d\n",nrpL1,nrpL2));
204
205 Double_t a[3]={xbeam,ybeam,0.};
206 Double_t b[3]={xbeam,ybeam,10.};
207 AliStrLine zeta(a,b,kTRUE);
208 Float_t bField=AliTracker::GetBz()/10.; //T
209 SetMeanPPtSelTracks(bField);
210
211 Int_t nolines = 0;
212 // Loop on modules of layer 1
213 for(Int_t modul1= firstL1; modul1<=lastL1;modul1++){ // Loop on modules of layer 1
214 if(!fUseModule[modul1]) continue;
215 UShort_t ladder=int(modul1/4)+1; // ladders are numbered starting from 1
216 branch->GetEvent(modul1);
217 Int_t nrecp1 = itsRec->GetEntries();
218 static TClonesArray prpl1("AliITSRecPoint",nrecp1);
219 prpl1.SetOwner();
220 for(Int_t j=0;j<nrecp1;j++){
221 AliITSRecPoint *recp = (AliITSRecPoint*)itsRec->At(j);
222 new(prpl1[j])AliITSRecPoint(*recp);
223 }
224 detTypeRec.ResetRecPoints();
225 for(Int_t j=0;j<nrecp1;j++){
226 AliITSRecPoint *recp1 = (AliITSRecPoint*)prpl1.At(j);
227 // Local coordinates of this recpoint
228 /*
229 lc[0]=recp1->GetDetLocalX();
230 lc[2]=recp1->GetDetLocalZ();
231 */
232 recp1->GetGlobalXYZ(gc1);
233 Double_t phi1 = TMath::ATan2(gc1[1]-ybeam,gc1[0]-xbeam);
234 if(phi1<0)phi1=2*TMath::Pi()+phi1;
235 for(Int_t ladl2=0 ; ladl2<fLadOnLay2*2+1;ladl2++){
236 for(Int_t k=0;k<4;k++){
237 Int_t ladmod=fLadders[ladder-1]+ladl2;
238 if(ladmod>AliITSgeomTGeo::GetNLadders(2)) ladmod=ladmod-AliITSgeomTGeo::GetNLadders(2);
239 Int_t modul2=AliITSgeomTGeo::GetModuleIndex(2,ladmod,k+1);
240 if(!fUseModule[modul2]) continue;
241 branch->GetEvent(modul2);
242 Int_t nrecp2 = itsRec->GetEntries();
243 for(Int_t j2=0;j2<nrecp2;j2++){
244 AliITSRecPoint *recp2 = (AliITSRecPoint*)itsRec->At(j2);
245 /*
246 lc2[0]=recp2->GetDetLocalX();
247 lc2[2]=recp2->GetDetLocalZ();
248 */
249 recp2->GetGlobalXYZ(gc2);
250 Double_t phi2 = TMath::ATan2(gc2[1]-ybeam,gc2[0]-xbeam);
251 if(phi2<0)phi2=2*TMath::Pi()+phi2;
252 Double_t diff = TMath::Abs(phi2-phi1);
253 if(diff>TMath::Pi())diff=2.*TMath::Pi()-diff;
254 if(diff>deltaPhi)continue;
255 AliStrLine line(gc1,gc2,kTRUE);
256 Double_t cp[3];
257 Int_t retcode = line.Cross(&zeta,cp);
258 if(retcode<0)continue;
259 Double_t dca = line.GetDCA(&zeta);
260 if(dca<0.) continue;
261 if(dca>deltaR)continue;
262 Double_t deltaZ=cp[2]-zvert;
263 if(TMath::Abs(deltaZ)>dZmax)continue;
264
265 if(nolines == 0){
266 if(fLines.GetEntriesFast()>0)fLines.Clear("C");
267 }
268 Float_t cov[6];
269 recp2->GetGlobalCov(cov);
270
271
272 Float_t rad1=TMath::Sqrt(gc1[0]*gc1[0]+gc1[1]*gc1[1]);
273 Float_t rad2=TMath::Sqrt(gc2[0]*gc2[0]+gc2[1]*gc2[1]);
274 Float_t factor=(rad1+rad2)/(rad2-rad1); //factor to account for error on tracklet direction
275
276 Float_t curvErr=0;
277 if(bField>0.00001){
278 Float_t curvRadius=fMeanPtSelTrk/(0.3*bField)*100; //cm
279 Float_t dRad=TMath::Sqrt(TMath::Power((gc1[0]-gc2[0]),2)+TMath::Power((gc1[1]-gc2[1]),2));
280 Float_t aux=dRad/2.+rad1;
281 curvErr=TMath::Sqrt(curvRadius*curvRadius-dRad*dRad/4.)-TMath::Sqrt(curvRadius*curvRadius-aux*aux); //cm
282 }
283
284 Float_t sigmasq[3];
285 sigmasq[0]=(cov[0]+curvErr*curvErr/2.)*factor*factor;
286 sigmasq[1]=(cov[3]+curvErr*curvErr/2.)*factor*factor;
287 sigmasq[2]=cov[5]*factor*factor;
288
289 // Multiple scattering
290 Float_t beta=1.;
291 Float_t beta2=beta*beta;
292 Float_t p2=fMeanPSelTrk*fMeanPSelTrk;
293 Float_t rBP=GetPipeRadius();
294 Float_t dBP=0.08/35.3; // 800 um of Be
295 Float_t dL1=0.01; //approx. 1% of radiation length
296 Float_t theta2BP=14.1*14.1/(beta2*p2*1e6)*TMath::Abs(dBP);
297 Float_t theta2L1=14.1*14.1/(beta2*p2*1e6)*TMath::Abs(dL1);
298 Float_t thetaBP=TMath::Sqrt(theta2BP);
299 Float_t thetaL1=TMath::Sqrt(theta2L1);
300// Float_t geomfac[3];
301// geomfac[0]=sin(phi1)*sin(phi1);
302// geomfac[1]=cos(phi1)*cos(phi1);
303// Float_t tgth=(gc2[2]-gc1[2])/(rad2-rad1);
304// geomfac[2]=1+tgth*tgth;
305 for(Int_t ico=0; ico<3;ico++){
306// printf("Error on coord. %d due to cov matrix+curvErr=%f\n",ico,sigmasq[ico]);
307// // sigmasq[ico]+=rad1*rad1*geomfac[ico]*theta2L1/2; // multiple scattering in layer 1
308// // sigmasq[ico]+=rBP*rBP*geomfac[ico]*theta2BP/2; // multiple scattering in beam pipe
309 sigmasq[ico]+=TMath::Power(rad1*TMath::Tan(thetaL1),2)/3.;
310 sigmasq[ico]+=TMath::Power(rBP*TMath::Tan(thetaBP),2)/3.;
311
312// printf("Multipl. scatt. contr %d = %f (LAY1), %f (BP)\n",ico,rad1*rad1*geomfac[ico]*theta2L1/2,rBP*rBP*geomfac[ico]*theta2BP/2);
313// printf("Total error on coord %d = %f\n",ico,sigmasq[ico]);
314 }
315 Float_t wmat[9]={1.,0.,0.,0.,1.,0.,0.,0.,1.};
316 if(sigmasq[0]!=0.) wmat[0]=1./sigmasq[0];
317 if(sigmasq[1]!=0.) wmat[4]=1./sigmasq[1];
318 if(sigmasq[2]!=0.) wmat[8]=1./sigmasq[2];
319 new(fLines[nolines++])AliStrLine(gc1,sigmasq,wmat,gc2,kTRUE);
320
321 }
322 detTypeRec.ResetRecPoints();
323 }
324 }
325 }
326 prpl1.Clear();
327 }
328 if(nolines == 0)return -2;
329 return nolines;
330}
331
332//______________________________________________________________________
333Int_t AliITSVertexer3D::Prepare3DVertex(Int_t optCuts){
334 // Finds the 3D vertex information using tracklets
335 Int_t retcode = -1;
336
337 Float_t xbeam=GetNominalPos()[0];
338 Float_t ybeam=GetNominalPos()[1];
339 Float_t zvert=0.;
340 Float_t deltaR=fCoarseMaxRCut;
341 Float_t dZmax=fZCutDiamond;
342 if(optCuts){
343 xbeam=fVert3D.GetXv();
344 ybeam=fVert3D.GetYv();
345 zvert=fVert3D.GetZv();
346 deltaR=fMaxRCut;
347 dZmax=fMaxZCut;
348 }
349
350 Int_t nbr=50;
351 Float_t rl=-fCoarseMaxRCut;
352 Float_t rh=fCoarseMaxRCut;
353 Int_t nbz=100;
354 Float_t zl=-fZCutDiamond;
355 Float_t zh=fZCutDiamond;
356 Float_t binsizer=(rh-rl)/nbr;
357 Float_t binsizez=(zh-zl)/nbz;
358 TH3F *h3d = new TH3F("h3d","xyz distribution",nbr,rl,rh,nbr,rl,rh,nbz,zl,zh);
359 Int_t nbrcs=25;
360 Int_t nbzcs=50;
361 TH3F *h3dcs = new TH3F("h3dcs","xyz distribution",nbrcs,rl,rh,nbrcs,rl,rh,nbzcs,zl,zh);
362
363 // cleanup of the TCLonesArray of tracklets (i.e. fakes are removed)
364 Int_t *validate = new Int_t [fLines.GetEntriesFast()];
365 for(Int_t i=0; i<fLines.GetEntriesFast();i++)validate[i]=0;
366 for(Int_t i=0; i<fLines.GetEntriesFast()-1;i++){
367 AliStrLine *l1 = (AliStrLine*)fLines.At(i);
368 for(Int_t j=i+1;j<fLines.GetEntriesFast();j++){
369 AliStrLine *l2 = (AliStrLine*)fLines.At(j);
370 Float_t dca=l1->GetDCA(l2);
371 if(dca > fDCAcut || dca<0.00001) continue;
372 Double_t point[3];
373 Int_t retc = l1->Cross(l2,point);
374 if(retc<0)continue;
375 Double_t deltaZ=point[2]-zvert;
376 if(TMath::Abs(deltaZ)>dZmax)continue;
377 Double_t rad=TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
378 if(rad>fCoarseMaxRCut)continue;
379 Double_t deltaX=point[0]-xbeam;
380 Double_t deltaY=point[1]-ybeam;
381 Double_t raddist=TMath::Sqrt(deltaX*deltaX+deltaY*deltaY);
382 if(raddist>deltaR)continue;
383 validate[i]=1;
384 validate[j]=1;
385 h3d->Fill(point[0],point[1],point[2]);
386 h3dcs->Fill(point[0],point[1],point[2]);
387 }
388 }
389
390
391
392 Int_t numbtracklets=0;
393 for(Int_t i=0; i<fLines.GetEntriesFast();i++)if(validate[i]>=1)numbtracklets++;
394 if(numbtracklets<2){delete [] validate; delete h3d; delete h3dcs; return retcode; }
395
396 for(Int_t i=0; i<fLines.GetEntriesFast();i++){
397 if(validate[i]<1)fLines.RemoveAt(i);
398 }
399 fLines.Compress();
400 AliDebug(1,Form("Number of tracklets (after compress)%d ",fLines.GetEntriesFast()));
401 delete [] validate;
402
403 // Find peaks in histos
404
405 Double_t peak[3]={0.,0.,0.};
406 Int_t ntrkl,ntimes;
407 FindPeaks(h3d,peak,ntrkl,ntimes);
408 delete h3d;
409
410 if(optCuts==0 && ntrkl<=2){
411 ntrkl=0;
412 ntimes=0;
413 FindPeaks(h3dcs,peak,ntrkl,ntimes);
414 binsizer=(rh-rl)/nbrcs;
415 binsizez=(zh-zl)/nbzcs;
416 if(ntrkl==1 || ntimes>1){delete h3dcs; return retcode;}
417 }
418 delete h3dcs;
419
420 // Second selection loop
421
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//________________________________________________________
452void 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//________________________________________________________
473void AliITSVertexer3D::FindPeaks(TH3F* histo, Double_t *peak, Int_t &nOfTracklets, Int_t &nOfTimes){
474 // Finds bin with max contents in 3D histo of tracket intersections
475 TAxis *xax = histo->GetXaxis();
476 TAxis *yax = histo->GetYaxis();
477 TAxis *zax = histo->GetZaxis();
478 peak[0]=0.;
479 peak[1]=0.;
480 peak[2]=0.;
481 nOfTracklets = 0;
482 nOfTimes=0;
483 for(Int_t i=xax->GetFirst();i<=xax->GetLast();i++){
484 Float_t xval = xax->GetBinCenter(i);
485 for(Int_t j=yax->GetFirst();j<=yax->GetLast();j++){
486 Float_t yval = yax->GetBinCenter(j);
487 for(Int_t k=zax->GetFirst();k<=zax->GetLast();k++){
488 Float_t zval = zax->GetBinCenter(k);
489 Int_t bc =(Int_t)histo->GetBinContent(i,j,k);
490 if(bc>nOfTracklets){
491 nOfTracklets = bc;
492 peak[2] = zval;
493 peak[1] = yval;
494 peak[0] = xval;
495 nOfTimes = 1;
496 }
497 if(bc==nOfTracklets){
498 nOfTimes++;
499 }
500 }
501 }
502 }
503
504}
505
506//________________________________________________________
507void AliITSVertexer3D::PrintStatus() const {
508 // Print current status
509 printf("=======================================================\n");
510 printf("Loose cut on Delta Phi %f\n",fCoarseDiffPhiCut);
511 printf("Cut on tracklet DCA to Z axis %f\n",fCoarseMaxRCut);
512 printf("Cut on tracklet DCA to beam axis %f\n",fMaxRCut);
513 printf("Cut on diamond (Z) %f\n",fZCutDiamond);
514 printf("Cut on DCA - tracklet to tracklet and to vertex %f\n",fDCAcut);
515 printf(" Max Phi difference: %f\n",fDiffPhiMax);
516 printf("=======================================================\n");
517}