1 /**************************************************************************
2 * Copyright(c) 2006-2008, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 #include "AliRunLoader.h"
17 #include "AliESDVertex.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
34 ////////////////////////////////////////////////////////////////
36 const Int_t AliITSVertexer3D::fgkMaxNumOfClDefault = 1000;
38 ClassImp(AliITSVertexer3D)
42 //______________________________________________________________________
43 AliITSVertexer3D::AliITSVertexer3D():
45 fLines("AliStrLine",1000),
47 fCoarseDiffPhiCut(0.),
59 fUsedCluster(kMaxCluPerMod*kNSPDMod),
62 fDiffPhiforPileup(0.),
66 fMaxNumOfCl(fgkMaxNumOfClDefault),
70 // Default constructor
71 SetCoarseDiffPhiCut();
84 SetDeltaPhiforPileup();
88 Double_t binsize=0.02; // default 200 micron
89 Int_t nbins=static_cast<Int_t>(1+2*fZCutDiamond/binsize);
90 fZHisto=new TH1F("hz","",nbins,-fZCutDiamond,-fZCutDiamond+binsize*nbins);
91 fGenerForDownScale=new TRandom3(987654321);
94 //______________________________________________________________________
95 AliITSVertexer3D::~AliITSVertexer3D() {
98 if(fZHisto) delete fZHisto;
99 if(fGenerForDownScale) delete fGenerForDownScale;
102 //______________________________________________________________________
103 void AliITSVertexer3D::ResetVert3D(){
109 fVert3D.SetDispersion(0.);
110 fVert3D.SetNContributors(0);
111 fUsedCluster.ResetAllBits(0);
113 //______________________________________________________________________
114 AliESDVertex* AliITSVertexer3D::FindVertexForCurrentEvent(TTree *itsClusterTree){
115 // Defines the AliESDVertex for the current event
117 AliDebug(1,"FindVertexForCurrentEvent - 3D - PROCESSING NEXT EVENT");
119 fCurrentVertex = NULL;
121 Int_t nolines = FindTracklets(itsClusterTree,0);
124 rc=Prepare3DVertex(0);
125 if(fVert3D.GetNContributors()>0){
127 nolines = FindTracklets(itsClusterTree,1);
129 rc=Prepare3DVertex(1);
130 if(fPileupAlgo == 2 && rc == 0) FindVertex3DIterative();
131 else if(fPileupAlgo!=2 && rc == 0) FindVertex3D(itsClusterTree);
132 if(rc!=0) fVert3D.SetNContributors(0); // exclude this vertex
138 AliITSVertexerZ vertz(GetNominalPos()[0],GetNominalPos()[1]);
139 vertz.SetDetTypeRec(GetDetTypeRec());
140 AliDebug(1,"Call Vertexer Z\n");
141 vertz.SetLowLimit(-fZCutDiamond);
142 vertz.SetHighLimit(fZCutDiamond);
143 AliESDVertex* vtxz = vertz.FindVertexForCurrentEvent(itsClusterTree);
145 Double_t position[3]={GetNominalPos()[0],GetNominalPos()[1],vtxz->GetZv()};
146 Double_t covmatrix[6];
147 vtxz->GetCovMatrix(covmatrix);
148 Double_t chi2=99999.;
149 Int_t nContr=vtxz->GetNContributors();
150 fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);
151 fCurrentVertex->SetDispersion(vtxz->GetDispersion());
152 fCurrentVertex->SetTitle("vertexer: Z");
153 fCurrentVertex->SetName("SPDVertexZ");
158 if(fComputeMultiplicity) FindMultiplicity(itsClusterTree);
159 return fCurrentVertex;
162 //______________________________________________________________________
163 void AliITSVertexer3D::FindVertex3D(TTree *itsClusterTree){
165 Double_t vRadius=TMath::Sqrt(fVert3D.GetXv()*fVert3D.GetXv()+fVert3D.GetYv()*fVert3D.GetYv());
166 if(vRadius<GetPipeRadius() && fVert3D.GetNContributors()>0){
167 Double_t position[3]={fVert3D.GetXv(),fVert3D.GetYv(),fVert3D.GetZv()};
168 Double_t covmatrix[6];
169 fVert3D.GetCovMatrix(covmatrix);
170 Double_t chi2=99999.;
171 Int_t nContr=fVert3D.GetNContributors();
172 fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);
173 fCurrentVertex->SetTitle("vertexer: 3D");
174 fCurrentVertex->SetName("SPDVertex3D");
175 fCurrentVertex->SetDispersion(fVert3D.GetDispersion());
179 case 0: PileupFromZ(); break;
180 case 1: FindOther3DVertices(itsClusterTree); break;
181 case 3: break; // no pileup algo
182 default: AliError("Wrong pileup algorithm"); break;
185 fVertArray = new AliESDVertex[1];
186 fVertArray[0]=(*fCurrentVertex);
191 //______________________________________________________________________
192 void AliITSVertexer3D::FindVertex3DIterative(){
195 Int_t nLines=fLines.GetEntriesFast();
196 Int_t maxPoints=nLines*(nLines-1)/2;
197 Double_t* xP=new Double_t[maxPoints];
198 Double_t* yP=new Double_t[maxPoints];
199 Double_t* zP=new Double_t[maxPoints];
200 Int_t* index1=new Int_t[maxPoints];
201 Int_t* index2=new Int_t[maxPoints];
202 Double_t xbeam=fVert3D.GetXv();
203 Double_t ybeam=fVert3D.GetYv();
206 for(Int_t ilin1=0; ilin1<nLines; ilin1++){
207 AliStrLine *l1 = (AliStrLine*)fLines.At(ilin1);
208 for(Int_t ilin2=ilin1+1; ilin2<nLines; ilin2++){
209 AliStrLine *l2 = (AliStrLine*)fLines.At(ilin2);
210 Double_t dca=l1->GetDCA(l2);
211 if(dca > fDCAcut || dca<0.00001) continue;
213 Int_t retc = l1->Cross(l2,point);
215 Double_t rad=TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
216 if(rad>fCoarseMaxRCut)continue;
217 Double_t distFromBeam=TMath::Sqrt((point[0]-xbeam)*(point[0]-xbeam)+(point[1]-ybeam)*(point[1]-ybeam));
218 if(distFromBeam>fMaxRCut2) continue;
222 index1[iPoint]=ilin1;
223 index2[iPoint]=ilin2;
227 Int_t npoints=iPoint++;
229 Short_t* mask=new Short_t[npoints];
230 for(Int_t ip=0;ip<npoints;ip++) mask[ip]=-1;
232 for(Int_t ip1=0;ip1<npoints;ip1++){
233 if(mask[ip1]==-1) mask[ip1]=index++;
234 for(Int_t ip2=ip1+1; ip2<npoints; ip2++){
235 if(mask[ip2]==mask[ip1] && mask[ip2]!=-1) continue;
236 Double_t dist2=(xP[ip1]-xP[ip2])*(xP[ip1]-xP[ip2]);
237 dist2+=(yP[ip1]-yP[ip2])*(yP[ip1]-yP[ip2]);
238 dist2+=(zP[ip1]-zP[ip2])*(zP[ip1]-zP[ip2]);
239 if(dist2<fCutOnPairs*fCutOnPairs){
240 if(mask[ip2]==-1) mask[ip2]=mask[ip1];
242 for(Int_t ip=0; ip<npoints;ip++){
243 if(mask[ip]==mask[ip2]) mask[ip]=mask[ip1];
251 // Count multiplicity of trackelts in clusters
252 UInt_t* isIndUsed=new UInt_t[index+1];
253 for(Int_t ind=0;ind<index+1;ind++) isIndUsed[ind]=0;
254 for(Int_t ip=0; ip<npoints;ip++){
259 // Count clusters/vertices and sort according to multiplicity
261 Int_t* sortedIndex=new Int_t[index+1];
262 for(Int_t ind1=0;ind1<index+1;ind1++){
263 if(isIndUsed[ind1]<=1) isIndUsed[ind1]=0;
266 if(ind1>0) cap=isIndUsed[sortedIndex[ind1-1]];
268 Int_t biggerindex=-1;
269 for(Int_t ind2=0;ind2<index+1;ind2++){
271 for(Int_t ind3=0; ind3<ind1; ind3++)
272 if(ind2==sortedIndex[ind3]) use=kFALSE;
273 if(use && isIndUsed[ind2]>bigger && isIndUsed[ind2]<=cap){
274 bigger=isIndUsed[ind2];
278 sortedIndex[ind1]=biggerindex;
280 AliDebug(3,Form("Number of clusters before merging = %d\n",nClusters));
282 // Assign lines to clusters/vertices and merge clusters which share 1 line
283 Int_t nClustersAfterMerge=nClusters;
284 Int_t* belongsTo=new Int_t[nLines];
285 for(Int_t ilin=0; ilin<nLines; ilin++) belongsTo[ilin]=-1;
286 for(Int_t iclu=0;iclu<nClusters;iclu++){
287 Int_t actualCluIndex=iclu;
288 for(Int_t ip=0; ip<npoints;ip++){
289 if(mask[ip]==sortedIndex[iclu]){
290 Int_t ind1=index1[ip];
291 if(belongsTo[ind1]==-1) belongsTo[ind1]=actualCluIndex;
292 else if(belongsTo[ind1]<actualCluIndex){
293 Int_t newCluIndex=belongsTo[ind1];
294 for(Int_t ilin=0; ilin<nLines; ilin++){
295 if(belongsTo[ilin]==actualCluIndex) belongsTo[ilin]=newCluIndex;
297 AliDebug(10,Form("Merged cluster %d with %d\n",actualCluIndex,newCluIndex));
298 actualCluIndex=newCluIndex;
299 nClustersAfterMerge--;
301 Int_t ind2=index2[ip];
302 if(belongsTo[ind2]==-1) belongsTo[ind2]=actualCluIndex;
303 else if(belongsTo[ind2]<actualCluIndex){
304 Int_t newCluIndex=belongsTo[ind2];
305 for(Int_t ilin=0; ilin<nLines; ilin++){
306 if(belongsTo[ilin]==actualCluIndex) belongsTo[ilin]=newCluIndex;
308 AliDebug(10,Form("Merged cluster %d with %d\n",actualCluIndex,newCluIndex));
309 actualCluIndex=newCluIndex;
310 nClustersAfterMerge--;
315 AliDebug(3,Form("Number of clusters after merging = %d\n",nClustersAfterMerge));
317 // Count lines associated to each cluster/vertex
318 UInt_t *cluSize=new UInt_t[nClusters];
319 for(Int_t iclu=0;iclu<nClusters;iclu++){
321 for(Int_t ilin=0; ilin<nLines; ilin++){
322 if(belongsTo[ilin]==iclu) cluSize[iclu]++;
326 // Count good vertices (>1 associated tracklet)
328 for(Int_t iclu=0;iclu<nClusters;iclu++){
329 AliDebug(3,Form("Vertex %d Size=%d\n",iclu,cluSize[iclu]));
330 if(cluSize[iclu]>1) nGoodVert++;
333 AliDebug(1,Form("Number of good vertices = %d\n",nGoodVert));
334 // Calculate vertex coordinates for each cluster
336 fVertArray = new AliESDVertex[nGoodVert];
338 for(Int_t iclu=0;iclu<nClusters;iclu++){
339 Int_t size=cluSize[iclu];
341 AliStrLine **arrlin = new AliStrLine*[size];
343 for(Int_t ilin=0; ilin<nLines; ilin++){
344 if(belongsTo[ilin]==iclu){
345 arrlin[nFilled++] = dynamic_cast<AliStrLine*>(fLines[ilin]);
348 AliDebug(3,Form("Vertex %d N associated tracklets = %d out of %d\n",iVert,size,nFilled));
350 fVertArray[iVert]=AliVertexerTracks::TrackletVertexFinder(arrlin,nFilled);
352 fVertArray[iVert].GetXYZ(peak);
353 AliStrLine **arrlin2 = new AliStrLine*[size];
355 for(Int_t i=0; i<nFilled;i++){
356 AliStrLine *l1 = arrlin[i];
357 if(l1->GetDistFromPoint(peak)< fDCAcut)
358 arrlin2[nFilled2++] = dynamic_cast<AliStrLine*>(l1);
361 AliDebug(3,Form("Vertex %d recalculated with %d tracklets\n",iVert,nFilled2));
362 fVertArray[iVert]=AliVertexerTracks::TrackletVertexFinder(arrlin2,nFilled2);
372 fNTrpuv = fVertArray[1].GetNContributors();
373 fZpuv = fVertArray[1].GetZv();
376 Double_t vRadius=TMath::Sqrt(fVertArray[0].GetXv()*fVertArray[0].GetXv()+fVertArray[0].GetYv()*fVertArray[0].GetYv());
377 if(vRadius<GetPipeRadius() && fVertArray[0].GetNContributors()>0){
378 Double_t position[3]={fVertArray[0].GetXv(),fVertArray[0].GetYv(),fVertArray[0].GetZv()};
379 Double_t covmatrix[6];
380 fVertArray[0].GetCovMatrix(covmatrix);
381 Double_t chi2=99999.;
382 Int_t nContr=fVertArray[0].GetNContributors();
383 fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);
384 fCurrentVertex->SetTitle("vertexer: 3D");
385 fCurrentVertex->SetName("SPDVertex3D");
386 fCurrentVertex->SetDispersion(fVertArray[0].GetDispersion());
394 delete [] sortedIndex;
401 //______________________________________________________________________
402 void AliITSVertexer3D::FindVertex3DIterativeMM(){
403 // Defines the AliESDVertex for the current event
404 Int_t numsor=fLines.GetEntriesFast()*(fLines.GetEntriesFast()-1)/2;
405 //cout<<"AliITSVertexer3D::FindVertexForCurentEvent: Number of tracklets selected for vertexing "<<fLines.GetEntriesFast()<<"; Number of pairs: "<<numsor<<endl;
406 AliITSSortTrkl srt(fLines,numsor,fCutOnPairs,fCoarseMaxRCut);
408 AliInfo(Form("Number of vertices: %d",srt.GetNumberOfClusters()));
410 fNoVertices = srt.GetNumberOfClusters();
411 //printf("fNoVertices = %d \n",fNoVertices);
413 fVertArray = new AliESDVertex[fNoVertices];
414 for(Int_t kk=0; kk<srt.GetNumberOfClusters(); kk++){
416 Int_t *labels = srt.GetTrackletsLab(kk,size);
418 Int_t *pairs = srt.GetClusters(kk);
419 Int_t nopai = srt.GetSizeOfCluster(kk);
420 cout<<"***** Vertex number "<<kk<<". Pairs: \n";
421 for(Int_t jj=0;jj<nopai;jj++){
422 cout<<pairs[jj]<<" - ";
423 if(jj>0 & jj%8==0)cout<<endl;
426 cout<<"***** Vertex number "<<kk<<". Labels: \n";
428 AliStrLine **tclo = new AliStrLine* [size];
429 for(Int_t jj=0;jj<size;jj++){
430 // cout<<labels[jj]<<" - ";
431 // if(jj>0 & jj%8==0)cout<<endl;
432 tclo[jj] = dynamic_cast<AliStrLine*>(fLines[labels[jj]]);
436 fVertArray[kk]=AliVertexerTracks::TrackletVertexFinder(tclo,size);
438 // fVertArray[kk].PrintStatus();
440 // at least one second vertex is present
442 fNTrpuv = fVertArray[kk].GetNContributors();
443 fZpuv = fVertArray[kk].GetZv();
446 Double_t vRadius=TMath::Sqrt(fVertArray[0].GetXv()*fVertArray[0].GetXv()+fVertArray[0].GetYv()*fVertArray[0].GetYv());
447 if(vRadius<GetPipeRadius() && fVertArray[0].GetNContributors()>0){
448 Double_t position[3]={fVertArray[0].GetXv(),fVertArray[0].GetYv(),fVertArray[0].GetZv()};
449 Double_t covmatrix[6];
450 fVertArray[0].GetCovMatrix(covmatrix);
451 Double_t chi2=99999.;
452 Int_t nContr=fVertArray[0].GetNContributors();
453 fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);
454 fCurrentVertex->SetTitle("vertexer: 3D");
455 fCurrentVertex->SetName("SPDVertex3D");
456 fCurrentVertex->SetDispersion(fVertArray[0].GetDispersion());
462 //______________________________________________________________________
463 Bool_t AliITSVertexer3D::DistBetweenVertices(AliESDVertex &a, AliESDVertex &b, Double_t test, Double_t &dist){
464 // method to compare the distance between vertices a and b with "test"
465 //it returns kTRUE is the distance is less or equal to test
466 dist = (a.GetX()-b.GetX()) * (a.GetX()-b.GetX());
467 dist += (a.GetY()-b.GetY()) * (a.GetY()-b.GetY());
468 dist += (a.GetZ()-b.GetZ()) * (a.GetZ()-b.GetZ());
469 dist = TMath::Sqrt(dist);
470 if(dist <= test)return kTRUE;
475 //______________________________________________________________________
476 Int_t AliITSVertexer3D::FindTracklets(TTree *itsClusterTree, Int_t optCuts){
477 // All the possible combinations between recpoints on layer 1and 2 are
478 // considered. Straight lines (=tracklets)are formed.
479 // The tracklets are processed in Prepare3DVertex
481 TClonesArray *itsRec = 0;
482 if(optCuts==0) fZHisto->Reset();
483 // gc1 are local and global coordinates for layer 1
484 Float_t gc1f[3]={0.,0.,0.};
485 Double_t gc1[3]={0.,0.,0.};
486 // gc2 are local and global coordinates for layer 2
487 Float_t gc2f[3]={0.,0.,0.};
488 Double_t gc2[3]={0.,0.,0.};
489 AliITSRecPointContainer* rpcont=AliITSRecPointContainer::Instance();
490 itsRec=rpcont->FetchClusters(0,itsClusterTree);
491 if(!rpcont->IsSPDActive()){
492 AliWarning("No SPD rec points found, 3D vertex not calculated");
496 // Set values for cuts
497 Double_t xbeam=GetNominalPos()[0];
498 Double_t ybeam=GetNominalPos()[1];
500 Double_t deltaPhi=fCoarseDiffPhiCut;
501 Double_t deltaR=fCoarseMaxRCut;
502 Double_t dZmax=fZCutDiamond;
504 xbeam=fVert3D.GetXv();
505 ybeam=fVert3D.GetYv();
506 zvert=fVert3D.GetZv();
507 deltaPhi = fDiffPhiMax;
510 if(fPileupAlgo == 2){
514 } else if(optCuts==2){
515 xbeam=fVert3D.GetXv();
516 ybeam=fVert3D.GetYv();
517 deltaPhi = fDiffPhiforPileup;
521 Int_t nrpL1 = 0; // number of rec points on layer 1
522 Int_t nrpL2 = 0; // number of rec points on layer 2
523 nrpL1=rpcont->GetNClustersInLayerFast(1);
524 nrpL2=rpcont->GetNClustersInLayerFast(2);
525 if(nrpL1 == 0 || nrpL2 == 0){
526 AliDebug(1,Form("No RecPoints in at least one SPD layer (%d %d)",nrpL1,nrpL2));
529 AliDebug(1,Form("RecPoints on Layer 1,2 = %d, %d\n",nrpL1,nrpL2));
531 Float_t factDownScal=1.;
532 Int_t origLaddersOnLayer2=fLadOnLay2;
534 if(nrpL1>fMaxNumOfCl || nrpL2>fMaxNumOfCl){
535 SetLaddersOnLayer2(2);
537 factDownScal=(Float_t)fMaxNumOfCl*(Float_t)fMaxNumOfCl/(Float_t)nrpL1/(Float_t)nrpL2;
539 factDownScal*=(fCoarseDiffPhiCut/fDiffPhiMax)*10;
542 SetLaddersOnLayer2(origLaddersOnLayer2);
545 else if(optCuts==2) return -1;
546 if(fDoDownScale)AliDebug(1,Form("Too many recpoints on SPD(%d %d ), downscale by %f",nrpL1,nrpL2,factDownScal));
549 Double_t a[3]={xbeam,ybeam,0.};
550 Double_t b[3]={xbeam,ybeam,10.};
551 AliStrLine zeta(a,b,kTRUE);
552 static Double_t bField=TMath::Abs(AliTracker::GetBz()/10.); //T
553 SetMeanPPtSelTracks(bField);
556 // Loop on modules of layer 1
557 Int_t firstL1 = TMath::Max(0,AliITSgeomTGeo::GetModuleIndex(1,1,1));
558 Int_t lastL1 = AliITSgeomTGeo::GetModuleIndex(2,1,1)-1;
559 for(Int_t modul1= firstL1; modul1<=lastL1;modul1++){ // Loop on modules of layer 1
560 if(!fUseModule[modul1]) continue;
562 UShort_t ladder=modul1/4+1; // ladders are numbered starting from 1
563 TClonesArray *prpl1=rpcont->UncheckedGetClusters(modul1);
564 Int_t nrecp1 = prpl1->GetEntries();
565 for(Int_t j=0;j<nrecp1;j++){
566 if(j>kMaxCluPerMod) continue;
567 UShort_t idClu1=modul1*kMaxCluPerMod+j;
568 if(fUsedCluster.TestBitNumber(idClu1)) continue;
569 if(fDoDownScale && fGenerForDownScale->Rndm()>factDownScal) continue;
570 AliITSRecPoint *recp1 = (AliITSRecPoint*)prpl1->At(j);
571 recp1->GetGlobalXYZ(gc1f);
572 for(Int_t ico=0;ico<3;ico++)gc1[ico]=gc1f[ico];
574 Double_t phi1 = TMath::ATan2(gc1[1]-ybeam,gc1[0]-xbeam);
575 if(phi1<0)phi1=2*TMath::Pi()+phi1;
576 for(Int_t ladl2=0 ; ladl2<fLadOnLay2*2+1;ladl2++){
577 for(Int_t k=0;k<4;k++){
578 Int_t ladmod=fLadders[ladder-1]+ladl2;
579 if(ladmod>AliITSgeomTGeo::GetNLadders(2)) ladmod=ladmod-AliITSgeomTGeo::GetNLadders(2);
580 Int_t modul2=AliITSgeomTGeo::GetModuleIndex(2,ladmod,k+1);
581 if(modul2<0)continue;
582 if(!fUseModule[modul2]) continue;
583 itsRec=rpcont->UncheckedGetClusters(modul2);
584 Int_t nrecp2 = itsRec->GetEntries();
585 for(Int_t j2=0;j2<nrecp2;j2++){
586 if(j2>kMaxCluPerMod) continue;
587 UShort_t idClu2=modul2*kMaxCluPerMod+j2;
588 if(fUsedCluster.TestBitNumber(idClu2)) continue;
590 AliITSRecPoint *recp2 = (AliITSRecPoint*)itsRec->At(j2);
591 recp2->GetGlobalXYZ(gc2f);
592 for(Int_t ico=0;ico<3;ico++)gc2[ico]=gc2f[ico];
593 Double_t phi2 = TMath::ATan2(gc2[1]-ybeam,gc2[0]-xbeam);
594 if(phi2<0)phi2=2*TMath::Pi()+phi2;
595 Double_t diff = TMath::Abs(phi2-phi1);
596 if(diff>TMath::Pi())diff=2.*TMath::Pi()-diff;
597 if(optCuts==0 && diff<fDiffPhiforPileup){
598 Double_t r1=TMath::Sqrt(gc1[0]*gc1[0]+gc1[1]*gc1[1]);
600 Double_t r2=TMath::Sqrt(gc2[0]*gc2[0]+gc2[1]*gc2[1]);
602 Double_t zr0=(r2*zc1-r1*zc2)/(r2-r1); //Z @ null radius
605 if(diff>deltaPhi)continue;
606 AliStrLine line(gc1,gc2,kTRUE);
608 Int_t retcode = line.Cross(&zeta,cp);
609 if(retcode<0)continue;
610 Double_t dca = line.GetDCA(&zeta);
612 if(dca>deltaR)continue;
613 Double_t deltaZ=cp[2]-zvert;
614 if(TMath::Abs(deltaZ)>dZmax)continue;
618 if(fLines.GetEntriesFast()>0)fLines.Clear("C");
621 recp2->GetGlobalCov(cov);
624 Double_t rad1=TMath::Sqrt(gc1[0]*gc1[0]+gc1[1]*gc1[1]);
625 Double_t rad2=TMath::Sqrt(gc2[0]*gc2[0]+gc2[1]*gc2[1]);
626 Double_t factor=(rad1+rad2)/(rad2-rad1); //factor to account for error on tracklet direction
630 Double_t curvRadius=fMeanPtSelTrk/(0.3*bField)*100; //cm
631 Double_t dRad=TMath::Sqrt((gc1[0]-gc2[0])*(gc1[0]-gc2[0])+(gc1[1]-gc2[1])*(gc1[1]-gc2[1]));
632 Double_t aux=dRad/2.+rad1;
633 curvErr=TMath::Sqrt(curvRadius*curvRadius-dRad*dRad/4.)-TMath::Sqrt(curvRadius*curvRadius-aux*aux); //cm
636 sigmasq[0]=(cov[0]+curvErr*curvErr/2.)*factor*factor;
637 sigmasq[1]=(cov[3]+curvErr*curvErr/2.)*factor*factor;
638 sigmasq[2]=cov[5]*factor*factor;
640 // Multiple scattering
641 Double_t pOverMass=fMeanPSelTrk/0.140;
642 Double_t beta2=pOverMass*pOverMass/(1+pOverMass*pOverMass);
643 Double_t p2=fMeanPSelTrk*fMeanPSelTrk;
644 Double_t rBP=GetPipeRadius();
645 Double_t dBP=0.08/35.3; // 800 um of Be
646 Double_t dL1=0.01; //approx. 1% of radiation length
647 Double_t theta2BP=14.1*14.1/(beta2*p2*1e6)*dBP;
648 Double_t theta2L1=14.1*14.1/(beta2*p2*1e6)*dL1;
649 Double_t rtantheta1=(rad2-rad1)*TMath::Tan(TMath::Sqrt(theta2L1));
650 Double_t rtanthetaBP=(rad1-rBP)*TMath::Tan(TMath::Sqrt(theta2BP));
651 for(Int_t ico=0; ico<3;ico++){
652 sigmasq[ico]+=rtantheta1*rtantheta1*factor*factor/3.;
653 sigmasq[ico]+=rtanthetaBP*rtanthetaBP*factor*factor/3.;
655 Double_t wmat[9]={1.,0.,0.,0.,1.,0.,0.,0.,1.};
656 if(sigmasq[0]!=0.) wmat[0]=1./sigmasq[0];
657 if(sigmasq[1]!=0.) wmat[4]=1./sigmasq[1];
658 if(sigmasq[2]!=0.) wmat[8]=1./sigmasq[2];
659 new(fLines[nolines++])AliStrLine(gc1,sigmasq,wmat,gc2,kTRUE,idClu1,idClu2);
668 SetLaddersOnLayer2(origLaddersOnLayer2);
672 if(nolines == 0)return -2;
676 //______________________________________________________________________
677 Int_t AliITSVertexer3D::Prepare3DVertex(Int_t optCuts){
678 // Finds the 3D vertex information using tracklets
681 Double_t xbeam=GetNominalPos()[0];
682 Double_t ybeam=GetNominalPos()[1];
684 Double_t deltaR=fCoarseMaxRCut;
685 Double_t dZmax=fZCutDiamond;
687 xbeam=fVert3D.GetXv();
688 ybeam=fVert3D.GetYv();
689 zvert=fVert3D.GetZv();
692 if(fPileupAlgo == 2){
696 }else if(optCuts==2){
697 xbeam=fVert3D.GetXv();
698 ybeam=fVert3D.GetYv();
702 Double_t origBinSizeR=fBinSizeR;
703 Double_t origBinSizeZ=fBinSizeZ;
709 Double_t rl=-fCoarseMaxRCut;
710 Double_t rh=fCoarseMaxRCut;
711 Double_t zl=-fZCutDiamond;
712 Double_t zh=fZCutDiamond;
713 Int_t nbr=(Int_t)((rh-rl)/fBinSizeR+0.0001);
714 Int_t nbz=(Int_t)((zh-zl)/fBinSizeZ+0.0001);
715 Int_t nbrcs=(Int_t)((rh-rl)/(fBinSizeR*2.)+0.0001);
716 Int_t nbzcs=(Int_t)((zh-zl)/(fBinSizeZ*2.)+0.0001);
718 TH3F *h3d = new TH3F("h3d","xyz distribution",nbr,rl,rh,nbr,rl,rh,nbz,zl,zh);
719 TH3F *h3dcs = new TH3F("h3dcs","xyz distribution",nbrcs,rl,rh,nbrcs,rl,rh,nbzcs,zl,zh);
721 // cleanup of the TCLonesArray of tracklets (i.e. fakes are removed)
722 Int_t *validate = new Int_t [fLines.GetEntriesFast()];
723 for(Int_t i=0; i<fLines.GetEntriesFast();i++)validate[i]=0;
724 for(Int_t i=0; i<fLines.GetEntriesFast()-1;i++){
725 AliStrLine *l1 = (AliStrLine*)fLines.At(i);
726 for(Int_t j=i+1;j<fLines.GetEntriesFast();j++){
727 AliStrLine *l2 = (AliStrLine*)fLines.At(j);
728 Double_t dca=l1->GetDCA(l2);
729 if(dca > fDCAcut || dca<0.00001) continue;
731 Int_t retc = l1->Cross(l2,point);
733 Double_t deltaZ=point[2]-zvert;
734 if(TMath::Abs(deltaZ)>dZmax)continue;
735 Double_t rad=TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
736 if(rad>fCoarseMaxRCut)continue;
737 Double_t deltaX=point[0]-xbeam;
738 Double_t deltaY=point[1]-ybeam;
739 Double_t raddist=TMath::Sqrt(deltaX*deltaX+deltaY*deltaY);
740 if(raddist>deltaR)continue;
743 h3d->Fill(point[0],point[1],point[2]);
744 h3dcs->Fill(point[0],point[1],point[2]);
750 Int_t numbtracklets=0;
751 for(Int_t i=0; i<fLines.GetEntriesFast();i++)if(validate[i]>=1)numbtracklets++;
759 for(Int_t i=0; i<fLines.GetEntriesFast();i++){
760 if(validate[i]<1)fLines.RemoveAt(i);
763 AliDebug(1,Form("Number of tracklets (after compress)%d ",fLines.GetEntriesFast()));
766 // Exit here if Pileup Algorithm 2 has been chosen during second loop
767 if(fPileupAlgo == 2 && optCuts==1){
773 // Find peaks in histos
775 Double_t peak[3]={0.,0.,0.};
777 FindPeaks(h3d,peak,ntrkl,ntimes);
779 Double_t binsizer=(rh-rl)/nbr;
780 Double_t binsizez=(zh-zl)/nbz;
781 if(optCuts==0 && (ntrkl<=2 || ntimes>1)){
784 FindPeaks(h3dcs,peak,ntrkl,ntimes);
785 binsizer=(rh-rl)/nbrcs;
786 binsizez=(zh-zl)/nbzcs;
787 if(ntrkl==1 || ntimes>1){delete h3dcs; return retcode;}
791 // Finer Histo in limited range in case of high mult.
795 Double_t xl=peak[0]-0.3;
796 Double_t xh=peak[0]+0.3;
797 Double_t yl=peak[1]-0.3;
798 Double_t yh=peak[1]+0.3;
801 Int_t nbxfs=(Int_t)((xh-xl)/fBinSizeR+0.0001);
802 Int_t nbyfs=(Int_t)((yh-yl)/fBinSizeR+0.0001);
803 Int_t nbzfs=(Int_t)((zh-zl)/fBinSizeZ+0.0001);
805 TH3F *h3dfs = new TH3F("h3dfs","xyz distribution",nbxfs,xl,xh,nbyfs,yl,yh,nbzfs,zl,zh);
806 for(Int_t i=0; i<fLines.GetEntriesFast()-1;i++){
807 AliStrLine *l1 = (AliStrLine*)fLines.At(i);
808 for(Int_t j=i+1;j<fLines.GetEntriesFast();j++){
809 AliStrLine *l2 = (AliStrLine*)fLines.At(j);
811 Int_t retc = l1->Cross(l2,point);
813 h3dfs->Fill(point[0],point[1],point[2]);
819 Double_t newpeak[3]={0.,0.,0.};
820 FindPeaks(h3dfs,newpeak,ntrkl,ntimes);
822 for(Int_t iCoo=0; iCoo<3; iCoo++) peak[iCoo]=newpeak[iCoo];
827 SetBinSizeR(origBinSizeR);
828 SetBinSizeZ(origBinSizeZ);
832 // Second selection loop
834 Double_t bs=(binsizer+binsizez)/2.;
835 for(Int_t i=0; i<fLines.GetEntriesFast();i++){
836 AliStrLine *l1 = (AliStrLine*)fLines.At(i);
837 if(l1->GetDistFromPoint(peak)>2.5*bs)fLines.RemoveAt(i);
840 AliDebug(1,Form("Number of tracklets (after 2nd compression) %d",fLines.GetEntriesFast()));
842 if(fLines.GetEntriesFast()>1){
844 // find a first candidate for the primary vertex
845 fVert3D=AliVertexerTracks::TrackletVertexFinder(&fLines,0);
846 // make a further selection on tracklets based on this first candidate
847 fVert3D.GetXYZ(peak);
848 AliDebug(1,Form("FIRST V candidate: %f ; %f ; %f",peak[0],peak[1],peak[2]));
849 Int_t *validate2 = new Int_t [fLines.GetEntriesFast()];
850 for(Int_t i=0; i<fLines.GetEntriesFast();i++) validate2[i]=1;
851 for(Int_t i=0; i<fLines.GetEntriesFast();i++){
852 if(validate2[i]==0) continue;
853 AliStrLine *l1 = (AliStrLine*)fLines.At(i);
854 if(l1->GetDistFromPoint(peak)> fDCAcut)fLines.RemoveAt(i);
855 if(optCuts==2){ // temporarily only for pileup
856 for(Int_t j=i+1; j<fLines.GetEntriesFast();j++){
857 AliStrLine *l2 = (AliStrLine*)fLines.At(j);
858 if(l1->GetDCA(l2)<0.00001){
859 Int_t delta1=(Int_t)l1->GetIdPoint(0)-(Int_t)l2->GetIdPoint(0);
860 Int_t delta2=(Int_t)l1->GetIdPoint(1)-(Int_t)l2->GetIdPoint(1);
861 Int_t deltamod1=(Int_t)l1->GetIdPoint(0)/kMaxCluPerMod
862 -(Int_t)l2->GetIdPoint(0)/kMaxCluPerMod;
863 Int_t deltamod2=(Int_t)l1->GetIdPoint(1)/kMaxCluPerMod
864 -(Int_t)l2->GetIdPoint(1)/kMaxCluPerMod;
865 // remove tracklets sharing a point
866 if( (delta1==0 && deltamod2==0) ||
867 (delta2==0 && deltamod1==0) ) validate2[j]=0;
872 for(Int_t i=0; i<fLines.GetEntriesFast();i++){
873 if(validate2[i]==0) fLines.RemoveAt(i);
877 AliDebug(1,Form("Number of tracklets (after 3rd compression) %d",fLines.GetEntriesFast()));
878 if(fLines.GetEntriesFast()>1){// this new tracklet selection is used
879 fVert3D=AliVertexerTracks::TrackletVertexFinder(&fLines,0);
885 //________________________________________________________
886 void AliITSVertexer3D::SetMeanPPtSelTracks(Double_t fieldTesla){
887 // Sets mean values of Pt based on the field
888 // for P (used in multiple scattering) the most probable value is used
889 if(TMath::Abs(fieldTesla-0.5)<0.01){
890 SetMeanPSelTracks(0.375);
891 SetMeanPtSelTracks(0.630);
892 }else if(TMath::Abs(fieldTesla-0.4)<0.01){
893 SetMeanPSelTracks(0.375);
894 SetMeanPtSelTracks(0.580);
895 }else if(TMath::Abs(fieldTesla-0.2)<0.01){
896 SetMeanPSelTracks(0.375);
897 SetMeanPtSelTracks(0.530);
898 }else if(fieldTesla<0.00001){
899 SetMeanPSelTracks(0.375);
900 SetMeanPtSelTracks(0.230);
903 SetMeanPtSelTracks();
907 //________________________________________________________
908 void AliITSVertexer3D::FindPeaks(TH3F* histo, Double_t *peak, Int_t &nOfTracklets, Int_t &nOfTimes){
909 // Finds bin with max contents in 3D histo of tracket intersections
910 TAxis *xax = histo->GetXaxis();
911 TAxis *yax = histo->GetYaxis();
912 TAxis *zax = histo->GetZaxis();
918 Int_t peakbin[3]={0,0,0};
919 Int_t peak2bin[3]={-1,-1,-1};
921 for(Int_t i=xax->GetFirst();i<=xax->GetLast();i++){
922 Double_t xval = xax->GetBinCenter(i);
923 for(Int_t j=yax->GetFirst();j<=yax->GetLast();j++){
924 Double_t yval = yax->GetBinCenter(j);
925 for(Int_t k=zax->GetFirst();k<=zax->GetLast();k++){
926 Double_t zval = zax->GetBinCenter(k);
927 Int_t bc =(Int_t)histo->GetBinContent(i,j,k);
942 }else if(bc==nOfTracklets){
943 if(TMath::Abs(i-peakbin[0])<=1 && TMath::Abs(j-peakbin[1])<=1 && TMath::Abs(k-peakbin[2])<=1){
956 if(peak2bin[0]>=-1 && bc2!=-1){ // two contiguous peak-cells with same contents
957 peak[0]=0.5*(xax->GetBinCenter(peakbin[0])+xax->GetBinCenter(peak2bin[0]));
958 peak[1]=0.5*(yax->GetBinCenter(peakbin[1])+yax->GetBinCenter(peak2bin[1]));
959 peak[2]=0.5*(zax->GetBinCenter(peakbin[2])+zax->GetBinCenter(peak2bin[2]));
964 //________________________________________________________
965 void AliITSVertexer3D::MarkUsedClusters(){
966 // Mark clusters of tracklets used in vertex claulation
967 for(Int_t i=0; i<fLines.GetEntriesFast();i++){
968 AliStrLine *lin = (AliStrLine*)fLines.At(i);
969 Int_t idClu1=lin->GetIdPoint(0);
970 Int_t idClu2=lin->GetIdPoint(1);
971 fUsedCluster.SetBitNumber(idClu1);
972 fUsedCluster.SetBitNumber(idClu2);
975 //________________________________________________________
976 Int_t AliITSVertexer3D::RemoveTracklets(){
977 // Remove trackelts close to first found vertex
978 Double_t vert[3]={fVert3D.GetXv(),fVert3D.GetYv(),fVert3D.GetZv()};
980 for(Int_t i=0; i<fLines.GetEntriesFast();i++){
981 AliStrLine *lin = (AliStrLine*)fLines.At(i);
982 if(lin->GetDistFromPoint(vert)<fDCAforPileup){
983 Int_t idClu1=lin->GetIdPoint(0);
984 Int_t idClu2=lin->GetIdPoint(1);
985 fUsedCluster.SetBitNumber(idClu1);
986 fUsedCluster.SetBitNumber(idClu2);
994 //________________________________________________________
995 void AliITSVertexer3D::FindOther3DVertices(TTree *itsClusterTree){
996 // pileup identification based on 3D vertexing with not used clusters
999 Int_t nolines = FindTracklets(itsClusterTree,2);
1001 Int_t nr=RemoveTracklets();
1004 Int_t rc=Prepare3DVertex(2);
1006 fVert3D=AliVertexerTracks::TrackletVertexFinder(&fLines,0);
1007 if(fVert3D.GetNContributors()>=fMinTrackletsForPilup){
1010 fVertArray = new AliESDVertex[2];
1011 fVertArray[0]=(*fCurrentVertex);
1012 fVertArray[1]=fVert3D;
1013 fZpuv=fVert3D.GetZv();
1014 fNTrpuv=fVert3D.GetNContributors();
1020 //______________________________________________________________________
1021 void AliITSVertexer3D::PileupFromZ(){
1022 // Calls the pileup algorithm of ALiITSVertexerZ
1023 Int_t binmin, binmax;
1024 Int_t nPeaks=AliITSVertexerZ::GetPeakRegion(fZHisto,binmin,binmax);
1025 if(nPeaks==2)AliWarning("2 peaks found");
1026 Int_t firstPeakCont=0;
1027 Double_t firstPeakPos=0.;
1028 for(Int_t i=binmin-1;i<=binmax+1;i++){
1029 firstPeakCont+=static_cast<Int_t>(fZHisto->GetBinContent(i));
1030 firstPeakPos+=fZHisto->GetBinContent(i)*fZHisto->GetBinCenter(i);
1032 if(firstPeakCont>0){
1033 firstPeakPos/=firstPeakCont;
1035 if(firstPeakCont>fMinTrackletsForPilup){
1037 ncontr2=AliITSVertexerZ::FindSecondPeak(fZHisto,binmin,binmax,secPeakPos);
1038 if(ncontr2>=fMinTrackletsForPilup){
1041 AliESDVertex secondVert(secPeakPos,0.1,ncontr2);
1042 fVertArray = new AliESDVertex[2];
1043 fVertArray[0]=(*fCurrentVertex);
1044 fVertArray[1]=secondVert;
1051 //________________________________________________________
1052 void AliITSVertexer3D::PrintStatus() const {
1053 // Print current status
1054 printf("========= First step selections =====================\n");
1055 printf("Cut on diamond (Z) %f\n",fZCutDiamond);
1056 printf("Loose cut on Delta Phi %f\n",fCoarseDiffPhiCut);
1057 printf("Loose cut on tracklet DCA to Z axis %f\n",fCoarseMaxRCut);
1058 printf("Cut on DCA - tracklet to tracklet and to vertex %f\n",fDCAcut);
1059 printf("========= Second step selections ====================\n");
1060 printf("Cut on tracklet-to-first-vertex Z distance %f\n",fMaxZCut);
1061 printf("Max Phi difference: %f\n",fDiffPhiMax);
1062 printf("Cut on tracklet DCA to beam axis %f\n",fMaxRCut);
1063 printf("Cut on tracklet DCA to beam axis (algo2) %f\n",fMaxRCut2);
1064 printf("========= Pileup selections =========================\n");
1065 printf("Pileup algo: %d\n",fPileupAlgo);
1066 printf("Min DCA to 1st vertex for pileup (algo 0 and 1): %f\n",fDCAforPileup);
1067 printf("Cut on distance between pair-vertices (algo 2): %f\n",fCutOnPairs);
1068 printf("Maximum number of clusters allowed on L1 or L2: %d\n",fMaxNumOfCl);
1069 printf("=======================================================\n");