]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSVertexer3D.cxx
bug fix
[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 "AliESDVertex.h"
17 #include "AliLog.h"
18 #include "AliStrLine.h"
19 #include "AliTracker.h"
20 #include "AliITSDetTypeRec.h"
21 #include "AliITSRecPoint.h"
22 #include "AliITSgeomTGeo.h"
23 #include "AliVertexerTracks.h"
24 #include "AliITSVertexer3D.h"
25 #include "AliITSVertexerZ.h"
26 #include "AliITSSortTrkl.h"
27 /////////////////////////////////////////////////////////////////
28 // this class implements a method to determine
29 // the 3 coordinates of the primary vertex
30 // for p-p collisions 
31 // It can be used successfully with Pb-Pb collisions
32 ////////////////////////////////////////////////////////////////
33
34 ClassImp(AliITSVertexer3D)
35
36 /* $Id$ */
37
38 //______________________________________________________________________
39 AliITSVertexer3D::AliITSVertexer3D():AliITSVertexer(),
40 fLines("AliStrLine",1000),
41 fVert3D(),
42 fCoarseDiffPhiCut(0.),
43 fFineDiffPhiCut(0.),
44 fCutOnPairs(0.),
45 fCoarseMaxRCut(0.),
46 fMaxRCut(0.),
47 fZCutDiamond(0.),
48 fMaxZCut(0.),
49 fDCAcut(0.),
50 fDiffPhiMax(0.),
51 fMeanPSelTrk(0.),
52 fMeanPtSelTrk(0.),
53 fUsedCluster(kMaxCluPerMod*kNSPDMod),
54 fZHisto(0),
55 fDCAforPileup(0.),
56 fPileupAlgo(0)
57 {
58   // Default constructor
59   SetCoarseDiffPhiCut();
60   SetFineDiffPhiCut();
61   SetCutOnPairs();
62   SetCoarseMaxRCut();
63   SetMaxRCut();
64   SetZCutDiamond();
65   SetMaxZCut();
66   SetDCACut();
67   SetDiffPhiMax();
68   SetMeanPSelTracks();
69   SetMeanPtSelTracks();
70   SetMinDCAforPileup();
71   SetPileupAlgo();
72   Float_t binsize=0.02; // default 200 micron
73   Int_t nbins=static_cast<Int_t>(1+2*fZCutDiamond/binsize);
74   fZHisto=new TH1F("hz","",nbins,-fZCutDiamond,-fZCutDiamond+binsize*nbins);
75 }
76
77 //______________________________________________________________________
78 AliITSVertexer3D::~AliITSVertexer3D() {
79   // Destructor
80   fLines.Clear("C");
81   if(fZHisto) delete fZHisto;
82 }
83
84 //______________________________________________________________________
85 void AliITSVertexer3D::ResetVert3D(){
86   //
87   ResetVertex();
88   fVert3D.SetXv(0.);
89   fVert3D.SetYv(0.);
90   fVert3D.SetZv(0.);
91   fVert3D.SetDispersion(0.);
92   fVert3D.SetNContributors(0);
93   fUsedCluster.ResetAllBits(0);
94 }
95 //______________________________________________________________________
96 AliESDVertex* AliITSVertexer3D::FindVertexForCurrentEvent(TTree *itsClusterTree){
97   // Defines the AliESDVertex for the current event
98   ResetVert3D();
99   AliDebug(1,"FindVertexForCurrentEvent - 3D - PROCESSING NEXT EVENT");
100   fLines.Clear("C");
101   fCurrentVertex = NULL;
102
103   Int_t nolines = FindTracklets(itsClusterTree,0);
104   if(nolines>=2){
105     Int_t rc=Prepare3DVertex(0);
106     if(fPileupAlgo == 2 && rc == 0) FindVertex3DIterative();
107     else if(fPileupAlgo<2 && rc == 0) FindVertex3D(itsClusterTree);
108   }
109
110   if(!fCurrentVertex){
111     AliITSVertexerZ vertz(GetNominalPos()[0],GetNominalPos()[1]);
112     vertz.SetDetTypeRec(GetDetTypeRec());
113     AliDebug(1,"Call Vertexer Z\n");
114     vertz.SetLowLimit(-fZCutDiamond);
115     vertz.SetHighLimit(fZCutDiamond);
116     AliESDVertex* vtxz = vertz.FindVertexForCurrentEvent(itsClusterTree);
117     if(vtxz){
118       Double_t position[3]={GetNominalPos()[0],GetNominalPos()[1],vtxz->GetZv()};
119       Double_t covmatrix[6];
120       vtxz->GetCovMatrix(covmatrix);
121       Double_t chi2=99999.;
122       Int_t    nContr=vtxz->GetNContributors();
123       fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);    
124       fCurrentVertex->SetTitle("vertexer: Z");
125       fCurrentVertex->SetName("SPDVertexZ");
126       delete vtxz;
127     }
128
129   }
130   FindMultiplicity(itsClusterTree);
131   return fCurrentVertex;
132 }  
133
134 //______________________________________________________________________
135 void AliITSVertexer3D::FindVertex3D(TTree *itsClusterTree){
136   // 3D algorithm
137   fVert3D=AliVertexerTracks::TrackletVertexFinder(&fLines,0);
138   /*  uncomment to debug
139       printf("Vertex found in first iteration:\n");
140       fVert3D.Print();
141       printf("Start second iteration\n");
142       end of debug lines  */
143   if(fVert3D.GetNContributors()>0){
144     fLines.Clear("C");
145     Int_t nolines = FindTracklets(itsClusterTree,1);
146     if(nolines>=2){
147       Int_t rc=Prepare3DVertex(1);
148       if(rc==0) fVert3D=AliVertexerTracks::TrackletVertexFinder(&fLines,0);
149     }
150   }
151   /*  uncomment to debug 
152       printf("Vertex found in second iteration:\n");
153       fVert3D.Print();
154       end of debug lines  */ 
155   
156   Float_t vRadius=TMath::Sqrt(fVert3D.GetXv()*fVert3D.GetXv()+fVert3D.GetYv()*fVert3D.GetYv());
157   if(vRadius<GetPipeRadius() && fVert3D.GetNContributors()>0){
158     Double_t position[3]={fVert3D.GetXv(),fVert3D.GetYv(),fVert3D.GetZv()};
159     Double_t covmatrix[6];
160     fVert3D.GetCovMatrix(covmatrix);
161     Double_t chi2=99999.;
162     Int_t    nContr=fVert3D.GetNContributors();
163     fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);    
164     fCurrentVertex->SetTitle("vertexer: 3D");
165     fCurrentVertex->SetName("SPDVertex3D");
166     fCurrentVertex->SetDispersion(fVert3D.GetDispersion());
167     fNoVertices=1;
168     
169     switch(fPileupAlgo){
170     case 0: PileupFromZ(); break;
171     case 1: FindOther3DVertices(itsClusterTree); break;
172     default: AliError("Wrong pileup algorithm"); break;
173     }
174     if(fNoVertices==1){
175       fVertArray = new AliESDVertex[1];
176       fVertArray[0]=(*fCurrentVertex);    
177     }
178   }
179 }
180
181 //______________________________________________________________________
182 void AliITSVertexer3D::FindVertex3DIterative(){
183   // Defines the AliESDVertex for the current event
184   Int_t numsor=fLines.GetEntriesFast()*(fLines.GetEntriesFast()-1)/2;
185   //cout<<"AliITSVertexer3D::FindVertexForCurentEvent: Number of tracklets selected for vertexing "<<fLines.GetEntriesFast()<<"; Number of pairs: "<<numsor<<endl;
186   AliITSSortTrkl srt(fLines,numsor,fCutOnPairs,fCoarseMaxRCut); 
187   srt.FindClusters();   
188   AliInfo(Form("Number of vertices: %d",srt.GetNumberOfClusters()));
189       
190   fNoVertices = srt.GetNumberOfClusters();
191   //printf("fNoVertices = %d \n",fNoVertices);
192   if(fNoVertices>0){
193     fVertArray = new AliESDVertex[fNoVertices];
194     for(Int_t kk=0; kk<srt.GetNumberOfClusters(); kk++){
195       Int_t size = 0;
196       Int_t *labels = srt.GetTrackletsLab(kk,size);
197       /*
198         Int_t *pairs = srt.GetClusters(kk);
199         Int_t nopai = srt.GetSizeOfCluster(kk);
200         cout<<"***** Vertex number "<<kk<<".  Pairs: \n";
201         for(Int_t jj=0;jj<nopai;jj++){
202         cout<<pairs[jj]<<" - ";
203         if(jj>0 & jj%8==0)cout<<endl;
204         }
205         cout<<endl;
206         cout<<"***** Vertex number "<<kk<<".  Labels: \n";
207       */
208       AliStrLine **tclo = new AliStrLine* [size];
209       for(Int_t jj=0;jj<size;jj++){
210         //          cout<<labels[jj]<<" - ";
211         //          if(jj>0 & jj%8==0)cout<<endl;
212         tclo[jj] = dynamic_cast<AliStrLine*>(fLines[labels[jj]]);
213       }
214       //          cout<<endl;
215       delete []labels;
216       fVertArray[kk]=AliVertexerTracks::TrackletVertexFinder(tclo,size);
217       delete [] tclo;
218       //          fVertArray[kk].PrintStatus();
219       if(kk == 1){
220         // at least one second vertex is present
221         fIsPileup = kTRUE;
222         fNTrpuv = fVertArray[kk].GetNContributors();
223         fZpuv = fVertArray[kk].GetZv();
224       }
225     }
226     Float_t vRadius=TMath::Sqrt(fVertArray[0].GetXv()*fVertArray[0].GetXv()+fVertArray[0].GetYv()*fVertArray[0].GetYv());
227     if(vRadius<GetPipeRadius() && fVertArray[0].GetNContributors()>0){
228       Double_t position[3]={fVertArray[0].GetXv(),fVertArray[0].GetYv(),fVertArray[0].GetZv()};
229       Double_t covmatrix[6];
230       fVertArray[0].GetCovMatrix(covmatrix);
231       Double_t chi2=99999.;
232       Int_t    nContr=fVertArray[0].GetNContributors();
233       fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);    
234       fCurrentVertex->SetTitle("vertexer: 3D");
235       fCurrentVertex->SetName("SPDVertex3D");
236       fCurrentVertex->SetDispersion(fVertArray[0].GetDispersion());  
237     }
238   }
239
240 }  
241
242 //______________________________________________________________________
243 Bool_t AliITSVertexer3D::DistBetweenVertices(AliESDVertex &a, AliESDVertex &b, Double_t test, Double_t &dist){
244   // method to compare the distance between vertices a and b with "test"
245   //it returns kTRUE is the distance is less or equal to test
246   dist = (a.GetX()-b.GetX()) * (a.GetX()-b.GetX());
247   dist +=  (a.GetY()-b.GetY()) * (a.GetY()-b.GetY());
248   dist +=  (a.GetZ()-b.GetZ()) * (a.GetZ()-b.GetZ());
249   dist = TMath::Sqrt(dist);
250   if(dist <= test)return kTRUE;
251   return kFALSE;
252 }
253
254
255 //______________________________________________________________________
256 Int_t AliITSVertexer3D::FindTracklets(TTree *itsClusterTree, Int_t optCuts){
257   // All the possible combinations between recpoints on layer 1and 2 are
258   // considered. Straight lines (=tracklets)are formed. 
259   // The tracklets are processed in Prepare3DVertex
260
261   if(!GetDetTypeRec())AliFatal("DetTypeRec pointer has not been set");
262
263   TTree *tR = itsClusterTree;
264   fDetTypeRec->ResetRecPoints();
265   fDetTypeRec->SetTreeAddressR(tR);
266   TClonesArray *itsRec  = 0;
267   if(optCuts==0) fZHisto->Reset();
268  // gc1 are local and global coordinates for layer 1
269   Float_t gc1[3]={0.,0.,0.};
270   // gc2 are local and global coordinates for layer 2
271   Float_t gc2[3]={0.,0.,0.};
272
273   itsRec = fDetTypeRec->RecPoints();
274   TBranch *branch = NULL;
275   branch = tR->GetBranch("ITSRecPoints");
276   if(!branch){
277     AliError("Null pointer for RecPoints branch");
278     return -1;
279   }
280
281   // Set values for cuts
282   Float_t xbeam=GetNominalPos()[0]; 
283   Float_t ybeam=GetNominalPos()[1];
284   Float_t zvert=0.;
285   Float_t deltaPhi=fCoarseDiffPhiCut;
286   Float_t deltaR=fCoarseMaxRCut;
287   Float_t dZmax=fZCutDiamond;
288   if(fPileupAlgo == 2){
289     deltaPhi=fFineDiffPhiCut;
290     deltaR=fMaxRCut;
291     if(optCuts != 0)AliWarning(Form("fPileupAlgo=2 AND optCuts=%d has been selected. It should be 0",optCuts));
292   } else if(optCuts==1){
293     xbeam=fVert3D.GetXv();
294     ybeam=fVert3D.GetYv();
295     zvert=fVert3D.GetZv();
296     deltaPhi = fDiffPhiMax; 
297     deltaR=fMaxRCut;
298     dZmax=fMaxZCut;
299   } else if(optCuts==2){
300     xbeam=fVert3D.GetXv();
301     ybeam=fVert3D.GetYv();
302     deltaPhi = fDiffPhiMax; 
303     deltaR=fMaxRCut;
304   }
305
306   Int_t nrpL1 = 0;    // number of rec points on layer 1
307   Int_t nrpL2 = 0;    // number of rec points on layer 2
308
309   // By default irstL1=0 and lastL1=79
310   Int_t firstL1 = AliITSgeomTGeo::GetModuleIndex(1,1,1);
311   Int_t lastL1 = AliITSgeomTGeo::GetModuleIndex(2,1,1)-1;
312   for(Int_t module= firstL1; module<=lastL1;module++){  // count number of recopints on layer 1
313     branch->GetEvent(module);
314     nrpL1+= itsRec->GetEntries();
315     fDetTypeRec->ResetRecPoints();
316   }
317   //By default firstL2=80 and lastL2=239
318   Int_t firstL2 = AliITSgeomTGeo::GetModuleIndex(2,1,1);
319   Int_t lastL2 = AliITSgeomTGeo::GetModuleIndex(3,1,1)-1;
320   for(Int_t module= firstL2; module<=lastL2;module++){  // count number of recopints on layer 2
321     branch->GetEvent(module);
322     nrpL2+= itsRec->GetEntries();
323     fDetTypeRec->ResetRecPoints();
324   }
325   if(nrpL1 == 0 || nrpL2 == 0){
326     return -1;
327   }
328   AliDebug(1,Form("RecPoints on Layer 1,2 = %d, %d\n",nrpL1,nrpL2));
329
330   Double_t a[3]={xbeam,ybeam,0.}; 
331   Double_t b[3]={xbeam,ybeam,10.};
332   AliStrLine zeta(a,b,kTRUE);
333   Float_t bField=AliTracker::GetBz()/10.; //T
334   SetMeanPPtSelTracks(bField);
335
336   Int_t nolines = 0;
337   // Loop on modules of layer 1
338   for(Int_t modul1= firstL1; modul1<=lastL1;modul1++){   // Loop on modules of layer 1
339     if(!fUseModule[modul1]) continue;
340     UShort_t ladder=int(modul1/4)+1; // ladders are numbered starting from 1
341     branch->GetEvent(modul1);
342     Int_t nrecp1 = itsRec->GetEntries();
343     static TClonesArray prpl1("AliITSRecPoint",nrecp1);
344     prpl1.SetOwner();
345     for(Int_t j=0;j<nrecp1;j++){
346       AliITSRecPoint *recp = (AliITSRecPoint*)itsRec->At(j);
347       new(prpl1[j])AliITSRecPoint(*recp);
348     }
349     fDetTypeRec->ResetRecPoints();
350     for(Int_t j=0;j<nrecp1;j++){
351       if(j>kMaxCluPerMod) continue;
352       UShort_t idClu1=modul1*kMaxCluPerMod+j;
353       if(fUsedCluster.TestBitNumber(idClu1)) continue;
354       AliITSRecPoint *recp1 = (AliITSRecPoint*)prpl1.At(j);
355       recp1->GetGlobalXYZ(gc1);
356       Double_t phi1 = TMath::ATan2(gc1[1]-ybeam,gc1[0]-xbeam);
357       if(phi1<0)phi1=2*TMath::Pi()+phi1;
358       for(Int_t ladl2=0 ; ladl2<fLadOnLay2*2+1;ladl2++){
359         for(Int_t k=0;k<4;k++){
360           Int_t ladmod=fLadders[ladder-1]+ladl2;
361           if(ladmod>AliITSgeomTGeo::GetNLadders(2)) ladmod=ladmod-AliITSgeomTGeo::GetNLadders(2);
362           Int_t modul2=AliITSgeomTGeo::GetModuleIndex(2,ladmod,k+1);
363           if(!fUseModule[modul2]) continue;
364           branch->GetEvent(modul2);
365           Int_t nrecp2 = itsRec->GetEntries();
366           for(Int_t j2=0;j2<nrecp2;j2++){
367             if(j2>kMaxCluPerMod) continue;
368             UShort_t idClu2=modul2*kMaxCluPerMod+j2;
369             if(fUsedCluster.TestBitNumber(idClu2)) continue;
370             AliITSRecPoint *recp2 = (AliITSRecPoint*)itsRec->At(j2);
371             recp2->GetGlobalXYZ(gc2);
372             Double_t phi2 = TMath::ATan2(gc2[1]-ybeam,gc2[0]-xbeam);
373             if(phi2<0)phi2=2*TMath::Pi()+phi2;
374             Double_t diff = TMath::Abs(phi2-phi1); 
375             if(diff>TMath::Pi())diff=2.*TMath::Pi()-diff; 
376             if(optCuts==0 && diff<fDiffPhiMax){
377               Float_t r1=TMath::Sqrt(gc1[0]*gc1[0]+gc1[1]*gc1[1]);
378               Float_t zc1=gc1[2];
379               Float_t r2=TMath::Sqrt(gc2[0]*gc2[0]+gc2[1]*gc2[1]);
380               Float_t zc2=gc2[2];
381               Float_t zr0=(r2*zc1-r1*zc2)/(r2-r1); //Z @ null radius
382               fZHisto->Fill(zr0);
383             }
384             if(diff>deltaPhi)continue;
385             AliStrLine line(gc1,gc2,kTRUE);
386             Double_t cp[3];
387             Int_t retcode = line.Cross(&zeta,cp);
388             if(retcode<0)continue;
389             Double_t dca = line.GetDCA(&zeta);
390             if(dca<0.) continue;
391             if(dca>deltaR)continue;
392             Double_t deltaZ=cp[2]-zvert;
393             if(TMath::Abs(deltaZ)>dZmax)continue;
394
395             if(nolines == 0){
396               if(fLines.GetEntriesFast()>0)fLines.Clear("C");
397             }
398             Float_t cov[6];
399             recp2->GetGlobalCov(cov);
400
401
402             Float_t rad1=TMath::Sqrt(gc1[0]*gc1[0]+gc1[1]*gc1[1]);
403             Float_t rad2=TMath::Sqrt(gc2[0]*gc2[0]+gc2[1]*gc2[1]);
404             Float_t factor=(rad1+rad2)/(rad2-rad1); //factor to account for error on tracklet direction 
405
406             Float_t curvErr=0;
407             if(bField>0.00001){
408               Float_t curvRadius=fMeanPtSelTrk/(0.3*bField)*100; //cm 
409               Float_t dRad=TMath::Sqrt(TMath::Power((gc1[0]-gc2[0]),2)+TMath::Power((gc1[1]-gc2[1]),2));
410               Float_t aux=dRad/2.+rad1;
411               curvErr=TMath::Sqrt(curvRadius*curvRadius-dRad*dRad/4.)-TMath::Sqrt(curvRadius*curvRadius-aux*aux); //cm
412             }
413             Float_t sigmasq[3];
414             sigmasq[0]=(cov[0]+curvErr*curvErr/2.)*factor*factor;
415             sigmasq[1]=(cov[3]+curvErr*curvErr/2.)*factor*factor;
416             sigmasq[2]=cov[5]*factor*factor;
417
418             // Multiple scattering
419             Float_t beta=1.;
420             Float_t beta2=beta*beta;
421             Float_t p2=fMeanPSelTrk*fMeanPSelTrk;
422             Float_t rBP=GetPipeRadius();
423             Float_t dBP=0.08/35.3; // 800 um of Be
424             Float_t dL1=0.01; //approx. 1% of radiation length  
425             Float_t theta2BP=14.1*14.1/(beta2*p2*1e6)*TMath::Abs(dBP);
426             Float_t theta2L1=14.1*14.1/(beta2*p2*1e6)*TMath::Abs(dL1);
427             Float_t thetaBP=TMath::Sqrt(theta2BP);
428             Float_t thetaL1=TMath::Sqrt(theta2L1);
429             for(Int_t ico=0; ico<3;ico++){    
430 //            printf("Error on coord. %d due to cov matrix+curvErr=%f\n",ico,sigmasq[ico]);
431 //            //              sigmasq[ico]+=rad1*rad1*geomfac[ico]*theta2L1/2; // multiple scattering in layer 1
432 //            //  sigmasq[ico]+=rBP*rBP*geomfac[ico]*theta2BP/2; // multiple scattering in beam pipe
433               sigmasq[ico]+=TMath::Power(rad1*TMath::Tan(thetaL1),2)/3.;
434               sigmasq[ico]+=TMath::Power(rBP*TMath::Tan(thetaBP),2)/3.;
435
436 //            printf("Multipl. scatt. contr %d = %f (LAY1), %f (BP)\n",ico,rad1*rad1*geomfac[ico]*theta2L1/2,rBP*rBP*geomfac[ico]*theta2BP/2);
437 //            printf("Total error on coord %d = %f\n",ico,sigmasq[ico]);
438             }
439             Float_t wmat[9]={1.,0.,0.,0.,1.,0.,0.,0.,1.};
440             if(sigmasq[0]!=0.) wmat[0]=1./sigmasq[0];
441             if(sigmasq[1]!=0.) wmat[4]=1./sigmasq[1];
442             if(sigmasq[2]!=0.) wmat[8]=1./sigmasq[2];
443             new(fLines[nolines++])AliStrLine(gc1,sigmasq,wmat,gc2,kTRUE,idClu1,idClu2);
444
445           }
446           fDetTypeRec->ResetRecPoints();
447         }
448       }
449     }
450     prpl1.Clear();
451   }
452   if(nolines == 0)return -2;
453   return nolines;
454 }
455
456 //______________________________________________________________________
457 Int_t  AliITSVertexer3D::Prepare3DVertex(Int_t optCuts){
458   // Finds the 3D vertex information using tracklets
459   Int_t retcode = -1;
460
461   Float_t xbeam=GetNominalPos()[0];
462   Float_t ybeam=GetNominalPos()[1];
463   Float_t zvert=0.;
464   Float_t deltaR=fCoarseMaxRCut;
465   if(fPileupAlgo == 2) {
466     deltaR=fMaxRCut;
467     if(optCuts!=0)AliWarning(Form("fPileupAlgo=2 AND optCuts=%d. It should be 0",optCuts));
468   }
469   Float_t dZmax=fZCutDiamond;
470   if(optCuts==1){
471     xbeam=fVert3D.GetXv();
472     ybeam=fVert3D.GetYv();
473     zvert=fVert3D.GetZv();
474     deltaR=fMaxRCut;
475     dZmax=fMaxZCut;
476   }else if(optCuts==2){
477     xbeam=fVert3D.GetXv();
478     ybeam=fVert3D.GetYv();
479     deltaR=fMaxRCut;
480   }
481
482   Int_t nbr=50;
483   Float_t rl=-fCoarseMaxRCut;
484   Float_t rh=fCoarseMaxRCut;
485   Int_t nbz=100;
486   Float_t zl=-fZCutDiamond;
487   Float_t zh=fZCutDiamond;
488   Float_t binsizer=(rh-rl)/nbr;
489   Float_t binsizez=(zh-zl)/nbz;
490   Int_t nbrcs=25;
491   Int_t nbzcs=50;
492   TH3F *h3d = NULL;
493   TH3F *h3dcs = NULL;
494   if(fPileupAlgo !=2){
495     h3d = new TH3F("h3d","xyz distribution",nbr,rl,rh,nbr,rl,rh,nbz,zl,zh);
496     h3dcs = new TH3F("h3dcs","xyz distribution",nbrcs,rl,rh,nbrcs,rl,rh,nbzcs,zl,zh);
497   }
498   // cleanup of the TCLonesArray of tracklets (i.e. fakes are removed)
499   Int_t *validate = new Int_t [fLines.GetEntriesFast()];
500   for(Int_t i=0; i<fLines.GetEntriesFast();i++)validate[i]=0;
501   for(Int_t i=0; i<fLines.GetEntriesFast()-1;i++){
502     AliStrLine *l1 = (AliStrLine*)fLines.At(i);
503     for(Int_t j=i+1;j<fLines.GetEntriesFast();j++){
504       AliStrLine *l2 = (AliStrLine*)fLines.At(j);
505       Float_t dca=l1->GetDCA(l2);
506       if(dca > fDCAcut || dca<0.00001) continue;
507       Double_t point[3];
508       Int_t retc = l1->Cross(l2,point);
509       if(retc<0)continue;
510       Double_t deltaZ=point[2]-zvert;
511      if(TMath::Abs(deltaZ)>dZmax)continue;
512       Double_t rad=TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
513       if(rad>fCoarseMaxRCut)continue;
514       Double_t deltaX=point[0]-xbeam;
515       Double_t deltaY=point[1]-ybeam;
516       Double_t raddist=TMath::Sqrt(deltaX*deltaX+deltaY*deltaY);
517       if(raddist>deltaR)continue;
518       validate[i]=1;
519       validate[j]=1;
520       if(fPileupAlgo != 2){
521         h3d->Fill(point[0],point[1],point[2]);
522         h3dcs->Fill(point[0],point[1],point[2]);
523       }
524     }
525   }
526
527
528
529   Int_t numbtracklets=0;
530   for(Int_t i=0; i<fLines.GetEntriesFast();i++)if(validate[i]>=1)numbtracklets++;
531   if(numbtracklets<2){
532     delete [] validate; 
533     if(fPileupAlgo != 2){
534       delete h3d; 
535       delete h3dcs; 
536     }
537     return retcode; 
538   }
539
540   for(Int_t i=0; i<fLines.GetEntriesFast();i++){
541     if(validate[i]<1)fLines.RemoveAt(i);
542   }
543   fLines.Compress();
544   AliDebug(1,Form("Number of tracklets (after compress)%d ",fLines.GetEntriesFast()));
545   delete [] validate;
546
547   // Exit here if Pileup Algorithm 2 has been chosen
548   if(fPileupAlgo == 2)return 0;
549   //         
550
551
552   //        Find peaks in histos
553
554   Double_t peak[3]={0.,0.,0.};
555   Int_t ntrkl,ntimes;
556   FindPeaks(h3d,peak,ntrkl,ntimes);  
557   delete h3d;
558
559   if(optCuts==0 && ntrkl<=2){
560     ntrkl=0;
561     ntimes=0;
562     FindPeaks(h3dcs,peak,ntrkl,ntimes);  
563     binsizer=(rh-rl)/nbrcs;
564     binsizez=(zh-zl)/nbzcs;
565     if(ntrkl==1 || ntimes>1){delete h3dcs; return retcode;}
566   }
567   delete h3dcs;
568
569   //         Second selection loop
570
571   Float_t bs=(binsizer+binsizez)/2.;
572   for(Int_t i=0; i<fLines.GetEntriesFast();i++){
573     AliStrLine *l1 = (AliStrLine*)fLines.At(i);
574     if(l1->GetDistFromPoint(peak)>2.5*bs)fLines.RemoveAt(i);
575   }
576   fLines.Compress();
577   AliDebug(1,Form("Number of tracklets (after 2nd compression) %d",fLines.GetEntriesFast()));
578
579   if(fLines.GetEntriesFast()>1){
580     //  find a first candidate for the primary vertex
581     fVert3D=AliVertexerTracks::TrackletVertexFinder(&fLines,0); 
582     // make a further selection on tracklets based on this first candidate
583     fVert3D.GetXYZ(peak);
584     AliDebug(1,Form("FIRST V candidate: %f ; %f ; %f",peak[0],peak[1],peak[2]));
585     for(Int_t i=0; i<fLines.GetEntriesFast();i++){
586       AliStrLine *l1 = (AliStrLine*)fLines.At(i);
587       if(l1->GetDistFromPoint(peak)> fDCAcut)fLines.RemoveAt(i);
588     }
589     fLines.Compress();
590     AliDebug(1,Form("Number of tracklets (after 3rd compression) %d",fLines.GetEntriesFast()));
591     if(fLines.GetEntriesFast()>1) retcode=0; // this new tracklet selection is used
592     else retcode =1; // the previous tracklet selection will be used
593   }
594   else {
595     retcode = 0;
596   }
597   return retcode;  
598 }
599
600 //________________________________________________________
601 void AliITSVertexer3D::SetMeanPPtSelTracks(Float_t fieldTesla){
602   // Sets mean values of P and Pt based on the field
603   if(TMath::Abs(fieldTesla-0.5)<0.01){
604     SetMeanPSelTracks(0.885);
605     SetMeanPtSelTracks(0.630);
606   }else if(TMath::Abs(fieldTesla-0.4)<0.01){
607     SetMeanPSelTracks(0.805);
608     SetMeanPtSelTracks(0.580);
609   }else if(TMath::Abs(fieldTesla-0.2)<0.01){
610     SetMeanPSelTracks(0.740);
611     SetMeanPtSelTracks(0.530);
612   }else if(fieldTesla<0.00001){
613     SetMeanPSelTracks(0.730);
614     SetMeanPtSelTracks(0.510);
615   }else{
616     SetMeanPSelTracks();
617     SetMeanPtSelTracks();
618   }
619 }
620
621 //________________________________________________________
622 void AliITSVertexer3D::FindPeaks(TH3F* histo, Double_t *peak, Int_t &nOfTracklets, Int_t &nOfTimes){
623   // Finds bin with max contents in 3D histo of tracket intersections
624   TAxis *xax = histo->GetXaxis();  
625   TAxis *yax = histo->GetYaxis();
626   TAxis *zax = histo->GetZaxis();
627   peak[0]=0.;
628   peak[1]=0.;
629   peak[2]=0.;
630   nOfTracklets = 0;
631   nOfTimes=0;
632   for(Int_t i=xax->GetFirst();i<=xax->GetLast();i++){
633     Float_t xval = xax->GetBinCenter(i);
634     for(Int_t j=yax->GetFirst();j<=yax->GetLast();j++){
635       Float_t yval = yax->GetBinCenter(j);
636       for(Int_t k=zax->GetFirst();k<=zax->GetLast();k++){
637         Float_t zval = zax->GetBinCenter(k);
638         Int_t bc =(Int_t)histo->GetBinContent(i,j,k);
639         if(bc>nOfTracklets){
640           nOfTracklets = bc;
641           peak[2] = zval;
642           peak[1] = yval;
643           peak[0] = xval;
644           nOfTimes = 1;
645         }else if(bc==nOfTracklets){
646           nOfTimes++;
647         }
648       }
649     }
650   }
651 }
652
653 //________________________________________________________
654 void AliITSVertexer3D::MarkUsedClusters(){
655   // Mark clusters of tracklets used in vertex claulation
656   for(Int_t i=0; i<fLines.GetEntriesFast();i++){
657     AliStrLine *lin = (AliStrLine*)fLines.At(i);
658     Int_t idClu1=lin->GetIdPoint(0);
659     Int_t idClu2=lin->GetIdPoint(1);
660     fUsedCluster.SetBitNumber(idClu1);
661     fUsedCluster.SetBitNumber(idClu2);
662   }
663 }
664 //________________________________________________________
665 Int_t AliITSVertexer3D::RemoveTracklets(){
666   // Remove trackelts close to first found vertex
667   Double_t vert[3]={fVert3D.GetXv(),fVert3D.GetYv(),fVert3D.GetZv()};
668   Int_t nRemoved=0;
669   for(Int_t i=0; i<fLines.GetEntriesFast();i++){
670     AliStrLine *lin = (AliStrLine*)fLines.At(i);
671     if(lin->GetDistFromPoint(vert)<fDCAforPileup){
672       Int_t idClu1=lin->GetIdPoint(0);
673       Int_t idClu2=lin->GetIdPoint(1);
674       fUsedCluster.SetBitNumber(idClu1);
675       fUsedCluster.SetBitNumber(idClu2);
676       fLines.RemoveAt(i);
677       ++nRemoved;
678     }
679   }
680   fLines.Compress();
681   return nRemoved;
682 }
683 //________________________________________________________
684 void AliITSVertexer3D::FindOther3DVertices(TTree *itsClusterTree){
685   // pileup identification based on 3D vertexing with not used clusters
686   MarkUsedClusters();
687   fLines.Clear("C");
688   Int_t nolines = FindTracklets(itsClusterTree,2);
689   if(nolines>=2){
690     Int_t nr=RemoveTracklets();
691     nolines-=nr;
692     if(nolines>=2){
693       Int_t rc=Prepare3DVertex(2);
694       if(rc==0){ 
695         fVert3D=AliVertexerTracks::TrackletVertexFinder(&fLines,0);
696         if(fVert3D.GetNContributors()>=fMinTrackletsForPilup){
697           fIsPileup=kTRUE;
698           fNoVertices=2;
699           fVertArray = new AliESDVertex[2];
700           fVertArray[0]=(*fCurrentVertex);
701           fVertArray[1]=fVert3D;
702           fZpuv=fVert3D.GetZv();
703           fNTrpuv=fVert3D.GetNContributors();
704         }
705       }
706     }
707   }
708 }
709 //______________________________________________________________________
710 void AliITSVertexer3D::PileupFromZ(){
711   // Calls the pileup algorithm of ALiITSVertexerZ
712   Int_t binmin, binmax;
713   Int_t nPeaks=AliITSVertexerZ::GetPeakRegion(fZHisto,binmin,binmax);   
714   if(nPeaks==2)AliWarning("2 peaks found");
715   Int_t firstPeakCont=0;
716   Float_t firstPeakPos=0.;
717   for(Int_t i=binmin-1;i<=binmax+1;i++){
718     firstPeakCont+=static_cast<Int_t>(fZHisto->GetBinContent(i));
719     firstPeakPos+=fZHisto->GetBinContent(i)*fZHisto->GetBinCenter(i);
720   }
721   if(firstPeakCont>0){ 
722     firstPeakPos/=firstPeakCont;
723     Int_t ncontr2=0;
724     if(firstPeakCont>fMinTrackletsForPilup){     
725       Float_t secPeakPos;
726       ncontr2=AliITSVertexerZ::FindSecondPeak(fZHisto,binmin,binmax,secPeakPos);
727       if(ncontr2>=fMinTrackletsForPilup){ 
728         fIsPileup=kTRUE;
729         fNoVertices=2;
730         AliESDVertex secondVert(secPeakPos,0.1,ncontr2);
731         fVertArray = new AliESDVertex[2];
732         fVertArray[0]=(*fCurrentVertex);
733         fVertArray[1]=secondVert;
734         fZpuv=secPeakPos;
735         fNTrpuv=ncontr2;
736       }
737     }
738   }
739 }
740 //________________________________________________________
741 void AliITSVertexer3D::PrintStatus() const {
742   // Print current status
743   printf("=======================================================\n");
744   printf("Loose cut on Delta Phi %f\n",fCoarseDiffPhiCut);
745   printf("Cut on tracklet DCA to Z axis %f\n",fCoarseMaxRCut);
746   printf("Cut on tracklet DCA to beam axis %f\n",fMaxRCut);
747   printf("Cut on diamond (Z) %f\n",fZCutDiamond);
748   printf("Cut on DCA - tracklet to tracklet and to vertex %f\n",fDCAcut);
749   printf("Max Phi difference: %f\n",fDiffPhiMax);
750   printf("Pileup algo: %d\n",fPileupAlgo);
751   printf("Min DCA to 1st vetrtex for pileup: %f\n",fDCAforPileup);
752   printf("=======================================================\n");
753 }