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