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