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 "AliESDVertex.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
31 // It can be used successfully with Pb-Pb collisions
32 ////////////////////////////////////////////////////////////////
34 ClassImp(AliITSVertexer3D)
38 //______________________________________________________________________
39 AliITSVertexer3D::AliITSVertexer3D():
41 fLines("AliStrLine",1000),
43 fCoarseDiffPhiCut(0.),
54 fUsedCluster(kMaxCluPerMod*kNSPDMod),
61 // Default constructor
62 SetCoarseDiffPhiCut();
77 Float_t binsize=0.02; // default 200 micron
78 Int_t nbins=static_cast<Int_t>(1+2*fZCutDiamond/binsize);
79 fZHisto=new TH1F("hz","",nbins,-fZCutDiamond,-fZCutDiamond+binsize*nbins);
82 //______________________________________________________________________
83 AliITSVertexer3D::~AliITSVertexer3D() {
86 if(fZHisto) delete fZHisto;
89 //______________________________________________________________________
90 void AliITSVertexer3D::ResetVert3D(){
96 fVert3D.SetDispersion(0.);
97 fVert3D.SetNContributors(0);
98 fUsedCluster.ResetAllBits(0);
100 //______________________________________________________________________
101 AliESDVertex* AliITSVertexer3D::FindVertexForCurrentEvent(TTree *itsClusterTree){
102 // Defines the AliESDVertex for the current event
104 AliDebug(1,"FindVertexForCurrentEvent - 3D - PROCESSING NEXT EVENT");
106 fCurrentVertex = NULL;
108 Int_t nolines = FindTracklets(itsClusterTree,0);
110 Int_t rc=Prepare3DVertex(0);
111 if(fPileupAlgo == 2 && rc == 0) FindVertex3DIterative();
112 else if(fPileupAlgo<2 && rc == 0) FindVertex3D(itsClusterTree);
116 AliITSVertexerZ vertz(GetNominalPos()[0],GetNominalPos()[1]);
117 vertz.SetDetTypeRec(GetDetTypeRec());
118 AliDebug(1,"Call Vertexer Z\n");
119 vertz.SetLowLimit(-fZCutDiamond);
120 vertz.SetHighLimit(fZCutDiamond);
121 AliESDVertex* vtxz = vertz.FindVertexForCurrentEvent(itsClusterTree);
123 Double_t position[3]={GetNominalPos()[0],GetNominalPos()[1],vtxz->GetZv()};
124 Double_t covmatrix[6];
125 vtxz->GetCovMatrix(covmatrix);
126 Double_t chi2=99999.;
127 Int_t nContr=vtxz->GetNContributors();
128 fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);
129 fCurrentVertex->SetTitle("vertexer: Z");
130 fCurrentVertex->SetName("SPDVertexZ");
135 FindMultiplicity(itsClusterTree);
136 return fCurrentVertex;
139 //______________________________________________________________________
140 void AliITSVertexer3D::FindVertex3D(TTree *itsClusterTree){
142 /* uncomment to debug
143 printf("Vertex found in first iteration:\n");
145 printf("Start second iteration\n");
146 end of debug lines */
147 if(fVert3D.GetNContributors()>0){
149 Int_t nolines = FindTracklets(itsClusterTree,1);
151 Int_t rc=Prepare3DVertex(1);
152 if(rc!=0) fVert3D.SetNContributors(0); // exclude this vertex
155 /* uncomment to debug
156 printf("Vertex found in second iteration:\n");
158 end of debug lines */
160 Float_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());
174 case 0: PileupFromZ(); break;
175 case 1: FindOther3DVertices(itsClusterTree); break;
176 default: AliError("Wrong pileup algorithm"); break;
179 fVertArray = new AliESDVertex[1];
180 fVertArray[0]=(*fCurrentVertex);
185 //______________________________________________________________________
186 void AliITSVertexer3D::FindVertex3DIterative(){
187 // Defines the AliESDVertex for the current event
188 Int_t numsor=fLines.GetEntriesFast()*(fLines.GetEntriesFast()-1)/2;
189 //cout<<"AliITSVertexer3D::FindVertexForCurentEvent: Number of tracklets selected for vertexing "<<fLines.GetEntriesFast()<<"; Number of pairs: "<<numsor<<endl;
190 AliITSSortTrkl srt(fLines,numsor,fCutOnPairs,fCoarseMaxRCut);
192 AliInfo(Form("Number of vertices: %d",srt.GetNumberOfClusters()));
194 fNoVertices = srt.GetNumberOfClusters();
195 //printf("fNoVertices = %d \n",fNoVertices);
197 fVertArray = new AliESDVertex[fNoVertices];
198 for(Int_t kk=0; kk<srt.GetNumberOfClusters(); kk++){
200 Int_t *labels = srt.GetTrackletsLab(kk,size);
202 Int_t *pairs = srt.GetClusters(kk);
203 Int_t nopai = srt.GetSizeOfCluster(kk);
204 cout<<"***** Vertex number "<<kk<<". Pairs: \n";
205 for(Int_t jj=0;jj<nopai;jj++){
206 cout<<pairs[jj]<<" - ";
207 if(jj>0 & jj%8==0)cout<<endl;
210 cout<<"***** Vertex number "<<kk<<". Labels: \n";
212 AliStrLine **tclo = new AliStrLine* [size];
213 for(Int_t jj=0;jj<size;jj++){
214 // cout<<labels[jj]<<" - ";
215 // if(jj>0 & jj%8==0)cout<<endl;
216 tclo[jj] = dynamic_cast<AliStrLine*>(fLines[labels[jj]]);
220 fVertArray[kk]=AliVertexerTracks::TrackletVertexFinder(tclo,size);
222 // fVertArray[kk].PrintStatus();
224 // at least one second vertex is present
226 fNTrpuv = fVertArray[kk].GetNContributors();
227 fZpuv = fVertArray[kk].GetZv();
230 Float_t vRadius=TMath::Sqrt(fVertArray[0].GetXv()*fVertArray[0].GetXv()+fVertArray[0].GetYv()*fVertArray[0].GetYv());
231 if(vRadius<GetPipeRadius() && fVertArray[0].GetNContributors()>0){
232 Double_t position[3]={fVertArray[0].GetXv(),fVertArray[0].GetYv(),fVertArray[0].GetZv()};
233 Double_t covmatrix[6];
234 fVertArray[0].GetCovMatrix(covmatrix);
235 Double_t chi2=99999.;
236 Int_t nContr=fVertArray[0].GetNContributors();
237 fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);
238 fCurrentVertex->SetTitle("vertexer: 3D");
239 fCurrentVertex->SetName("SPDVertex3D");
240 fCurrentVertex->SetDispersion(fVertArray[0].GetDispersion());
246 //______________________________________________________________________
247 Bool_t AliITSVertexer3D::DistBetweenVertices(AliESDVertex &a, AliESDVertex &b, Double_t test, Double_t &dist){
248 // method to compare the distance between vertices a and b with "test"
249 //it returns kTRUE is the distance is less or equal to test
250 dist = (a.GetX()-b.GetX()) * (a.GetX()-b.GetX());
251 dist += (a.GetY()-b.GetY()) * (a.GetY()-b.GetY());
252 dist += (a.GetZ()-b.GetZ()) * (a.GetZ()-b.GetZ());
253 dist = TMath::Sqrt(dist);
254 if(dist <= test)return kTRUE;
259 //______________________________________________________________________
260 Int_t AliITSVertexer3D::FindTracklets(TTree *itsClusterTree, Int_t optCuts){
261 // All the possible combinations between recpoints on layer 1and 2 are
262 // considered. Straight lines (=tracklets)are formed.
263 // The tracklets are processed in Prepare3DVertex
265 if(!GetDetTypeRec())AliFatal("DetTypeRec pointer has not been set");
267 TTree *tR = itsClusterTree;
268 fDetTypeRec->ResetRecPoints();
269 fDetTypeRec->SetTreeAddressR(tR);
270 TClonesArray *itsRec = 0;
271 if(optCuts==0) fZHisto->Reset();
272 // gc1 are local and global coordinates for layer 1
273 Float_t gc1[3]={0.,0.,0.};
274 // gc2 are local and global coordinates for layer 2
275 Float_t gc2[3]={0.,0.,0.};
277 itsRec = fDetTypeRec->RecPoints();
278 TBranch *branch = NULL;
279 branch = tR->GetBranch("ITSRecPoints");
281 AliError("Null pointer for RecPoints branch");
285 // Set values for cuts
286 Float_t xbeam=GetNominalPos()[0];
287 Float_t ybeam=GetNominalPos()[1];
289 Float_t deltaPhi=fCoarseDiffPhiCut;
290 Float_t deltaR=fCoarseMaxRCut;
291 Float_t dZmax=fZCutDiamond;
292 if(fPileupAlgo == 2){
293 deltaPhi=fFineDiffPhiCut;
295 if(optCuts != 0)AliWarning(Form("fPileupAlgo=2 AND optCuts=%d has been selected. It should be 0",optCuts));
296 } else if(optCuts==1){
297 xbeam=fVert3D.GetXv();
298 ybeam=fVert3D.GetYv();
299 zvert=fVert3D.GetZv();
300 deltaPhi = fDiffPhiMax;
303 } else if(optCuts==2){
304 xbeam=fVert3D.GetXv();
305 ybeam=fVert3D.GetYv();
306 deltaPhi = fDiffPhiMax;
310 Int_t nrpL1 = 0; // number of rec points on layer 1
311 Int_t nrpL2 = 0; // number of rec points on layer 2
313 // By default irstL1=0 and lastL1=79
314 Int_t firstL1 = AliITSgeomTGeo::GetModuleIndex(1,1,1);
315 Int_t lastL1 = AliITSgeomTGeo::GetModuleIndex(2,1,1)-1;
316 for(Int_t module= firstL1; module<=lastL1;module++){ // count number of recopints on layer 1
317 branch->GetEvent(module);
318 nrpL1+= itsRec->GetEntries();
319 fDetTypeRec->ResetRecPoints();
321 //By default firstL2=80 and lastL2=239
322 Int_t firstL2 = AliITSgeomTGeo::GetModuleIndex(2,1,1);
323 Int_t lastL2 = AliITSgeomTGeo::GetModuleIndex(3,1,1)-1;
324 for(Int_t module= firstL2; module<=lastL2;module++){ // count number of recopints on layer 2
325 branch->GetEvent(module);
326 nrpL2+= itsRec->GetEntries();
327 fDetTypeRec->ResetRecPoints();
329 if(nrpL1 == 0 || nrpL2 == 0){
332 AliDebug(1,Form("RecPoints on Layer 1,2 = %d, %d\n",nrpL1,nrpL2));
334 Double_t a[3]={xbeam,ybeam,0.};
335 Double_t b[3]={xbeam,ybeam,10.};
336 AliStrLine zeta(a,b,kTRUE);
337 Float_t bField=AliTracker::GetBz()/10.; //T
338 SetMeanPPtSelTracks(bField);
341 // Loop on modules of layer 1
342 for(Int_t modul1= firstL1; modul1<=lastL1;modul1++){ // Loop on modules of layer 1
343 if(!fUseModule[modul1]) continue;
344 UShort_t ladder=int(modul1/4)+1; // ladders are numbered starting from 1
345 branch->GetEvent(modul1);
346 Int_t nrecp1 = itsRec->GetEntries();
347 static TClonesArray prpl1("AliITSRecPoint",nrecp1);
349 for(Int_t j=0;j<nrecp1;j++){
350 AliITSRecPoint *recp = (AliITSRecPoint*)itsRec->At(j);
351 new(prpl1[j])AliITSRecPoint(*recp);
353 fDetTypeRec->ResetRecPoints();
354 for(Int_t j=0;j<nrecp1;j++){
355 if(j>kMaxCluPerMod) continue;
356 UShort_t idClu1=modul1*kMaxCluPerMod+j;
357 if(fUsedCluster.TestBitNumber(idClu1)) continue;
358 AliITSRecPoint *recp1 = (AliITSRecPoint*)prpl1.At(j);
359 recp1->GetGlobalXYZ(gc1);
360 Double_t phi1 = TMath::ATan2(gc1[1]-ybeam,gc1[0]-xbeam);
361 if(phi1<0)phi1=2*TMath::Pi()+phi1;
362 for(Int_t ladl2=0 ; ladl2<fLadOnLay2*2+1;ladl2++){
363 for(Int_t k=0;k<4;k++){
364 Int_t ladmod=fLadders[ladder-1]+ladl2;
365 if(ladmod>AliITSgeomTGeo::GetNLadders(2)) ladmod=ladmod-AliITSgeomTGeo::GetNLadders(2);
366 Int_t modul2=AliITSgeomTGeo::GetModuleIndex(2,ladmod,k+1);
367 if(!fUseModule[modul2]) continue;
368 branch->GetEvent(modul2);
369 Int_t nrecp2 = itsRec->GetEntries();
370 for(Int_t j2=0;j2<nrecp2;j2++){
371 if(j2>kMaxCluPerMod) continue;
372 UShort_t idClu2=modul2*kMaxCluPerMod+j2;
373 if(fUsedCluster.TestBitNumber(idClu2)) continue;
374 AliITSRecPoint *recp2 = (AliITSRecPoint*)itsRec->At(j2);
375 recp2->GetGlobalXYZ(gc2);
376 Double_t phi2 = TMath::ATan2(gc2[1]-ybeam,gc2[0]-xbeam);
377 if(phi2<0)phi2=2*TMath::Pi()+phi2;
378 Double_t diff = TMath::Abs(phi2-phi1);
379 if(diff>TMath::Pi())diff=2.*TMath::Pi()-diff;
380 if(optCuts==0 && diff<fDiffPhiMax){
381 Float_t r1=TMath::Sqrt(gc1[0]*gc1[0]+gc1[1]*gc1[1]);
383 Float_t r2=TMath::Sqrt(gc2[0]*gc2[0]+gc2[1]*gc2[1]);
385 Float_t zr0=(r2*zc1-r1*zc2)/(r2-r1); //Z @ null radius
388 if(diff>deltaPhi)continue;
389 AliStrLine line(gc1,gc2,kTRUE);
391 Int_t retcode = line.Cross(&zeta,cp);
392 if(retcode<0)continue;
393 Double_t dca = line.GetDCA(&zeta);
395 if(dca>deltaR)continue;
396 Double_t deltaZ=cp[2]-zvert;
397 if(TMath::Abs(deltaZ)>dZmax)continue;
400 if(fLines.GetEntriesFast()>0)fLines.Clear("C");
403 recp2->GetGlobalCov(cov);
406 Float_t rad1=TMath::Sqrt(gc1[0]*gc1[0]+gc1[1]*gc1[1]);
407 Float_t rad2=TMath::Sqrt(gc2[0]*gc2[0]+gc2[1]*gc2[1]);
408 Float_t factor=(rad1+rad2)/(rad2-rad1); //factor to account for error on tracklet direction
412 Float_t curvRadius=fMeanPtSelTrk/(0.3*bField)*100; //cm
413 Float_t dRad=TMath::Sqrt(TMath::Power((gc1[0]-gc2[0]),2)+TMath::Power((gc1[1]-gc2[1]),2));
414 Float_t aux=dRad/2.+rad1;
415 curvErr=TMath::Sqrt(curvRadius*curvRadius-dRad*dRad/4.)-TMath::Sqrt(curvRadius*curvRadius-aux*aux); //cm
418 sigmasq[0]=(cov[0]+curvErr*curvErr/2.)*factor*factor;
419 sigmasq[1]=(cov[3]+curvErr*curvErr/2.)*factor*factor;
420 sigmasq[2]=cov[5]*factor*factor;
422 // Multiple scattering
424 Float_t beta2=beta*beta;
425 Float_t p2=fMeanPSelTrk*fMeanPSelTrk;
426 Float_t rBP=GetPipeRadius();
427 Float_t dBP=0.08/35.3; // 800 um of Be
428 Float_t dL1=0.01; //approx. 1% of radiation length
429 Float_t theta2BP=14.1*14.1/(beta2*p2*1e6)*TMath::Abs(dBP);
430 Float_t theta2L1=14.1*14.1/(beta2*p2*1e6)*TMath::Abs(dL1);
431 Float_t thetaBP=TMath::Sqrt(theta2BP);
432 Float_t thetaL1=TMath::Sqrt(theta2L1);
433 for(Int_t ico=0; ico<3;ico++){
434 // printf("Error on coord. %d due to cov matrix+curvErr=%f\n",ico,sigmasq[ico]);
435 // // sigmasq[ico]+=rad1*rad1*geomfac[ico]*theta2L1/2; // multiple scattering in layer 1
436 // // sigmasq[ico]+=rBP*rBP*geomfac[ico]*theta2BP/2; // multiple scattering in beam pipe
437 sigmasq[ico]+=TMath::Power(rad1*TMath::Tan(thetaL1),2)/3.;
438 sigmasq[ico]+=TMath::Power(rBP*TMath::Tan(thetaBP),2)/3.;
440 // printf("Multipl. scatt. contr %d = %f (LAY1), %f (BP)\n",ico,rad1*rad1*geomfac[ico]*theta2L1/2,rBP*rBP*geomfac[ico]*theta2BP/2);
441 // printf("Total error on coord %d = %f\n",ico,sigmasq[ico]);
443 Float_t wmat[9]={1.,0.,0.,0.,1.,0.,0.,0.,1.};
444 if(sigmasq[0]!=0.) wmat[0]=1./sigmasq[0];
445 if(sigmasq[1]!=0.) wmat[4]=1./sigmasq[1];
446 if(sigmasq[2]!=0.) wmat[8]=1./sigmasq[2];
447 new(fLines[nolines++])AliStrLine(gc1,sigmasq,wmat,gc2,kTRUE,idClu1,idClu2);
450 fDetTypeRec->ResetRecPoints();
456 if(nolines == 0)return -2;
460 //______________________________________________________________________
461 Int_t AliITSVertexer3D::Prepare3DVertex(Int_t optCuts){
462 // Finds the 3D vertex information using tracklets
465 Float_t xbeam=GetNominalPos()[0];
466 Float_t ybeam=GetNominalPos()[1];
468 Float_t deltaR=fCoarseMaxRCut;
469 if(fPileupAlgo == 2) {
471 if(optCuts!=0)AliWarning(Form("fPileupAlgo=2 AND optCuts=%d. It should be 0",optCuts));
473 Float_t dZmax=fZCutDiamond;
475 xbeam=fVert3D.GetXv();
476 ybeam=fVert3D.GetYv();
477 zvert=fVert3D.GetZv();
480 }else if(optCuts==2){
481 xbeam=fVert3D.GetXv();
482 ybeam=fVert3D.GetYv();
486 Float_t rl=-fCoarseMaxRCut;
487 Float_t rh=fCoarseMaxRCut;
488 Float_t zl=-fZCutDiamond;
489 Float_t zh=fZCutDiamond;
490 Int_t nbr=(Int_t)((rh-rl)/fBinSizeR+0.0001);
491 Int_t nbz=(Int_t)((zh-zl)/fBinSizeZ+0.0001);
492 Int_t nbrcs=(Int_t)((rh-rl)/(fBinSizeR*2.)+0.0001);
493 Int_t nbzcs=(Int_t)((zh-zl)/(fBinSizeZ*2.)+0.0001);
498 h3d = new TH3F("h3d","xyz distribution",nbr,rl,rh,nbr,rl,rh,nbz,zl,zh);
499 h3dcs = new TH3F("h3dcs","xyz distribution",nbrcs,rl,rh,nbrcs,rl,rh,nbzcs,zl,zh);
501 // cleanup of the TCLonesArray of tracklets (i.e. fakes are removed)
502 Int_t *validate = new Int_t [fLines.GetEntriesFast()];
503 for(Int_t i=0; i<fLines.GetEntriesFast();i++)validate[i]=0;
504 for(Int_t i=0; i<fLines.GetEntriesFast()-1;i++){
505 AliStrLine *l1 = (AliStrLine*)fLines.At(i);
506 for(Int_t j=i+1;j<fLines.GetEntriesFast();j++){
507 AliStrLine *l2 = (AliStrLine*)fLines.At(j);
508 Float_t dca=l1->GetDCA(l2);
509 if(dca > fDCAcut || dca<0.00001) continue;
511 Int_t retc = l1->Cross(l2,point);
513 Double_t deltaZ=point[2]-zvert;
514 if(TMath::Abs(deltaZ)>dZmax)continue;
515 Double_t rad=TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
516 if(rad>fCoarseMaxRCut)continue;
517 Double_t deltaX=point[0]-xbeam;
518 Double_t deltaY=point[1]-ybeam;
519 Double_t raddist=TMath::Sqrt(deltaX*deltaX+deltaY*deltaY);
520 if(raddist>deltaR)continue;
523 if(fPileupAlgo != 2){
524 h3d->Fill(point[0],point[1],point[2]);
525 h3dcs->Fill(point[0],point[1],point[2]);
532 Int_t numbtracklets=0;
533 for(Int_t i=0; i<fLines.GetEntriesFast();i++)if(validate[i]>=1)numbtracklets++;
536 if(fPileupAlgo != 2){
543 for(Int_t i=0; i<fLines.GetEntriesFast();i++){
544 if(validate[i]<1)fLines.RemoveAt(i);
547 AliDebug(1,Form("Number of tracklets (after compress)%d ",fLines.GetEntriesFast()));
550 // Exit here if Pileup Algorithm 2 has been chosen
551 if(fPileupAlgo == 2)return 0;
555 // Find peaks in histos
557 Double_t peak[3]={0.,0.,0.};
559 FindPeaks(h3d,peak,ntrkl,ntimes);
561 Float_t binsizer=(rh-rl)/nbr;
562 Float_t binsizez=(zh-zl)/nbz;
563 if(optCuts==0 && (ntrkl<=2 || ntimes>1)){
566 FindPeaks(h3dcs,peak,ntrkl,ntimes);
567 binsizer=(rh-rl)/nbrcs;
568 binsizez=(zh-zl)/nbzcs;
569 if(ntrkl==1 || ntimes>1){delete h3dcs; return retcode;}
573 // Second selection loop
575 Float_t bs=(binsizer+binsizez)/2.;
576 for(Int_t i=0; i<fLines.GetEntriesFast();i++){
577 AliStrLine *l1 = (AliStrLine*)fLines.At(i);
578 if(l1->GetDistFromPoint(peak)>2.5*bs)fLines.RemoveAt(i);
581 AliDebug(1,Form("Number of tracklets (after 2nd compression) %d",fLines.GetEntriesFast()));
583 if(fLines.GetEntriesFast()>1){
585 // find a first candidate for the primary vertex
586 fVert3D=AliVertexerTracks::TrackletVertexFinder(&fLines,0);
587 // make a further selection on tracklets based on this first candidate
588 fVert3D.GetXYZ(peak);
589 AliDebug(1,Form("FIRST V candidate: %f ; %f ; %f",peak[0],peak[1],peak[2]));
590 for(Int_t i=0; i<fLines.GetEntriesFast();i++){
591 AliStrLine *l1 = (AliStrLine*)fLines.At(i);
592 if(l1->GetDistFromPoint(peak)> fDCAcut)fLines.RemoveAt(i);
595 AliDebug(1,Form("Number of tracklets (after 3rd compression) %d",fLines.GetEntriesFast()));
596 if(fLines.GetEntriesFast()>1){// this new tracklet selection is used
597 fVert3D=AliVertexerTracks::TrackletVertexFinder(&fLines,0);
603 //________________________________________________________
604 void AliITSVertexer3D::SetMeanPPtSelTracks(Float_t fieldTesla){
605 // Sets mean values of P and Pt based on the field
606 if(TMath::Abs(fieldTesla-0.5)<0.01){
607 SetMeanPSelTracks(0.885);
608 SetMeanPtSelTracks(0.630);
609 }else if(TMath::Abs(fieldTesla-0.4)<0.01){
610 SetMeanPSelTracks(0.805);
611 SetMeanPtSelTracks(0.580);
612 }else if(TMath::Abs(fieldTesla-0.2)<0.01){
613 SetMeanPSelTracks(0.740);
614 SetMeanPtSelTracks(0.530);
615 }else if(fieldTesla<0.00001){
616 SetMeanPSelTracks(0.730);
617 SetMeanPtSelTracks(0.510);
620 SetMeanPtSelTracks();
624 //________________________________________________________
625 void AliITSVertexer3D::FindPeaks(TH3F* histo, Double_t *peak, Int_t &nOfTracklets, Int_t &nOfTimes){
626 // Finds bin with max contents in 3D histo of tracket intersections
627 TAxis *xax = histo->GetXaxis();
628 TAxis *yax = histo->GetYaxis();
629 TAxis *zax = histo->GetZaxis();
635 Int_t peakbin[3]={0,0,0};
636 Int_t peak2bin[3]={-1,-1,-1};
638 for(Int_t i=xax->GetFirst();i<=xax->GetLast();i++){
639 Float_t xval = xax->GetBinCenter(i);
640 for(Int_t j=yax->GetFirst();j<=yax->GetLast();j++){
641 Float_t yval = yax->GetBinCenter(j);
642 for(Int_t k=zax->GetFirst();k<=zax->GetLast();k++){
643 Float_t zval = zax->GetBinCenter(k);
644 Int_t bc =(Int_t)histo->GetBinContent(i,j,k);
659 }else if(bc==nOfTracklets){
660 if(TMath::Abs(i-peakbin[0])<=1 && TMath::Abs(j-peakbin[1])<=1 && TMath::Abs(k-peakbin[2])<=1){
673 if(peak2bin[0]>=-1 && bc2!=-1){ // two contiguous peak-cells with same contents
674 peak[0]=0.5*(xax->GetBinCenter(peakbin[0])+xax->GetBinCenter(peak2bin[0]));
675 peak[1]=0.5*(yax->GetBinCenter(peakbin[1])+yax->GetBinCenter(peak2bin[1]));
676 peak[2]=0.5*(zax->GetBinCenter(peakbin[2])+zax->GetBinCenter(peak2bin[2]));
681 //________________________________________________________
682 void AliITSVertexer3D::MarkUsedClusters(){
683 // Mark clusters of tracklets used in vertex claulation
684 for(Int_t i=0; i<fLines.GetEntriesFast();i++){
685 AliStrLine *lin = (AliStrLine*)fLines.At(i);
686 Int_t idClu1=lin->GetIdPoint(0);
687 Int_t idClu2=lin->GetIdPoint(1);
688 fUsedCluster.SetBitNumber(idClu1);
689 fUsedCluster.SetBitNumber(idClu2);
692 //________________________________________________________
693 Int_t AliITSVertexer3D::RemoveTracklets(){
694 // Remove trackelts close to first found vertex
695 Double_t vert[3]={fVert3D.GetXv(),fVert3D.GetYv(),fVert3D.GetZv()};
697 for(Int_t i=0; i<fLines.GetEntriesFast();i++){
698 AliStrLine *lin = (AliStrLine*)fLines.At(i);
699 if(lin->GetDistFromPoint(vert)<fDCAforPileup){
700 Int_t idClu1=lin->GetIdPoint(0);
701 Int_t idClu2=lin->GetIdPoint(1);
702 fUsedCluster.SetBitNumber(idClu1);
703 fUsedCluster.SetBitNumber(idClu2);
711 //________________________________________________________
712 void AliITSVertexer3D::FindOther3DVertices(TTree *itsClusterTree){
713 // pileup identification based on 3D vertexing with not used clusters
716 Int_t nolines = FindTracklets(itsClusterTree,2);
718 Int_t nr=RemoveTracklets();
721 Int_t rc=Prepare3DVertex(2);
723 fVert3D=AliVertexerTracks::TrackletVertexFinder(&fLines,0);
724 if(fVert3D.GetNContributors()>=fMinTrackletsForPilup){
727 fVertArray = new AliESDVertex[2];
728 fVertArray[0]=(*fCurrentVertex);
729 fVertArray[1]=fVert3D;
730 fZpuv=fVert3D.GetZv();
731 fNTrpuv=fVert3D.GetNContributors();
737 //______________________________________________________________________
738 void AliITSVertexer3D::PileupFromZ(){
739 // Calls the pileup algorithm of ALiITSVertexerZ
740 Int_t binmin, binmax;
741 Int_t nPeaks=AliITSVertexerZ::GetPeakRegion(fZHisto,binmin,binmax);
742 if(nPeaks==2)AliWarning("2 peaks found");
743 Int_t firstPeakCont=0;
744 Float_t firstPeakPos=0.;
745 for(Int_t i=binmin-1;i<=binmax+1;i++){
746 firstPeakCont+=static_cast<Int_t>(fZHisto->GetBinContent(i));
747 firstPeakPos+=fZHisto->GetBinContent(i)*fZHisto->GetBinCenter(i);
750 firstPeakPos/=firstPeakCont;
752 if(firstPeakCont>fMinTrackletsForPilup){
754 ncontr2=AliITSVertexerZ::FindSecondPeak(fZHisto,binmin,binmax,secPeakPos);
755 if(ncontr2>=fMinTrackletsForPilup){
758 AliESDVertex secondVert(secPeakPos,0.1,ncontr2);
759 fVertArray = new AliESDVertex[2];
760 fVertArray[0]=(*fCurrentVertex);
761 fVertArray[1]=secondVert;
768 //________________________________________________________
769 void AliITSVertexer3D::PrintStatus() const {
770 // Print current status
771 printf("=======================================================\n");
772 printf("Loose cut on Delta Phi %f\n",fCoarseDiffPhiCut);
773 printf("Cut on tracklet DCA to Z axis %f\n",fCoarseMaxRCut);
774 printf("Cut on tracklet DCA to beam axis %f\n",fMaxRCut);
775 printf("Cut on diamond (Z) %f\n",fZCutDiamond);
776 printf("Cut on DCA - tracklet to tracklet and to vertex %f\n",fDCAcut);
777 printf("Max Phi difference: %f\n",fDiffPhiMax);
778 printf("Pileup algo: %d\n",fPileupAlgo);
779 printf("Min DCA to 1st vetrtex for pileup: %f\n",fDCAforPileup);
780 printf("=======================================================\n");