]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSVertexer3D.cxx
bccbc54ce4701d1f45d9375ed0b855365cb5fcbc
[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 "AliRunLoader.h"
17 #include "AliESDVertex.h"
18 #include "AliLog.h"
19 #include "AliStrLine.h"
20 #include "AliTracker.h"
21 #include "AliITSDetTypeRec.h"
22 #include "AliITSRecPoint.h"
23 #include "AliITSRecPointContainer.h"
24 #include "AliITSgeomTGeo.h"
25 #include "AliVertexerTracks.h"
26 #include "AliITSVertexer3D.h"
27 #include "AliITSVertexerZ.h"
28 #include "AliITSSortTrkl.h"
29 /////////////////////////////////////////////////////////////////
30 // this class implements a method to determine
31 // the 3 coordinates of the primary vertex
32 // optimized for 
33 // p-p collisions
34 ////////////////////////////////////////////////////////////////
35
36 const Int_t    AliITSVertexer3D::fgkMaxNumOfClDefault = 500;
37
38 ClassImp(AliITSVertexer3D)
39
40 /* $Id$ */
41
42 //______________________________________________________________________
43 AliITSVertexer3D::AliITSVertexer3D():
44   AliITSVertexer(),
45   fLines("AliStrLine",1000),
46   fVert3D(),
47   fCoarseDiffPhiCut(0.),
48   fFineDiffPhiCut(0.),
49   fCutOnPairs(0.),
50   fCoarseMaxRCut(0.),
51   fMaxRCut(0.),
52   fMaxRCut2(0.),
53   fZCutDiamond(0.),
54   fMaxZCut(0.),
55   fDCAcut(0.),
56   fDiffPhiMax(0.),
57   fMeanPSelTrk(0.),
58   fMeanPtSelTrk(0.),
59   fUsedCluster(kMaxCluPerMod*kNSPDMod),
60   fZHisto(0),
61   fDCAforPileup(0.),
62   fDiffPhiforPileup(0.),
63   fBinSizeR(0.),
64   fBinSizeZ(0.),
65   fPileupAlgo(0),
66   fMaxNumOfCl(fgkMaxNumOfClDefault)
67 {
68   // Default constructor
69   SetCoarseDiffPhiCut();
70   SetFineDiffPhiCut();
71   SetCutOnPairs();
72   SetCoarseMaxRCut();
73   SetMaxRCut();
74   SetMaxRCutAlgo2();
75   SetZCutDiamond();
76   SetMaxZCut();
77   SetDCACut();
78   SetDiffPhiMax();
79   SetMeanPSelTracks();
80   SetMeanPtSelTracks();
81   SetMinDCAforPileup();
82   SetDeltaPhiforPileup();
83   SetPileupAlgo();
84   SetBinSizeR();
85   SetBinSizeZ();
86   Double_t binsize=0.02; // default 200 micron
87   Int_t nbins=static_cast<Int_t>(1+2*fZCutDiamond/binsize);
88   fZHisto=new TH1F("hz","",nbins,-fZCutDiamond,-fZCutDiamond+binsize*nbins);
89 }
90
91 //______________________________________________________________________
92 AliITSVertexer3D::~AliITSVertexer3D() {
93   // Destructor
94   fLines.Clear("C");
95   if(fZHisto) delete fZHisto;
96 }
97
98 //______________________________________________________________________
99 void AliITSVertexer3D::ResetVert3D(){
100   //
101   ResetVertex();
102   fVert3D.SetXv(0.);
103   fVert3D.SetYv(0.);
104   fVert3D.SetZv(0.);
105   fVert3D.SetDispersion(0.);
106   fVert3D.SetNContributors(0);
107   fUsedCluster.ResetAllBits(0);
108 }
109 //______________________________________________________________________
110 AliESDVertex* AliITSVertexer3D::FindVertexForCurrentEvent(TTree *itsClusterTree){
111   // Defines the AliESDVertex for the current event
112   ResetVert3D();
113   AliDebug(1,"FindVertexForCurrentEvent - 3D - PROCESSING NEXT EVENT");
114   fLines.Clear("C");
115   fCurrentVertex = NULL;
116
117   Int_t nolines = FindTracklets(itsClusterTree,0);
118   Int_t rc;
119   if(nolines>=2){
120     rc=Prepare3DVertex(0);
121     if(fVert3D.GetNContributors()>0){
122       fLines.Clear("C");
123       nolines = FindTracklets(itsClusterTree,1);
124       if(nolines>=2){
125         rc=Prepare3DVertex(1);
126         if(fPileupAlgo == 2 && rc == 0) FindVertex3DIterative();
127         else if(fPileupAlgo<2 && rc == 0) FindVertex3D(itsClusterTree);
128         if(rc!=0) fVert3D.SetNContributors(0); // exclude this vertex      
129       }
130     }
131   }
132
133   if(!fCurrentVertex){
134     AliITSVertexerZ vertz(GetNominalPos()[0],GetNominalPos()[1]);
135     vertz.SetDetTypeRec(GetDetTypeRec());
136     AliDebug(1,"Call Vertexer Z\n");
137     vertz.SetLowLimit(-fZCutDiamond);
138     vertz.SetHighLimit(fZCutDiamond);
139     AliESDVertex* vtxz = vertz.FindVertexForCurrentEvent(itsClusterTree);
140     if(vtxz){
141       Double_t position[3]={GetNominalPos()[0],GetNominalPos()[1],vtxz->GetZv()};
142       Double_t covmatrix[6];
143       vtxz->GetCovMatrix(covmatrix);
144       Double_t chi2=99999.;
145       Int_t    nContr=vtxz->GetNContributors();
146       fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);    
147       fCurrentVertex->SetDispersion(vtxz->GetDispersion());
148       fCurrentVertex->SetTitle("vertexer: Z");
149       fCurrentVertex->SetName("SPDVertexZ");
150       delete vtxz;
151     }
152
153   }
154   //RS
155   //  FindMultiplicity(itsClusterTree);
156   return fCurrentVertex;
157 }  
158
159 //______________________________________________________________________
160 void AliITSVertexer3D::FindVertex3D(TTree *itsClusterTree){
161  
162   Double_t vRadius=TMath::Sqrt(fVert3D.GetXv()*fVert3D.GetXv()+fVert3D.GetYv()*fVert3D.GetYv());
163   if(vRadius<GetPipeRadius() && fVert3D.GetNContributors()>0){
164     Double_t position[3]={fVert3D.GetXv(),fVert3D.GetYv(),fVert3D.GetZv()};
165     Double_t covmatrix[6];
166     fVert3D.GetCovMatrix(covmatrix);
167     Double_t chi2=99999.;
168     Int_t    nContr=fVert3D.GetNContributors();
169     fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);    
170     fCurrentVertex->SetTitle("vertexer: 3D");
171     fCurrentVertex->SetName("SPDVertex3D");
172     fCurrentVertex->SetDispersion(fVert3D.GetDispersion());
173     fNoVertices=1;
174     
175     switch(fPileupAlgo){
176     case 0: PileupFromZ(); break;
177     case 1: FindOther3DVertices(itsClusterTree); break;
178     default: AliError("Wrong pileup algorithm"); break;
179     }
180     if(fNoVertices==1){
181       fVertArray = new AliESDVertex[1];
182       fVertArray[0]=(*fCurrentVertex);    
183     }
184   }
185 }
186
187 //______________________________________________________________________
188 void AliITSVertexer3D::FindVertex3DIterative(){
189   //
190
191   Int_t nLines=fLines.GetEntriesFast();
192   Int_t maxPoints=nLines*(nLines-1)/2;
193   Double_t* xP=new Double_t[maxPoints];
194   Double_t* yP=new Double_t[maxPoints];
195   Double_t* zP=new Double_t[maxPoints];
196   Int_t* index1=new Int_t[maxPoints];
197   Int_t* index2=new Int_t[maxPoints];
198   Double_t xbeam=fVert3D.GetXv();
199   Double_t ybeam=fVert3D.GetYv();
200
201   Int_t iPoint=0;
202   for(Int_t ilin1=0; ilin1<nLines; ilin1++){
203     AliStrLine *l1 = (AliStrLine*)fLines.At(ilin1);
204     for(Int_t ilin2=ilin1+1; ilin2<nLines; ilin2++){
205       AliStrLine *l2 = (AliStrLine*)fLines.At(ilin2);
206       Double_t dca=l1->GetDCA(l2);
207       if(dca > fDCAcut || dca<0.00001) continue;
208       Double_t point[3];
209       Int_t retc = l1->Cross(l2,point);
210       if(retc<0)continue;
211       Double_t rad=TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
212       if(rad>fCoarseMaxRCut)continue;
213       Double_t distFromBeam=TMath::Sqrt((point[0]-xbeam)*(point[0]-xbeam)+(point[1]-ybeam)*(point[1]-ybeam));
214       if(distFromBeam>fMaxRCut2) continue;
215       xP[iPoint]=point[0];
216       yP[iPoint]=point[1];
217       zP[iPoint]=point[2];
218       index1[iPoint]=ilin1;
219       index2[iPoint]=ilin2;
220       iPoint++;
221     }
222   }
223   Int_t npoints=iPoint++;
224   Int_t index=0;
225   Short_t* mask=new Short_t[npoints];
226   for(Int_t ip=0;ip<npoints;ip++) mask[ip]=-1;
227  
228   for(Int_t ip1=0;ip1<npoints;ip1++){
229     if(mask[ip1]==-1) mask[ip1]=index++;
230     for(Int_t ip2=ip1+1; ip2<npoints; ip2++){
231       if(mask[ip2]==mask[ip1] && mask[ip2]!=-1) continue;
232       Double_t dist2=(xP[ip1]-xP[ip2])*(xP[ip1]-xP[ip2]);
233       dist2+=(yP[ip1]-yP[ip2])*(yP[ip1]-yP[ip2]);
234       dist2+=(zP[ip1]-zP[ip2])*(zP[ip1]-zP[ip2]);
235       if(dist2<fCutOnPairs*fCutOnPairs){ 
236         if(mask[ip2]==-1) mask[ip2]=mask[ip1];
237         else{
238           for(Int_t ip=0; ip<npoints;ip++){
239             if(mask[ip]==mask[ip2]) mask[ip]=mask[ip1];
240           }
241         }
242       }
243     }
244   }
245
246
247   // Count multiplicity of trackelts in clusters
248   UInt_t* isIndUsed=new UInt_t[index+1];
249   for(Int_t ind=0;ind<index+1;ind++) isIndUsed[ind]=0;
250   for(Int_t ip=0; ip<npoints;ip++){
251     Int_t ind=mask[ip];
252     isIndUsed[ind]++;
253   }
254
255   // Count clusters/vertices and sort according to multiplicity
256   Int_t nClusters=0;
257   Int_t* sortedIndex=new Int_t[index+1];
258   for(Int_t ind1=0;ind1<index+1;ind1++){
259     if(isIndUsed[ind1]<=1) isIndUsed[ind1]=0;
260     else nClusters++;
261     UInt_t cap=9999999;
262     if(ind1>0) cap=isIndUsed[sortedIndex[ind1-1]];
263     UInt_t bigger=0;
264     Int_t biggerindex=-1;
265     for(Int_t ind2=0;ind2<index+1;ind2++){
266       Bool_t use=kTRUE;
267       for(Int_t ind3=0; ind3<ind1; ind3++)
268         if(ind2==sortedIndex[ind3]) use=kFALSE;
269       if(use && isIndUsed[ind2]>bigger && isIndUsed[ind2]<=cap){
270         bigger=isIndUsed[ind2];
271         biggerindex=ind2;
272       }
273     }
274     sortedIndex[ind1]=biggerindex;    
275   }
276   AliDebug(3,Form("Number of clusters before merging = %d\n",nClusters));
277
278   // Assign lines to clusters/vertices and merge clusters which share 1 line
279   Int_t nClustersAfterMerge=nClusters;
280   Int_t* belongsTo=new Int_t[nLines];
281   for(Int_t ilin=0; ilin<nLines; ilin++) belongsTo[ilin]=-1;
282   for(Int_t iclu=0;iclu<nClusters;iclu++){
283     Int_t actualCluIndex=iclu;
284     for(Int_t ip=0; ip<npoints;ip++){
285       if(mask[ip]==sortedIndex[iclu]){
286         Int_t ind1=index1[ip];
287         if(belongsTo[ind1]==-1) belongsTo[ind1]=actualCluIndex;
288         else if(belongsTo[ind1]<actualCluIndex){
289           Int_t newCluIndex=belongsTo[ind1];
290           for(Int_t ilin=0; ilin<nLines; ilin++){
291             if(belongsTo[ilin]==actualCluIndex) belongsTo[ilin]=newCluIndex;
292           }
293           AliDebug(10,Form("Merged cluster %d with %d\n",actualCluIndex,newCluIndex));
294           actualCluIndex=newCluIndex;
295           nClustersAfterMerge--;
296         }
297         Int_t ind2=index2[ip];      
298         if(belongsTo[ind2]==-1) belongsTo[ind2]=actualCluIndex;
299         else if(belongsTo[ind2]<actualCluIndex){
300           Int_t newCluIndex=belongsTo[ind2];
301           for(Int_t ilin=0; ilin<nLines; ilin++){
302             if(belongsTo[ilin]==actualCluIndex) belongsTo[ilin]=newCluIndex;
303           }
304           AliDebug(10,Form("Merged cluster %d with %d\n",actualCluIndex,newCluIndex));
305           actualCluIndex=newCluIndex;
306           nClustersAfterMerge--;
307         }
308       }
309     }
310   }
311   AliDebug(3,Form("Number of clusters after merging = %d\n",nClustersAfterMerge));
312   
313   // Count lines associated to each cluster/vertex
314   UInt_t *cluSize=new UInt_t[nClusters];
315   for(Int_t iclu=0;iclu<nClusters;iclu++){ 
316     cluSize[iclu]=0;
317     for(Int_t ilin=0; ilin<nLines; ilin++){
318       if(belongsTo[ilin]==iclu) cluSize[iclu]++;
319     }
320   }
321
322   // Count good vertices (>1 associated tracklet)
323   UInt_t nGoodVert=0;
324   for(Int_t iclu=0;iclu<nClusters;iclu++){ 
325     AliDebug(3,Form("Vertex %d Size=%d\n",iclu,cluSize[iclu]));
326     if(cluSize[iclu]>1) nGoodVert++;
327   }
328     
329   AliDebug(1,Form("Number of good vertices = %d\n",nGoodVert));
330   // Calculate vertex coordinates for each cluster
331   if(nGoodVert>0){
332     fVertArray = new AliESDVertex[nGoodVert];
333     Int_t iVert=0;
334     for(Int_t iclu=0;iclu<nClusters;iclu++){
335       Int_t size=cluSize[iclu];
336       if(size>1){
337         AliStrLine **arrlin = new AliStrLine*[size];
338         Int_t nFilled=0;
339         for(Int_t ilin=0; ilin<nLines; ilin++){
340           if(belongsTo[ilin]==iclu){
341             arrlin[nFilled++] = dynamic_cast<AliStrLine*>(fLines[ilin]);
342           }
343         }      
344         AliDebug(3,Form("Vertex %d  N associated tracklets = %d out of %d\n",iVert,size,nFilled));
345
346         fVertArray[iVert]=AliVertexerTracks::TrackletVertexFinder(arrlin,nFilled);
347         Double_t peak[3];
348         fVertArray[iVert].GetXYZ(peak);
349         AliStrLine **arrlin2 = new AliStrLine*[size];
350         Int_t nFilled2=0;       
351         for(Int_t i=0; i<nFilled;i++){
352           AliStrLine *l1 = arrlin[i];     
353           if(l1->GetDistFromPoint(peak)< fDCAcut)
354             arrlin2[nFilled2++] = dynamic_cast<AliStrLine*>(l1);
355         }
356         if(nFilled2>1){
357           AliDebug(3,Form("Vertex %d  recalculated with %d tracklets\n",iVert,nFilled2));
358           fVertArray[iVert]=AliVertexerTracks::TrackletVertexFinder(arrlin2,nFilled2);
359         }
360         delete [] arrlin;
361         delete [] arrlin2;
362         ++iVert;
363       }
364     }
365     
366     if(nGoodVert > 1){
367       fIsPileup = kTRUE;
368       fNTrpuv = fVertArray[1].GetNContributors();
369       fZpuv = fVertArray[1].GetZv();
370     }
371     
372     Double_t vRadius=TMath::Sqrt(fVertArray[0].GetXv()*fVertArray[0].GetXv()+fVertArray[0].GetYv()*fVertArray[0].GetYv());
373     if(vRadius<GetPipeRadius() && fVertArray[0].GetNContributors()>0){
374       Double_t position[3]={fVertArray[0].GetXv(),fVertArray[0].GetYv(),fVertArray[0].GetZv()};
375       Double_t covmatrix[6];
376       fVertArray[0].GetCovMatrix(covmatrix);
377       Double_t chi2=99999.;
378       Int_t    nContr=fVertArray[0].GetNContributors();
379       fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);    
380       fCurrentVertex->SetTitle("vertexer: 3D");
381       fCurrentVertex->SetName("SPDVertex3D");
382       fCurrentVertex->SetDispersion(fVertArray[0].GetDispersion());  
383     }
384   }
385
386   delete [] index1;
387   delete [] index2;
388   delete [] mask;
389   delete [] isIndUsed;
390   delete [] sortedIndex;
391   delete [] belongsTo;
392   delete [] cluSize;
393   delete [] xP;
394   delete [] yP;
395   delete [] zP;
396 }
397 //______________________________________________________________________
398 void AliITSVertexer3D::FindVertex3DIterativeMM(){
399   // Defines the AliESDVertex for the current event
400   Int_t numsor=fLines.GetEntriesFast()*(fLines.GetEntriesFast()-1)/2;
401   //cout<<"AliITSVertexer3D::FindVertexForCurentEvent: Number of tracklets selected for vertexing "<<fLines.GetEntriesFast()<<"; Number of pairs: "<<numsor<<endl;
402   AliITSSortTrkl srt(fLines,numsor,fCutOnPairs,fCoarseMaxRCut); 
403   srt.FindClusters();   
404   AliInfo(Form("Number of vertices: %d",srt.GetNumberOfClusters()));
405       
406   fNoVertices = srt.GetNumberOfClusters();
407   //printf("fNoVertices = %d \n",fNoVertices);
408   if(fNoVertices>0){
409     fVertArray = new AliESDVertex[fNoVertices];
410     for(Int_t kk=0; kk<srt.GetNumberOfClusters(); kk++){
411       Int_t size = 0;
412       Int_t *labels = srt.GetTrackletsLab(kk,size);
413       /*
414         Int_t *pairs = srt.GetClusters(kk);
415         Int_t nopai = srt.GetSizeOfCluster(kk);
416         cout<<"***** Vertex number "<<kk<<".  Pairs: \n";
417         for(Int_t jj=0;jj<nopai;jj++){
418         cout<<pairs[jj]<<" - ";
419         if(jj>0 & jj%8==0)cout<<endl;
420         }
421         cout<<endl;
422         cout<<"***** Vertex number "<<kk<<".  Labels: \n";
423       */
424       AliStrLine **tclo = new AliStrLine* [size];
425       for(Int_t jj=0;jj<size;jj++){
426         //          cout<<labels[jj]<<" - ";
427         //          if(jj>0 & jj%8==0)cout<<endl;
428         tclo[jj] = dynamic_cast<AliStrLine*>(fLines[labels[jj]]);
429       }
430       //          cout<<endl;
431       delete []labels;
432       fVertArray[kk]=AliVertexerTracks::TrackletVertexFinder(tclo,size);
433       delete [] tclo;
434       //          fVertArray[kk].PrintStatus();
435       if(kk == 1){
436         // at least one second vertex is present
437         fIsPileup = kTRUE;
438         fNTrpuv = fVertArray[kk].GetNContributors();
439         fZpuv = fVertArray[kk].GetZv();
440       }
441     }
442     Double_t vRadius=TMath::Sqrt(fVertArray[0].GetXv()*fVertArray[0].GetXv()+fVertArray[0].GetYv()*fVertArray[0].GetYv());
443     if(vRadius<GetPipeRadius() && fVertArray[0].GetNContributors()>0){
444       Double_t position[3]={fVertArray[0].GetXv(),fVertArray[0].GetYv(),fVertArray[0].GetZv()};
445       Double_t covmatrix[6];
446       fVertArray[0].GetCovMatrix(covmatrix);
447       Double_t chi2=99999.;
448       Int_t    nContr=fVertArray[0].GetNContributors();
449       fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);    
450       fCurrentVertex->SetTitle("vertexer: 3D");
451       fCurrentVertex->SetName("SPDVertex3D");
452       fCurrentVertex->SetDispersion(fVertArray[0].GetDispersion());  
453     }
454   }
455
456 }  
457
458 //______________________________________________________________________
459 Bool_t AliITSVertexer3D::DistBetweenVertices(AliESDVertex &a, AliESDVertex &b, Double_t test, Double_t &dist){
460   // method to compare the distance between vertices a and b with "test"
461   //it returns kTRUE is the distance is less or equal to test
462   dist = (a.GetX()-b.GetX()) * (a.GetX()-b.GetX());
463   dist +=  (a.GetY()-b.GetY()) * (a.GetY()-b.GetY());
464   dist +=  (a.GetZ()-b.GetZ()) * (a.GetZ()-b.GetZ());
465   dist = TMath::Sqrt(dist);
466   if(dist <= test)return kTRUE;
467   return kFALSE;
468 }
469
470
471 //______________________________________________________________________
472 Int_t AliITSVertexer3D::FindTracklets(TTree *itsClusterTree, Int_t optCuts){
473   // All the possible combinations between recpoints on layer 1and 2 are
474   // considered. Straight lines (=tracklets)are formed. 
475   // The tracklets are processed in Prepare3DVertex
476
477   TClonesArray *itsRec  = 0;
478   if(optCuts==0) fZHisto->Reset();
479  // gc1 are local and global coordinates for layer 1
480   Float_t gc1f[3]={0.,0.,0.};
481   Double_t gc1[3]={0.,0.,0.};
482   // gc2 are local and global coordinates for layer 2
483   Float_t gc2f[3]={0.,0.,0.};
484   Double_t gc2[3]={0.,0.,0.};
485   AliITSRecPointContainer* rpcont=AliITSRecPointContainer::Instance();
486   itsRec=rpcont->FetchClusters(0,itsClusterTree);
487   if(!rpcont->IsSPDActive()){
488     AliWarning("No SPD rec points found, 3D vertex not calculated");
489     return -1;
490   }
491
492   // Set values for cuts
493   Double_t xbeam=GetNominalPos()[0]; 
494   Double_t ybeam=GetNominalPos()[1];
495   Double_t zvert=0.;
496   Double_t deltaPhi=fCoarseDiffPhiCut;
497   Double_t deltaR=fCoarseMaxRCut;
498   Double_t dZmax=fZCutDiamond;
499   if(optCuts==1){
500     xbeam=fVert3D.GetXv();
501     ybeam=fVert3D.GetYv();
502     zvert=fVert3D.GetZv();
503     deltaPhi = fDiffPhiMax; 
504     deltaR=fMaxRCut;
505     dZmax=fMaxZCut;
506     if(fPileupAlgo == 2){
507       dZmax=fZCutDiamond;
508       deltaR=fMaxRCut2;
509     }
510   } else if(optCuts==2){
511     xbeam=fVert3D.GetXv();
512     ybeam=fVert3D.GetYv();
513     deltaPhi = fDiffPhiforPileup;
514     deltaR=fMaxRCut;
515   }
516
517   Int_t nrpL1 = 0;    // number of rec points on layer 1
518   Int_t nrpL2 = 0;    // number of rec points on layer 2
519   nrpL1=rpcont->GetNClustersInLayerFast(1);
520   nrpL2=rpcont->GetNClustersInLayerFast(2);
521   if(nrpL1 == 0 || nrpL2 == 0){
522     AliDebug(1,Form("No RecPoints in at least one SPD layer (%d %d)",nrpL1,nrpL2));
523     return -1;
524   }
525   AliDebug(1,Form("RecPoints on Layer 1,2 = %d, %d\n",nrpL1,nrpL2));
526   if(nrpL1>fMaxNumOfCl || nrpL2>fMaxNumOfCl){
527     AliWarning(Form("Too many recpoints on SPD(%d %d ), call vertexerZ",nrpL1,nrpL2));
528     return -1;
529   }
530
531   Double_t a[3]={xbeam,ybeam,0.}; 
532   Double_t b[3]={xbeam,ybeam,10.};
533   AliStrLine zeta(a,b,kTRUE);
534   static Double_t bField=TMath::Abs(AliTracker::GetBz()/10.); //T
535   SetMeanPPtSelTracks(bField);
536
537   Int_t nolines = 0;
538   // Loop on modules of layer 1
539   Int_t firstL1 = AliITSgeomTGeo::GetModuleIndex(1,1,1);
540   Int_t lastL1 = AliITSgeomTGeo::GetModuleIndex(2,1,1)-1;
541   for(Int_t modul1= firstL1; modul1<=lastL1;modul1++){   // Loop on modules of layer 1
542     if(!fUseModule[modul1]) continue;
543     UShort_t ladder=modul1/4+1; // ladders are numbered starting from 1
544     TClonesArray *prpl1=rpcont->UncheckedGetClusters(modul1);
545     Int_t nrecp1 = prpl1->GetEntries();
546     for(Int_t j=0;j<nrecp1;j++){
547       if(j>kMaxCluPerMod) continue;
548       UShort_t idClu1=modul1*kMaxCluPerMod+j;
549       if(fUsedCluster.TestBitNumber(idClu1)) continue;
550       AliITSRecPoint *recp1 = (AliITSRecPoint*)prpl1->At(j);
551       recp1->GetGlobalXYZ(gc1f);
552       for(Int_t ico=0;ico<3;ico++)gc1[ico]=gc1f[ico];
553
554       Double_t phi1 = TMath::ATan2(gc1[1]-ybeam,gc1[0]-xbeam);
555       if(phi1<0)phi1=2*TMath::Pi()+phi1;
556       for(Int_t ladl2=0 ; ladl2<fLadOnLay2*2+1;ladl2++){
557         for(Int_t k=0;k<4;k++){
558           Int_t ladmod=fLadders[ladder-1]+ladl2;
559           if(ladmod>AliITSgeomTGeo::GetNLadders(2)) ladmod=ladmod-AliITSgeomTGeo::GetNLadders(2);
560           Int_t modul2=AliITSgeomTGeo::GetModuleIndex(2,ladmod,k+1);
561           if(!fUseModule[modul2]) continue;
562           itsRec=rpcont->UncheckedGetClusters(modul2);
563           Int_t nrecp2 = itsRec->GetEntries();
564           for(Int_t j2=0;j2<nrecp2;j2++){
565             if(j2>kMaxCluPerMod) continue;
566             UShort_t idClu2=modul2*kMaxCluPerMod+j2;
567             if(fUsedCluster.TestBitNumber(idClu2)) continue;
568             AliITSRecPoint *recp2 = (AliITSRecPoint*)itsRec->At(j2);
569             recp2->GetGlobalXYZ(gc2f);
570             for(Int_t ico=0;ico<3;ico++)gc2[ico]=gc2f[ico];
571             Double_t phi2 = TMath::ATan2(gc2[1]-ybeam,gc2[0]-xbeam);
572             if(phi2<0)phi2=2*TMath::Pi()+phi2;
573             Double_t diff = TMath::Abs(phi2-phi1); 
574             if(diff>TMath::Pi())diff=2.*TMath::Pi()-diff; 
575             if(optCuts==0 && diff<fDiffPhiforPileup){
576               Double_t r1=TMath::Sqrt(gc1[0]*gc1[0]+gc1[1]*gc1[1]);
577               Double_t zc1=gc1[2];
578               Double_t r2=TMath::Sqrt(gc2[0]*gc2[0]+gc2[1]*gc2[1]);
579               Double_t zc2=gc2[2];
580               Double_t zr0=(r2*zc1-r1*zc2)/(r2-r1); //Z @ null radius
581               fZHisto->Fill(zr0);
582             }
583             if(diff>deltaPhi)continue;
584             AliStrLine line(gc1,gc2,kTRUE);
585             Double_t cp[3];
586             Int_t retcode = line.Cross(&zeta,cp);
587             if(retcode<0)continue;
588             Double_t dca = line.GetDCA(&zeta);
589             if(dca<0.) continue;
590             if(dca>deltaR)continue;
591             Double_t deltaZ=cp[2]-zvert;
592             if(TMath::Abs(deltaZ)>dZmax)continue;
593
594             if(nolines == 0){
595               if(fLines.GetEntriesFast()>0)fLines.Clear("C");
596             }
597             Float_t cov[6];
598             recp2->GetGlobalCov(cov);
599
600
601             Double_t rad1=TMath::Sqrt(gc1[0]*gc1[0]+gc1[1]*gc1[1]);
602             Double_t rad2=TMath::Sqrt(gc2[0]*gc2[0]+gc2[1]*gc2[1]);
603             Double_t factor=(rad1+rad2)/(rad2-rad1); //factor to account for error on tracklet direction 
604
605             Double_t curvErr=0;
606             if(bField>0.00001){
607               Double_t curvRadius=fMeanPtSelTrk/(0.3*bField)*100; //cm 
608               Double_t dRad=TMath::Sqrt((gc1[0]-gc2[0])*(gc1[0]-gc2[0])+(gc1[1]-gc2[1])*(gc1[1]-gc2[1]));
609               Double_t aux=dRad/2.+rad1;
610               curvErr=TMath::Sqrt(curvRadius*curvRadius-dRad*dRad/4.)-TMath::Sqrt(curvRadius*curvRadius-aux*aux); //cm
611             }
612             Double_t sigmasq[3];
613             sigmasq[0]=(cov[0]+curvErr*curvErr/2.)*factor*factor;
614             sigmasq[1]=(cov[3]+curvErr*curvErr/2.)*factor*factor;
615             sigmasq[2]=cov[5]*factor*factor;
616
617             // Multiple scattering
618             Double_t pOverMass=fMeanPSelTrk/0.140;
619             Double_t beta2=pOverMass*pOverMass/(1+pOverMass*pOverMass);
620             Double_t p2=fMeanPSelTrk*fMeanPSelTrk;
621             Double_t rBP=GetPipeRadius();
622             Double_t dBP=0.08/35.3; // 800 um of Be
623             Double_t dL1=0.01; //approx. 1% of radiation length  
624             Double_t theta2BP=14.1*14.1/(beta2*p2*1e6)*dBP;
625             Double_t theta2L1=14.1*14.1/(beta2*p2*1e6)*dL1;
626             Double_t rtantheta1=(rad2-rad1)*TMath::Tan(TMath::Sqrt(theta2L1));
627             Double_t rtanthetaBP=(rad1-rBP)*TMath::Tan(TMath::Sqrt(theta2BP));
628             for(Int_t ico=0; ico<3;ico++){    
629               sigmasq[ico]+=rtantheta1*rtantheta1*factor*factor/3.;
630               sigmasq[ico]+=rtanthetaBP*rtanthetaBP*factor*factor/3.;
631             }
632             Double_t wmat[9]={1.,0.,0.,0.,1.,0.,0.,0.,1.};
633             if(sigmasq[0]!=0.) wmat[0]=1./sigmasq[0];
634             if(sigmasq[1]!=0.) wmat[4]=1./sigmasq[1];
635             if(sigmasq[2]!=0.) wmat[8]=1./sigmasq[2];
636             new(fLines[nolines++])AliStrLine(gc1,sigmasq,wmat,gc2,kTRUE,idClu1,idClu2);
637
638           }
639         }
640       }
641     }
642   }
643   if(nolines == 0)return -2;
644   return nolines;
645 }
646
647 //______________________________________________________________________
648 Int_t  AliITSVertexer3D::Prepare3DVertex(Int_t optCuts){
649   // Finds the 3D vertex information using tracklets
650   Int_t retcode = -1;
651
652   Double_t xbeam=GetNominalPos()[0];
653   Double_t ybeam=GetNominalPos()[1];
654   Double_t zvert=0.;
655   Double_t deltaR=fCoarseMaxRCut;
656   Double_t dZmax=fZCutDiamond;
657   if(optCuts==1){
658     xbeam=fVert3D.GetXv();
659     ybeam=fVert3D.GetYv();
660     zvert=fVert3D.GetZv();
661     deltaR=fMaxRCut;
662     dZmax=fMaxZCut;
663     if(fPileupAlgo == 2){ 
664       dZmax=fZCutDiamond;
665       deltaR=fMaxRCut2;
666     }
667   }else if(optCuts==2){
668     xbeam=fVert3D.GetXv();
669     ybeam=fVert3D.GetYv();
670     deltaR=fMaxRCut;
671   }
672
673   Double_t rl=-fCoarseMaxRCut;
674   Double_t rh=fCoarseMaxRCut;
675   Double_t zl=-fZCutDiamond;
676   Double_t zh=fZCutDiamond;
677   Int_t nbr=(Int_t)((rh-rl)/fBinSizeR+0.0001);
678   Int_t nbz=(Int_t)((zh-zl)/fBinSizeZ+0.0001);
679   Int_t nbrcs=(Int_t)((rh-rl)/(fBinSizeR*2.)+0.0001);
680   Int_t nbzcs=(Int_t)((zh-zl)/(fBinSizeZ*2.)+0.0001);
681
682   TH3F *h3d = new TH3F("h3d","xyz distribution",nbr,rl,rh,nbr,rl,rh,nbz,zl,zh);
683   TH3F *h3dcs = new TH3F("h3dcs","xyz distribution",nbrcs,rl,rh,nbrcs,rl,rh,nbzcs,zl,zh);
684
685   // cleanup of the TCLonesArray of tracklets (i.e. fakes are removed)
686   Int_t *validate = new Int_t [fLines.GetEntriesFast()];
687   for(Int_t i=0; i<fLines.GetEntriesFast();i++)validate[i]=0;
688   for(Int_t i=0; i<fLines.GetEntriesFast()-1;i++){
689     AliStrLine *l1 = (AliStrLine*)fLines.At(i);
690     for(Int_t j=i+1;j<fLines.GetEntriesFast();j++){
691       AliStrLine *l2 = (AliStrLine*)fLines.At(j);
692       Double_t dca=l1->GetDCA(l2);
693       if(dca > fDCAcut || dca<0.00001) continue;
694       Double_t point[3];
695       Int_t retc = l1->Cross(l2,point);
696       if(retc<0)continue;
697       Double_t deltaZ=point[2]-zvert;
698      if(TMath::Abs(deltaZ)>dZmax)continue;
699       Double_t rad=TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
700       if(rad>fCoarseMaxRCut)continue;
701       Double_t deltaX=point[0]-xbeam;
702       Double_t deltaY=point[1]-ybeam;
703       Double_t raddist=TMath::Sqrt(deltaX*deltaX+deltaY*deltaY);
704       if(raddist>deltaR)continue;
705       validate[i]=1;
706       validate[j]=1;
707       h3d->Fill(point[0],point[1],point[2]);
708       h3dcs->Fill(point[0],point[1],point[2]);
709     }
710   }
711
712
713
714   Int_t numbtracklets=0;
715   for(Int_t i=0; i<fLines.GetEntriesFast();i++)if(validate[i]>=1)numbtracklets++;
716   if(numbtracklets<2){
717     delete [] validate; 
718     delete h3d; 
719     delete h3dcs; 
720     return retcode; 
721   }
722
723   for(Int_t i=0; i<fLines.GetEntriesFast();i++){
724     if(validate[i]<1)fLines.RemoveAt(i);
725   }
726   fLines.Compress();
727   AliDebug(1,Form("Number of tracklets (after compress)%d ",fLines.GetEntriesFast()));
728   delete [] validate;
729
730   // Exit here if Pileup Algorithm 2 has been chosen during second loop
731   if(fPileupAlgo == 2 && optCuts==1){
732     delete h3d; 
733     delete h3dcs;     
734     return 0;
735   }
736
737   //        Find peaks in histos
738
739   Double_t peak[3]={0.,0.,0.};
740   Int_t ntrkl,ntimes;
741   FindPeaks(h3d,peak,ntrkl,ntimes);  
742   delete h3d;
743   Double_t binsizer=(rh-rl)/nbr;
744   Double_t binsizez=(zh-zl)/nbz;
745   if(optCuts==0 && (ntrkl<=2 || ntimes>1)){
746     ntrkl=0;
747     ntimes=0;
748     FindPeaks(h3dcs,peak,ntrkl,ntimes);  
749     binsizer=(rh-rl)/nbrcs;
750     binsizez=(zh-zl)/nbzcs;
751     if(ntrkl==1 || ntimes>1){delete h3dcs; return retcode;}
752   }
753   delete h3dcs;
754
755   //         Second selection loop
756
757   Double_t bs=(binsizer+binsizez)/2.;
758   for(Int_t i=0; i<fLines.GetEntriesFast();i++){
759     AliStrLine *l1 = (AliStrLine*)fLines.At(i);
760     if(l1->GetDistFromPoint(peak)>2.5*bs)fLines.RemoveAt(i);
761   }
762   fLines.Compress();
763   AliDebug(1,Form("Number of tracklets (after 2nd compression) %d",fLines.GetEntriesFast()));
764
765   if(fLines.GetEntriesFast()>1){
766     retcode=0;
767     //  find a first candidate for the primary vertex
768     fVert3D=AliVertexerTracks::TrackletVertexFinder(&fLines,0); 
769     // make a further selection on tracklets based on this first candidate
770     fVert3D.GetXYZ(peak);
771     AliDebug(1,Form("FIRST V candidate: %f ; %f ; %f",peak[0],peak[1],peak[2]));
772     Int_t *validate2 = new Int_t [fLines.GetEntriesFast()];
773     for(Int_t i=0; i<fLines.GetEntriesFast();i++) validate2[i]=1; 
774     for(Int_t i=0; i<fLines.GetEntriesFast();i++){
775       if(validate2[i]==0) continue; 
776       AliStrLine *l1 = (AliStrLine*)fLines.At(i);
777       if(l1->GetDistFromPoint(peak)> fDCAcut)fLines.RemoveAt(i);
778       if(optCuts==2){ // temporarily only for pileup
779         for(Int_t j=i+1; j<fLines.GetEntriesFast();j++){
780           AliStrLine *l2 = (AliStrLine*)fLines.At(j);
781           if(l1->GetDCA(l2)<0.00001){ 
782             Int_t delta1=(Int_t)l1->GetIdPoint(0)-(Int_t)l2->GetIdPoint(0);
783             Int_t delta2=(Int_t)l1->GetIdPoint(1)-(Int_t)l2->GetIdPoint(1);
784             Int_t deltamod1=(Int_t)l1->GetIdPoint(0)/kMaxCluPerMod
785               -(Int_t)l2->GetIdPoint(0)/kMaxCluPerMod;
786             Int_t deltamod2=(Int_t)l1->GetIdPoint(1)/kMaxCluPerMod
787               -(Int_t)l2->GetIdPoint(1)/kMaxCluPerMod;
788             // remove tracklets sharing a point
789             if( (delta1==0 && deltamod2==0)  || 
790                 (delta2==0 && deltamod1==0)  ) validate2[j]=0; 
791           }
792         }
793       }
794     }
795     for(Int_t i=0; i<fLines.GetEntriesFast();i++){
796       if(validate2[i]==0)  fLines.RemoveAt(i);
797     }
798     delete [] validate2;
799     fLines.Compress();
800     AliDebug(1,Form("Number of tracklets (after 3rd compression) %d",fLines.GetEntriesFast()));
801     if(fLines.GetEntriesFast()>1){// this new tracklet selection is used
802       fVert3D=AliVertexerTracks::TrackletVertexFinder(&fLines,0);
803     }
804   }
805   return retcode;  
806 }
807
808 //________________________________________________________
809 void AliITSVertexer3D::SetMeanPPtSelTracks(Double_t fieldTesla){
810   // Sets mean values of Pt based on the field
811   // for P (used in multiple scattering) the most probable value is used
812   if(TMath::Abs(fieldTesla-0.5)<0.01){
813     SetMeanPSelTracks(0.375);
814     SetMeanPtSelTracks(0.630);
815   }else if(TMath::Abs(fieldTesla-0.4)<0.01){
816     SetMeanPSelTracks(0.375);
817     SetMeanPtSelTracks(0.580);
818   }else if(TMath::Abs(fieldTesla-0.2)<0.01){
819     SetMeanPSelTracks(0.375);
820     SetMeanPtSelTracks(0.530);
821   }else if(fieldTesla<0.00001){
822     SetMeanPSelTracks(0.375);
823     SetMeanPtSelTracks(0.230);
824   }else{
825     SetMeanPSelTracks();
826     SetMeanPtSelTracks();
827   }
828 }
829
830 //________________________________________________________
831 void AliITSVertexer3D::FindPeaks(TH3F* histo, Double_t *peak, Int_t &nOfTracklets, Int_t &nOfTimes){
832   // Finds bin with max contents in 3D histo of tracket intersections
833   TAxis *xax = histo->GetXaxis();  
834   TAxis *yax = histo->GetYaxis();
835   TAxis *zax = histo->GetZaxis();
836   peak[0]=0.;
837   peak[1]=0.;
838   peak[2]=0.;
839   nOfTracklets = 0;
840   nOfTimes=0;
841   Int_t peakbin[3]={0,0,0};
842   Int_t peak2bin[3]={-1,-1,-1};
843   Int_t bc2=-1;
844   for(Int_t i=xax->GetFirst();i<=xax->GetLast();i++){
845     Double_t xval = xax->GetBinCenter(i);
846     for(Int_t j=yax->GetFirst();j<=yax->GetLast();j++){
847       Double_t yval = yax->GetBinCenter(j);
848       for(Int_t k=zax->GetFirst();k<=zax->GetLast();k++){
849         Double_t zval = zax->GetBinCenter(k);
850         Int_t bc =(Int_t)histo->GetBinContent(i,j,k);
851         if(bc==0) continue;
852         if(bc>nOfTracklets){
853           nOfTracklets=bc;
854           peak[2] = zval;
855           peak[1] = yval;
856           peak[0] = xval;
857           peakbin[2] = k;
858           peakbin[1] = j;
859           peakbin[0] = i;
860           peak2bin[2] = -1;
861           peak2bin[1] = -1;
862           peak2bin[0] = -1;
863           bc2=-1;
864           nOfTimes = 1;
865         }else if(bc==nOfTracklets){
866           if(TMath::Abs(i-peakbin[0])<=1 && TMath::Abs(j-peakbin[1])<=1 && TMath::Abs(k-peakbin[2])<=1){
867             peak2bin[2] = k;
868             peak2bin[1] = j;
869             peak2bin[0] = i;
870             bc2=bc;
871             nOfTimes = 1;
872           }else{
873             nOfTimes++;
874           }
875         }
876       }
877     }
878   }
879   if(peak2bin[0]>=-1 && bc2!=-1){ // two contiguous peak-cells with same contents
880     peak[0]=0.5*(xax->GetBinCenter(peakbin[0])+xax->GetBinCenter(peak2bin[0]));
881     peak[1]=0.5*(yax->GetBinCenter(peakbin[1])+yax->GetBinCenter(peak2bin[1]));
882     peak[2]=0.5*(zax->GetBinCenter(peakbin[2])+zax->GetBinCenter(peak2bin[2]));
883     nOfTracklets+=bc2;
884     nOfTimes=1;
885   }
886 }
887 //________________________________________________________
888 void AliITSVertexer3D::MarkUsedClusters(){
889   // Mark clusters of tracklets used in vertex claulation
890   for(Int_t i=0; i<fLines.GetEntriesFast();i++){
891     AliStrLine *lin = (AliStrLine*)fLines.At(i);
892     Int_t idClu1=lin->GetIdPoint(0);
893     Int_t idClu2=lin->GetIdPoint(1);
894     fUsedCluster.SetBitNumber(idClu1);
895     fUsedCluster.SetBitNumber(idClu2);
896   }
897 }
898 //________________________________________________________
899 Int_t AliITSVertexer3D::RemoveTracklets(){
900   // Remove trackelts close to first found vertex
901   Double_t vert[3]={fVert3D.GetXv(),fVert3D.GetYv(),fVert3D.GetZv()};
902   Int_t nRemoved=0;
903   for(Int_t i=0; i<fLines.GetEntriesFast();i++){
904     AliStrLine *lin = (AliStrLine*)fLines.At(i);
905     if(lin->GetDistFromPoint(vert)<fDCAforPileup){
906       Int_t idClu1=lin->GetIdPoint(0);
907       Int_t idClu2=lin->GetIdPoint(1);
908       fUsedCluster.SetBitNumber(idClu1);
909       fUsedCluster.SetBitNumber(idClu2);
910       fLines.RemoveAt(i);
911       ++nRemoved;
912     }
913   }
914   fLines.Compress();
915   return nRemoved;
916 }
917 //________________________________________________________
918 void AliITSVertexer3D::FindOther3DVertices(TTree *itsClusterTree){
919   // pileup identification based on 3D vertexing with not used clusters
920   MarkUsedClusters();
921   fLines.Clear("C");
922   Int_t nolines = FindTracklets(itsClusterTree,2);
923   if(nolines>=2){
924     Int_t nr=RemoveTracklets();
925     nolines-=nr;
926     if(nolines>=2){
927       Int_t rc=Prepare3DVertex(2);
928       if(rc==0){ 
929         fVert3D=AliVertexerTracks::TrackletVertexFinder(&fLines,0);
930         if(fVert3D.GetNContributors()>=fMinTrackletsForPilup){
931           fIsPileup=kTRUE;
932           fNoVertices=2;
933           fVertArray = new AliESDVertex[2];
934           fVertArray[0]=(*fCurrentVertex);
935           fVertArray[1]=fVert3D;
936           fZpuv=fVert3D.GetZv();
937           fNTrpuv=fVert3D.GetNContributors();
938         }
939       }
940     }
941   }
942 }
943 //______________________________________________________________________
944 void AliITSVertexer3D::PileupFromZ(){
945   // Calls the pileup algorithm of ALiITSVertexerZ
946   Int_t binmin, binmax;
947   Int_t nPeaks=AliITSVertexerZ::GetPeakRegion(fZHisto,binmin,binmax);   
948   if(nPeaks==2)AliWarning("2 peaks found");
949   Int_t firstPeakCont=0;
950   Double_t firstPeakPos=0.;
951   for(Int_t i=binmin-1;i<=binmax+1;i++){
952     firstPeakCont+=static_cast<Int_t>(fZHisto->GetBinContent(i));
953     firstPeakPos+=fZHisto->GetBinContent(i)*fZHisto->GetBinCenter(i);
954   }
955   if(firstPeakCont>0){ 
956     firstPeakPos/=firstPeakCont;
957     Int_t ncontr2=0;
958     if(firstPeakCont>fMinTrackletsForPilup){     
959       Float_t secPeakPos;
960       ncontr2=AliITSVertexerZ::FindSecondPeak(fZHisto,binmin,binmax,secPeakPos);
961       if(ncontr2>=fMinTrackletsForPilup){ 
962         fIsPileup=kTRUE;
963         fNoVertices=2;
964         AliESDVertex secondVert(secPeakPos,0.1,ncontr2);
965         fVertArray = new AliESDVertex[2];
966         fVertArray[0]=(*fCurrentVertex);
967         fVertArray[1]=secondVert;
968         fZpuv=secPeakPos;
969         fNTrpuv=ncontr2;
970       }
971     }
972   }
973 }
974 //________________________________________________________
975 void AliITSVertexer3D::PrintStatus() const {
976   // Print current status
977   printf("========= First step selections =====================\n");
978   printf("Cut on diamond (Z) %f\n",fZCutDiamond);
979   printf("Loose cut on Delta Phi %f\n",fCoarseDiffPhiCut);
980   printf("Loose cut on tracklet DCA to Z axis %f\n",fCoarseMaxRCut);
981   printf("Cut on DCA - tracklet to tracklet and to vertex %f\n",fDCAcut);
982   printf("========= Second step selections ====================\n");
983   printf("Cut on tracklet-to-first-vertex Z distance %f\n",fMaxZCut);
984   printf("Max Phi difference: %f\n",fDiffPhiMax);
985   printf("Cut on tracklet DCA to beam axis %f\n",fMaxRCut);
986   printf("Cut on tracklet DCA to beam axis (algo2) %f\n",fMaxRCut2);
987   printf("========= Pileup selections =========================\n");
988   printf("Pileup algo: %d\n",fPileupAlgo);
989   printf("Min DCA to 1st vertex for pileup (algo 0 and 1): %f\n",fDCAforPileup);
990   printf("Cut on distance between pair-vertices  (algo 2): %f\n",fCutOnPairs);
991   printf("Maximum number of clusters allowed on L1 or L2: %d\n",fMaxNumOfCl);
992   printf("=======================================================\n");
993 }
994