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