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