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
33 // It can be used successfully with Pb-Pb collisions
34 ////////////////////////////////////////////////////////////////
36 ClassImp(AliITSVertexer3D)
40 //______________________________________________________________________
41 AliITSVertexer3D::AliITSVertexer3D():
43 fLines("AliStrLine",1000),
45 fCoarseDiffPhiCut(0.),
57 fUsedCluster(kMaxCluPerMod*kNSPDMod),
60 fDiffPhiforPileup(0.),
65 // Default constructor
66 SetCoarseDiffPhiCut();
79 SetDeltaPhiforPileup();
83 Double_t binsize=0.02; // default 200 micron
84 Int_t nbins=static_cast<Int_t>(1+2*fZCutDiamond/binsize);
85 fZHisto=new TH1F("hz","",nbins,-fZCutDiamond,-fZCutDiamond+binsize*nbins);
88 //______________________________________________________________________
89 AliITSVertexer3D::~AliITSVertexer3D() {
92 if(fZHisto) delete fZHisto;
95 //______________________________________________________________________
96 void AliITSVertexer3D::ResetVert3D(){
102 fVert3D.SetDispersion(0.);
103 fVert3D.SetNContributors(0);
104 fUsedCluster.ResetAllBits(0);
106 //______________________________________________________________________
107 AliESDVertex* AliITSVertexer3D::FindVertexForCurrentEvent(TTree *itsClusterTree){
108 // Defines the AliESDVertex for the current event
110 AliDebug(1,"FindVertexForCurrentEvent - 3D - PROCESSING NEXT EVENT");
112 fCurrentVertex = NULL;
114 Int_t nolines = FindTracklets(itsClusterTree,0);
117 rc=Prepare3DVertex(0);
118 if(fVert3D.GetNContributors()>0){
120 nolines = FindTracklets(itsClusterTree,1);
122 rc=Prepare3DVertex(1);
123 if(fPileupAlgo == 2 && rc == 0) FindVertex3DIterative();
124 else if(fPileupAlgo<2 && rc == 0) FindVertex3D(itsClusterTree);
125 if(rc!=0) fVert3D.SetNContributors(0); // exclude this vertex
131 AliITSVertexerZ vertz(GetNominalPos()[0],GetNominalPos()[1]);
132 vertz.SetDetTypeRec(GetDetTypeRec());
133 AliDebug(1,"Call Vertexer Z\n");
134 vertz.SetLowLimit(-fZCutDiamond);
135 vertz.SetHighLimit(fZCutDiamond);
136 AliESDVertex* vtxz = vertz.FindVertexForCurrentEvent(itsClusterTree);
138 Double_t position[3]={GetNominalPos()[0],GetNominalPos()[1],vtxz->GetZv()};
139 Double_t covmatrix[6];
140 vtxz->GetCovMatrix(covmatrix);
141 Double_t chi2=99999.;
142 Int_t nContr=vtxz->GetNContributors();
143 fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);
144 fCurrentVertex->SetTitle("vertexer: Z");
145 fCurrentVertex->SetName("SPDVertexZ");
150 FindMultiplicity(itsClusterTree);
151 return fCurrentVertex;
154 //______________________________________________________________________
155 void AliITSVertexer3D::FindVertex3D(TTree *itsClusterTree){
157 Double_t vRadius=TMath::Sqrt(fVert3D.GetXv()*fVert3D.GetXv()+fVert3D.GetYv()*fVert3D.GetYv());
158 if(vRadius<GetPipeRadius() && fVert3D.GetNContributors()>0){
159 Double_t position[3]={fVert3D.GetXv(),fVert3D.GetYv(),fVert3D.GetZv()};
160 Double_t covmatrix[6];
161 fVert3D.GetCovMatrix(covmatrix);
162 Double_t chi2=99999.;
163 Int_t nContr=fVert3D.GetNContributors();
164 fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);
165 fCurrentVertex->SetTitle("vertexer: 3D");
166 fCurrentVertex->SetName("SPDVertex3D");
167 fCurrentVertex->SetDispersion(fVert3D.GetDispersion());
171 case 0: PileupFromZ(); break;
172 case 1: FindOther3DVertices(itsClusterTree); break;
173 default: AliError("Wrong pileup algorithm"); break;
176 fVertArray = new AliESDVertex[1];
177 fVertArray[0]=(*fCurrentVertex);
182 //______________________________________________________________________
183 void AliITSVertexer3D::FindVertex3DIterative(){
186 Int_t nLines=fLines.GetEntriesFast();
187 Int_t maxPoints=nLines*(nLines-1)/2;
188 Double_t* xP=new Double_t[maxPoints];
189 Double_t* yP=new Double_t[maxPoints];
190 Double_t* zP=new Double_t[maxPoints];
191 Int_t* index1=new Int_t[maxPoints];
192 Int_t* index2=new Int_t[maxPoints];
193 Double_t xbeam=fVert3D.GetXv();
194 Double_t ybeam=fVert3D.GetYv();
197 for(Int_t ilin1=0; ilin1<nLines; ilin1++){
198 AliStrLine *l1 = (AliStrLine*)fLines.At(ilin1);
199 for(Int_t ilin2=ilin1+1; ilin2<nLines; ilin2++){
200 AliStrLine *l2 = (AliStrLine*)fLines.At(ilin2);
201 Double_t dca=l1->GetDCA(l2);
202 if(dca > fDCAcut || dca<0.00001) continue;
204 Int_t retc = l1->Cross(l2,point);
206 Double_t rad=TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
207 if(rad>fCoarseMaxRCut)continue;
208 Double_t distFromBeam=TMath::Sqrt((point[0]-xbeam)*(point[0]-xbeam)+(point[1]-ybeam)*(point[1]-ybeam));
209 if(distFromBeam>fMaxRCut2) continue;
213 index1[iPoint]=ilin1;
214 index2[iPoint]=ilin2;
218 Int_t npoints=iPoint++;
220 Short_t* mask=new Short_t[npoints];
221 for(Int_t ip=0;ip<npoints;ip++) mask[ip]=-1;
223 for(Int_t ip1=0;ip1<npoints;ip1++){
224 if(mask[ip1]==-1) mask[ip1]=index++;
225 for(Int_t ip2=ip1+1; ip2<npoints; ip2++){
226 if(mask[ip2]==mask[ip1] && mask[ip2]!=-1) continue;
227 Double_t dist2=(xP[ip1]-xP[ip2])*(xP[ip1]-xP[ip2]);
228 dist2+=(yP[ip1]-yP[ip2])*(yP[ip1]-yP[ip2]);
229 dist2+=(zP[ip1]-zP[ip2])*(zP[ip1]-zP[ip2]);
230 if(dist2<fCutOnPairs*fCutOnPairs){
231 if(mask[ip2]==-1) mask[ip2]=mask[ip1];
233 for(Int_t ip=0; ip<npoints;ip++){
234 if(mask[ip]==mask[ip2]) mask[ip]=mask[ip1];
242 // Count multiplicity of trackelts in clusters
243 UInt_t* isIndUsed=new UInt_t[index+1];
244 for(Int_t ind=0;ind<index+1;ind++) isIndUsed[ind]=0;
245 for(Int_t ip=0; ip<npoints;ip++){
250 // Count clusters/vertices and sort according to multiplicity
252 Int_t* sortedIndex=new Int_t[index+1];
253 for(Int_t ind1=0;ind1<index+1;ind1++){
254 if(isIndUsed[ind1]<=1) isIndUsed[ind1]=0;
257 if(ind1>0) cap=isIndUsed[sortedIndex[ind1-1]];
259 Int_t biggerindex=-1;
260 for(Int_t ind2=0;ind2<index+1;ind2++){
262 for(Int_t ind3=0; ind3<ind1; ind3++)
263 if(ind2==sortedIndex[ind3]) use=kFALSE;
264 if(use && isIndUsed[ind2]>bigger && isIndUsed[ind2]<=cap){
265 bigger=isIndUsed[ind2];
269 sortedIndex[ind1]=biggerindex;
271 AliDebug(3,Form("Number of clusters before merging = %d\n",nClusters));
273 // Assign lines to clusters/vertices and merge clusters which share 1 line
274 Int_t nClustersAfterMerge=nClusters;
275 Int_t* belongsTo=new Int_t[nLines];
276 for(Int_t ilin=0; ilin<nLines; ilin++) belongsTo[ilin]=-1;
277 for(Int_t iclu=0;iclu<nClusters;iclu++){
278 Int_t actualCluIndex=iclu;
279 for(Int_t ip=0; ip<npoints;ip++){
280 if(mask[ip]==sortedIndex[iclu]){
281 Int_t ind1=index1[ip];
282 if(belongsTo[ind1]==-1) belongsTo[ind1]=actualCluIndex;
283 else if(belongsTo[ind1]<actualCluIndex){
284 Int_t newCluIndex=belongsTo[ind1];
285 for(Int_t ilin=0; ilin<nLines; ilin++){
286 if(belongsTo[ilin]==actualCluIndex) belongsTo[ilin]=newCluIndex;
288 AliDebug(10,Form("Merged cluster %d with %d\n",actualCluIndex,newCluIndex));
289 actualCluIndex=newCluIndex;
290 nClustersAfterMerge--;
292 Int_t ind2=index2[ip];
293 if(belongsTo[ind2]==-1) belongsTo[ind2]=actualCluIndex;
294 else if(belongsTo[ind2]<actualCluIndex){
295 Int_t newCluIndex=belongsTo[ind2];
296 for(Int_t ilin=0; ilin<nLines; ilin++){
297 if(belongsTo[ilin]==actualCluIndex) belongsTo[ilin]=newCluIndex;
299 AliDebug(10,Form("Merged cluster %d with %d\n",actualCluIndex,newCluIndex));
300 actualCluIndex=newCluIndex;
301 nClustersAfterMerge--;
306 AliDebug(3,Form("Number of clusters after merging = %d\n",nClustersAfterMerge));
308 // Count lines associated to each cluster/vertex
309 UInt_t *cluSize=new UInt_t[nClusters];
310 for(Int_t iclu=0;iclu<nClusters;iclu++){
312 for(Int_t ilin=0; ilin<nLines; ilin++){
313 if(belongsTo[ilin]==iclu) cluSize[iclu]++;
317 // Count good vertices (>1 associated tracklet)
319 for(Int_t iclu=0;iclu<nClusters;iclu++){
320 AliDebug(3,Form("Vertex %d Size=%d\n",iclu,cluSize[iclu]));
321 if(cluSize[iclu]>1) nGoodVert++;
324 AliDebug(1,Form("Number of good vertices = %d\n",nGoodVert));
325 // Calculate vertex coordinates for each cluster
327 fVertArray = new AliESDVertex[nGoodVert];
329 for(Int_t iclu=0;iclu<nClusters;iclu++){
330 Int_t size=cluSize[iclu];
332 AliStrLine **arrlin = new AliStrLine*[size];
334 for(Int_t ilin=0; ilin<nLines; ilin++){
335 if(belongsTo[ilin]==iclu){
336 arrlin[nFilled++] = dynamic_cast<AliStrLine*>(fLines[ilin]);
339 AliDebug(3,Form("Vertex %d N associated tracklets = %d out of %d\n",iVert,size,nFilled));
341 fVertArray[iVert]=AliVertexerTracks::TrackletVertexFinder(arrlin,nFilled);
343 fVertArray[iVert].GetXYZ(peak);
344 AliStrLine **arrlin2 = new AliStrLine*[size];
346 for(Int_t i=0; i<nFilled;i++){
347 AliStrLine *l1 = arrlin[i];
348 if(l1->GetDistFromPoint(peak)< fDCAcut)
349 arrlin2[nFilled2++] = dynamic_cast<AliStrLine*>(l1);
352 AliDebug(3,Form("Vertex %d recalculated with %d tracklets\n",iVert,nFilled2));
353 fVertArray[iVert]=AliVertexerTracks::TrackletVertexFinder(arrlin2,nFilled2);
363 fNTrpuv = fVertArray[1].GetNContributors();
364 fZpuv = fVertArray[1].GetZv();
367 Double_t vRadius=TMath::Sqrt(fVertArray[0].GetXv()*fVertArray[0].GetXv()+fVertArray[0].GetYv()*fVertArray[0].GetYv());
368 if(vRadius<GetPipeRadius() && fVertArray[0].GetNContributors()>0){
369 Double_t position[3]={fVertArray[0].GetXv(),fVertArray[0].GetYv(),fVertArray[0].GetZv()};
370 Double_t covmatrix[6];
371 fVertArray[0].GetCovMatrix(covmatrix);
372 Double_t chi2=99999.;
373 Int_t nContr=fVertArray[0].GetNContributors();
374 fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);
375 fCurrentVertex->SetTitle("vertexer: 3D");
376 fCurrentVertex->SetName("SPDVertex3D");
377 fCurrentVertex->SetDispersion(fVertArray[0].GetDispersion());
385 delete [] sortedIndex;
392 //______________________________________________________________________
393 void AliITSVertexer3D::FindVertex3DIterativeMM(){
394 // Defines the AliESDVertex for the current event
395 Int_t numsor=fLines.GetEntriesFast()*(fLines.GetEntriesFast()-1)/2;
396 //cout<<"AliITSVertexer3D::FindVertexForCurentEvent: Number of tracklets selected for vertexing "<<fLines.GetEntriesFast()<<"; Number of pairs: "<<numsor<<endl;
397 AliITSSortTrkl srt(fLines,numsor,fCutOnPairs,fCoarseMaxRCut);
399 AliInfo(Form("Number of vertices: %d",srt.GetNumberOfClusters()));
401 fNoVertices = srt.GetNumberOfClusters();
402 //printf("fNoVertices = %d \n",fNoVertices);
404 fVertArray = new AliESDVertex[fNoVertices];
405 for(Int_t kk=0; kk<srt.GetNumberOfClusters(); kk++){
407 Int_t *labels = srt.GetTrackletsLab(kk,size);
409 Int_t *pairs = srt.GetClusters(kk);
410 Int_t nopai = srt.GetSizeOfCluster(kk);
411 cout<<"***** Vertex number "<<kk<<". Pairs: \n";
412 for(Int_t jj=0;jj<nopai;jj++){
413 cout<<pairs[jj]<<" - ";
414 if(jj>0 & jj%8==0)cout<<endl;
417 cout<<"***** Vertex number "<<kk<<". Labels: \n";
419 AliStrLine **tclo = new AliStrLine* [size];
420 for(Int_t jj=0;jj<size;jj++){
421 // cout<<labels[jj]<<" - ";
422 // if(jj>0 & jj%8==0)cout<<endl;
423 tclo[jj] = dynamic_cast<AliStrLine*>(fLines[labels[jj]]);
427 fVertArray[kk]=AliVertexerTracks::TrackletVertexFinder(tclo,size);
429 // fVertArray[kk].PrintStatus();
431 // at least one second vertex is present
433 fNTrpuv = fVertArray[kk].GetNContributors();
434 fZpuv = fVertArray[kk].GetZv();
437 Double_t vRadius=TMath::Sqrt(fVertArray[0].GetXv()*fVertArray[0].GetXv()+fVertArray[0].GetYv()*fVertArray[0].GetYv());
438 if(vRadius<GetPipeRadius() && fVertArray[0].GetNContributors()>0){
439 Double_t position[3]={fVertArray[0].GetXv(),fVertArray[0].GetYv(),fVertArray[0].GetZv()};
440 Double_t covmatrix[6];
441 fVertArray[0].GetCovMatrix(covmatrix);
442 Double_t chi2=99999.;
443 Int_t nContr=fVertArray[0].GetNContributors();
444 fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nContr);
445 fCurrentVertex->SetTitle("vertexer: 3D");
446 fCurrentVertex->SetName("SPDVertex3D");
447 fCurrentVertex->SetDispersion(fVertArray[0].GetDispersion());
453 //______________________________________________________________________
454 Bool_t AliITSVertexer3D::DistBetweenVertices(AliESDVertex &a, AliESDVertex &b, Double_t test, Double_t &dist){
455 // method to compare the distance between vertices a and b with "test"
456 //it returns kTRUE is the distance is less or equal to test
457 dist = (a.GetX()-b.GetX()) * (a.GetX()-b.GetX());
458 dist += (a.GetY()-b.GetY()) * (a.GetY()-b.GetY());
459 dist += (a.GetZ()-b.GetZ()) * (a.GetZ()-b.GetZ());
460 dist = TMath::Sqrt(dist);
461 if(dist <= test)return kTRUE;
466 //______________________________________________________________________
467 Int_t AliITSVertexer3D::FindTracklets(TTree *itsClusterTree, Int_t optCuts){
468 // All the possible combinations between recpoints on layer 1and 2 are
469 // considered. Straight lines (=tracklets)are formed.
470 // The tracklets are processed in Prepare3DVertex
472 TClonesArray *itsRec = 0;
473 if(optCuts==0) fZHisto->Reset();
474 // gc1 are local and global coordinates for layer 1
475 Float_t gc1f[3]={0.,0.,0.};
476 Double_t gc1[3]={0.,0.,0.};
477 // gc2 are local and global coordinates for layer 2
478 Float_t gc2f[3]={0.,0.,0.};
479 Double_t gc2[3]={0.,0.,0.};
480 AliITSRecPointContainer* rpcont=AliITSRecPointContainer::Instance();
481 itsRec=rpcont->FetchClusters(0,itsClusterTree);
482 if(!rpcont->IsSPDActive()){
483 AliWarning("No SPD rec points found, 3D vertex not calculated");
487 // Set values for cuts
488 Double_t xbeam=GetNominalPos()[0];
489 Double_t ybeam=GetNominalPos()[1];
491 Double_t deltaPhi=fCoarseDiffPhiCut;
492 Double_t deltaR=fCoarseMaxRCut;
493 Double_t dZmax=fZCutDiamond;
495 xbeam=fVert3D.GetXv();
496 ybeam=fVert3D.GetYv();
497 zvert=fVert3D.GetZv();
498 deltaPhi = fDiffPhiMax;
501 if(fPileupAlgo == 2){
505 } else if(optCuts==2){
506 xbeam=fVert3D.GetXv();
507 ybeam=fVert3D.GetYv();
508 deltaPhi = fDiffPhiforPileup;
512 Int_t nrpL1 = 0; // number of rec points on layer 1
513 Int_t nrpL2 = 0; // number of rec points on layer 2
515 // By default irstL1=0 and lastL1=79
516 Int_t firstL1 = AliITSgeomTGeo::GetModuleIndex(1,1,1);
517 Int_t lastL1 = AliITSgeomTGeo::GetModuleIndex(2,1,1)-1;
518 for(Int_t module= firstL1; module<=lastL1;module++){ // count number of recopints on layer 1
519 itsRec=rpcont->UncheckedGetClusters(module);
520 nrpL1+= itsRec->GetEntries();
522 //By default firstL2=80 and lastL2=239
523 Int_t firstL2 = AliITSgeomTGeo::GetModuleIndex(2,1,1);
524 Int_t lastL2 = AliITSgeomTGeo::GetModuleIndex(3,1,1)-1;
525 for(Int_t module= firstL2; module<=lastL2;module++){ // count number of recopints on layer 2
526 itsRec=rpcont->UncheckedGetClusters(module);
527 nrpL2+= itsRec->GetEntries();
529 if(nrpL1 == 0 || nrpL2 == 0){
530 AliDebug(1,Form("No RecPoints in at least one SPD layer (%d %d)",nrpL1,nrpL2));
533 AliDebug(1,Form("RecPoints on Layer 1,2 = %d, %d\n",nrpL1,nrpL2));
535 Double_t a[3]={xbeam,ybeam,0.};
536 Double_t b[3]={xbeam,ybeam,10.};
537 AliStrLine zeta(a,b,kTRUE);
538 static Double_t bField=TMath::Abs(AliTracker::GetBz()/10.); //T
539 SetMeanPPtSelTracks(bField);
542 // Loop on modules of layer 1
543 for(Int_t modul1= firstL1; modul1<=lastL1;modul1++){ // Loop on modules of layer 1
544 if(!fUseModule[modul1]) continue;
545 UShort_t ladder=int(modul1/4)+1; // ladders are numbered starting from 1
546 TClonesArray *prpl1=rpcont->UncheckedGetClusters(modul1);
547 Int_t nrecp1 = prpl1->GetEntries();
548 for(Int_t j=0;j<nrecp1;j++){
549 if(j>kMaxCluPerMod) continue;
550 UShort_t idClu1=modul1*kMaxCluPerMod+j;
551 if(fUsedCluster.TestBitNumber(idClu1)) continue;
552 AliITSRecPoint *recp1 = (AliITSRecPoint*)prpl1->At(j);
553 recp1->GetGlobalXYZ(gc1f);
554 for(Int_t ico=0;ico<3;ico++)gc1[ico]=gc1f[ico];
556 Double_t phi1 = TMath::ATan2(gc1[1]-ybeam,gc1[0]-xbeam);
557 if(phi1<0)phi1=2*TMath::Pi()+phi1;
558 for(Int_t ladl2=0 ; ladl2<fLadOnLay2*2+1;ladl2++){
559 for(Int_t k=0;k<4;k++){
560 Int_t ladmod=fLadders[ladder-1]+ladl2;
561 if(ladmod>AliITSgeomTGeo::GetNLadders(2)) ladmod=ladmod-AliITSgeomTGeo::GetNLadders(2);
562 Int_t modul2=AliITSgeomTGeo::GetModuleIndex(2,ladmod,k+1);
563 if(!fUseModule[modul2]) continue;
564 itsRec=rpcont->UncheckedGetClusters(modul2);
565 Int_t nrecp2 = itsRec->GetEntries();
566 for(Int_t j2=0;j2<nrecp2;j2++){
567 if(j2>kMaxCluPerMod) continue;
568 UShort_t idClu2=modul2*kMaxCluPerMod+j2;
569 if(fUsedCluster.TestBitNumber(idClu2)) continue;
570 AliITSRecPoint *recp2 = (AliITSRecPoint*)itsRec->At(j2);
571 recp2->GetGlobalXYZ(gc2f);
572 for(Int_t ico=0;ico<3;ico++)gc2[ico]=gc2f[ico];
573 Double_t phi2 = TMath::ATan2(gc2[1]-ybeam,gc2[0]-xbeam);
574 if(phi2<0)phi2=2*TMath::Pi()+phi2;
575 Double_t diff = TMath::Abs(phi2-phi1);
576 if(diff>TMath::Pi())diff=2.*TMath::Pi()-diff;
577 if(optCuts==0 && diff<fDiffPhiforPileup){
578 Double_t r1=TMath::Sqrt(gc1[0]*gc1[0]+gc1[1]*gc1[1]);
580 Double_t r2=TMath::Sqrt(gc2[0]*gc2[0]+gc2[1]*gc2[1]);
582 Double_t zr0=(r2*zc1-r1*zc2)/(r2-r1); //Z @ null radius
585 if(diff>deltaPhi)continue;
586 AliStrLine line(gc1,gc2,kTRUE);
588 Int_t retcode = line.Cross(&zeta,cp);
589 if(retcode<0)continue;
590 Double_t dca = line.GetDCA(&zeta);
592 if(dca>deltaR)continue;
593 Double_t deltaZ=cp[2]-zvert;
594 if(TMath::Abs(deltaZ)>dZmax)continue;
597 if(fLines.GetEntriesFast()>0)fLines.Clear("C");
600 recp2->GetGlobalCov(cov);
603 Double_t rad1=TMath::Sqrt(gc1[0]*gc1[0]+gc1[1]*gc1[1]);
604 Double_t rad2=TMath::Sqrt(gc2[0]*gc2[0]+gc2[1]*gc2[1]);
605 Double_t factor=(rad1+rad2)/(rad2-rad1); //factor to account for error on tracklet direction
609 Double_t curvRadius=fMeanPtSelTrk/(0.3*bField)*100; //cm
610 Double_t dRad=TMath::Sqrt((gc1[0]-gc2[0])*(gc1[0]-gc2[0])+(gc1[1]-gc2[1])*(gc1[1]-gc2[1]));
611 Double_t aux=dRad/2.+rad1;
612 curvErr=TMath::Sqrt(curvRadius*curvRadius-dRad*dRad/4.)-TMath::Sqrt(curvRadius*curvRadius-aux*aux); //cm
615 sigmasq[0]=(cov[0]+curvErr*curvErr/2.)*factor*factor;
616 sigmasq[1]=(cov[3]+curvErr*curvErr/2.)*factor*factor;
617 sigmasq[2]=cov[5]*factor*factor;
619 // Multiple scattering
620 Double_t pOverMass=fMeanPSelTrk/0.140;
621 Double_t beta2=pOverMass*pOverMass/(1+pOverMass*pOverMass);
622 Double_t p2=fMeanPSelTrk*fMeanPSelTrk;
623 Double_t rBP=GetPipeRadius();
624 Double_t dBP=0.08/35.3; // 800 um of Be
625 Double_t dL1=0.01; //approx. 1% of radiation length
626 Double_t theta2BP=14.1*14.1/(beta2*p2*1e6)*dBP;
627 Double_t theta2L1=14.1*14.1/(beta2*p2*1e6)*dL1;
628 Double_t rtantheta1=(rad2-rad1)*TMath::Tan(TMath::Sqrt(theta2L1));
629 Double_t rtanthetaBP=(rad1-rBP)*TMath::Tan(TMath::Sqrt(theta2BP));
630 for(Int_t ico=0; ico<3;ico++){
631 sigmasq[ico]+=rtantheta1*rtantheta1*factor*factor/3.;
632 sigmasq[ico]+=rtanthetaBP*rtanthetaBP*factor*factor/3.;
634 Double_t wmat[9]={1.,0.,0.,0.,1.,0.,0.,0.,1.};
635 if(sigmasq[0]!=0.) wmat[0]=1./sigmasq[0];
636 if(sigmasq[1]!=0.) wmat[4]=1./sigmasq[1];
637 if(sigmasq[2]!=0.) wmat[8]=1./sigmasq[2];
638 new(fLines[nolines++])AliStrLine(gc1,sigmasq,wmat,gc2,kTRUE,idClu1,idClu2);
645 if(nolines == 0)return -2;
649 //______________________________________________________________________
650 Int_t AliITSVertexer3D::Prepare3DVertex(Int_t optCuts){
651 // Finds the 3D vertex information using tracklets
654 Double_t xbeam=GetNominalPos()[0];
655 Double_t ybeam=GetNominalPos()[1];
657 Double_t deltaR=fCoarseMaxRCut;
658 Double_t dZmax=fZCutDiamond;
660 xbeam=fVert3D.GetXv();
661 ybeam=fVert3D.GetYv();
662 zvert=fVert3D.GetZv();
665 if(fPileupAlgo == 2){
669 }else if(optCuts==2){
670 xbeam=fVert3D.GetXv();
671 ybeam=fVert3D.GetYv();
675 Double_t rl=-fCoarseMaxRCut;
676 Double_t rh=fCoarseMaxRCut;
677 Double_t zl=-fZCutDiamond;
678 Double_t zh=fZCutDiamond;
679 Int_t nbr=(Int_t)((rh-rl)/fBinSizeR+0.0001);
680 Int_t nbz=(Int_t)((zh-zl)/fBinSizeZ+0.0001);
681 Int_t nbrcs=(Int_t)((rh-rl)/(fBinSizeR*2.)+0.0001);
682 Int_t nbzcs=(Int_t)((zh-zl)/(fBinSizeZ*2.)+0.0001);
684 TH3F *h3d = new TH3F("h3d","xyz distribution",nbr,rl,rh,nbr,rl,rh,nbz,zl,zh);
685 TH3F *h3dcs = new TH3F("h3dcs","xyz distribution",nbrcs,rl,rh,nbrcs,rl,rh,nbzcs,zl,zh);
687 // cleanup of the TCLonesArray of tracklets (i.e. fakes are removed)
688 Int_t *validate = new Int_t [fLines.GetEntriesFast()];
689 for(Int_t i=0; i<fLines.GetEntriesFast();i++)validate[i]=0;
690 for(Int_t i=0; i<fLines.GetEntriesFast()-1;i++){
691 AliStrLine *l1 = (AliStrLine*)fLines.At(i);
692 for(Int_t j=i+1;j<fLines.GetEntriesFast();j++){
693 AliStrLine *l2 = (AliStrLine*)fLines.At(j);
694 Double_t dca=l1->GetDCA(l2);
695 if(dca > fDCAcut || dca<0.00001) continue;
697 Int_t retc = l1->Cross(l2,point);
699 Double_t deltaZ=point[2]-zvert;
700 if(TMath::Abs(deltaZ)>dZmax)continue;
701 Double_t rad=TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
702 if(rad>fCoarseMaxRCut)continue;
703 Double_t deltaX=point[0]-xbeam;
704 Double_t deltaY=point[1]-ybeam;
705 Double_t raddist=TMath::Sqrt(deltaX*deltaX+deltaY*deltaY);
706 if(raddist>deltaR)continue;
709 h3d->Fill(point[0],point[1],point[2]);
710 h3dcs->Fill(point[0],point[1],point[2]);
716 Int_t numbtracklets=0;
717 for(Int_t i=0; i<fLines.GetEntriesFast();i++)if(validate[i]>=1)numbtracklets++;
725 for(Int_t i=0; i<fLines.GetEntriesFast();i++){
726 if(validate[i]<1)fLines.RemoveAt(i);
729 AliDebug(1,Form("Number of tracklets (after compress)%d ",fLines.GetEntriesFast()));
732 // Exit here if Pileup Algorithm 2 has been chosen during second loop
733 if(fPileupAlgo == 2 && optCuts==1){
739 // Find peaks in histos
741 Double_t peak[3]={0.,0.,0.};
743 FindPeaks(h3d,peak,ntrkl,ntimes);
745 Double_t binsizer=(rh-rl)/nbr;
746 Double_t binsizez=(zh-zl)/nbz;
747 if(optCuts==0 && (ntrkl<=2 || ntimes>1)){
750 FindPeaks(h3dcs,peak,ntrkl,ntimes);
751 binsizer=(rh-rl)/nbrcs;
752 binsizez=(zh-zl)/nbzcs;
753 if(ntrkl==1 || ntimes>1){delete h3dcs; return retcode;}
757 // Second selection loop
759 Double_t bs=(binsizer+binsizez)/2.;
760 for(Int_t i=0; i<fLines.GetEntriesFast();i++){
761 AliStrLine *l1 = (AliStrLine*)fLines.At(i);
762 if(l1->GetDistFromPoint(peak)>2.5*bs)fLines.RemoveAt(i);
765 AliDebug(1,Form("Number of tracklets (after 2nd compression) %d",fLines.GetEntriesFast()));
767 if(fLines.GetEntriesFast()>1){
769 // find a first candidate for the primary vertex
770 fVert3D=AliVertexerTracks::TrackletVertexFinder(&fLines,0);
771 // make a further selection on tracklets based on this first candidate
772 fVert3D.GetXYZ(peak);
773 AliDebug(1,Form("FIRST V candidate: %f ; %f ; %f",peak[0],peak[1],peak[2]));
774 for(Int_t i=0; i<fLines.GetEntriesFast();i++){
775 AliStrLine *l1 = (AliStrLine*)fLines.At(i);
776 if(l1->GetDistFromPoint(peak)> fDCAcut)fLines.RemoveAt(i);
779 AliDebug(1,Form("Number of tracklets (after 3rd compression) %d",fLines.GetEntriesFast()));
780 if(fLines.GetEntriesFast()>1){// this new tracklet selection is used
781 fVert3D=AliVertexerTracks::TrackletVertexFinder(&fLines,0);
787 //________________________________________________________
788 void AliITSVertexer3D::SetMeanPPtSelTracks(Double_t fieldTesla){
789 // Sets mean values of Pt based on the field
790 // for P (used in multiple scattering) the most probable value is used
791 if(TMath::Abs(fieldTesla-0.5)<0.01){
792 SetMeanPSelTracks(0.375);
793 SetMeanPtSelTracks(0.630);
794 }else if(TMath::Abs(fieldTesla-0.4)<0.01){
795 SetMeanPSelTracks(0.375);
796 SetMeanPtSelTracks(0.580);
797 }else if(TMath::Abs(fieldTesla-0.2)<0.01){
798 SetMeanPSelTracks(0.375);
799 SetMeanPtSelTracks(0.530);
800 }else if(fieldTesla<0.00001){
801 SetMeanPSelTracks(0.375);
802 SetMeanPtSelTracks(0.230);
805 SetMeanPtSelTracks();
809 //________________________________________________________
810 void AliITSVertexer3D::FindPeaks(TH3F* histo, Double_t *peak, Int_t &nOfTracklets, Int_t &nOfTimes){
811 // Finds bin with max contents in 3D histo of tracket intersections
812 TAxis *xax = histo->GetXaxis();
813 TAxis *yax = histo->GetYaxis();
814 TAxis *zax = histo->GetZaxis();
820 Int_t peakbin[3]={0,0,0};
821 Int_t peak2bin[3]={-1,-1,-1};
823 for(Int_t i=xax->GetFirst();i<=xax->GetLast();i++){
824 Double_t xval = xax->GetBinCenter(i);
825 for(Int_t j=yax->GetFirst();j<=yax->GetLast();j++){
826 Double_t yval = yax->GetBinCenter(j);
827 for(Int_t k=zax->GetFirst();k<=zax->GetLast();k++){
828 Double_t zval = zax->GetBinCenter(k);
829 Int_t bc =(Int_t)histo->GetBinContent(i,j,k);
844 }else if(bc==nOfTracklets){
845 if(TMath::Abs(i-peakbin[0])<=1 && TMath::Abs(j-peakbin[1])<=1 && TMath::Abs(k-peakbin[2])<=1){
858 if(peak2bin[0]>=-1 && bc2!=-1){ // two contiguous peak-cells with same contents
859 peak[0]=0.5*(xax->GetBinCenter(peakbin[0])+xax->GetBinCenter(peak2bin[0]));
860 peak[1]=0.5*(yax->GetBinCenter(peakbin[1])+yax->GetBinCenter(peak2bin[1]));
861 peak[2]=0.5*(zax->GetBinCenter(peakbin[2])+zax->GetBinCenter(peak2bin[2]));
866 //________________________________________________________
867 void AliITSVertexer3D::MarkUsedClusters(){
868 // Mark clusters of tracklets used in vertex claulation
869 for(Int_t i=0; i<fLines.GetEntriesFast();i++){
870 AliStrLine *lin = (AliStrLine*)fLines.At(i);
871 Int_t idClu1=lin->GetIdPoint(0);
872 Int_t idClu2=lin->GetIdPoint(1);
873 fUsedCluster.SetBitNumber(idClu1);
874 fUsedCluster.SetBitNumber(idClu2);
877 //________________________________________________________
878 Int_t AliITSVertexer3D::RemoveTracklets(){
879 // Remove trackelts close to first found vertex
880 Double_t vert[3]={fVert3D.GetXv(),fVert3D.GetYv(),fVert3D.GetZv()};
882 for(Int_t i=0; i<fLines.GetEntriesFast();i++){
883 AliStrLine *lin = (AliStrLine*)fLines.At(i);
884 if(lin->GetDistFromPoint(vert)<fDCAforPileup){
885 Int_t idClu1=lin->GetIdPoint(0);
886 Int_t idClu2=lin->GetIdPoint(1);
887 fUsedCluster.SetBitNumber(idClu1);
888 fUsedCluster.SetBitNumber(idClu2);
896 //________________________________________________________
897 void AliITSVertexer3D::FindOther3DVertices(TTree *itsClusterTree){
898 // pileup identification based on 3D vertexing with not used clusters
901 Int_t nolines = FindTracklets(itsClusterTree,2);
903 Int_t nr=RemoveTracklets();
906 Int_t rc=Prepare3DVertex(2);
908 fVert3D=AliVertexerTracks::TrackletVertexFinder(&fLines,0);
909 if(fVert3D.GetNContributors()>=fMinTrackletsForPilup){
912 fVertArray = new AliESDVertex[2];
913 fVertArray[0]=(*fCurrentVertex);
914 fVertArray[1]=fVert3D;
915 fZpuv=fVert3D.GetZv();
916 fNTrpuv=fVert3D.GetNContributors();
922 //______________________________________________________________________
923 void AliITSVertexer3D::PileupFromZ(){
924 // Calls the pileup algorithm of ALiITSVertexerZ
925 Int_t binmin, binmax;
926 Int_t nPeaks=AliITSVertexerZ::GetPeakRegion(fZHisto,binmin,binmax);
927 if(nPeaks==2)AliWarning("2 peaks found");
928 Int_t firstPeakCont=0;
929 Double_t firstPeakPos=0.;
930 for(Int_t i=binmin-1;i<=binmax+1;i++){
931 firstPeakCont+=static_cast<Int_t>(fZHisto->GetBinContent(i));
932 firstPeakPos+=fZHisto->GetBinContent(i)*fZHisto->GetBinCenter(i);
935 firstPeakPos/=firstPeakCont;
937 if(firstPeakCont>fMinTrackletsForPilup){
939 ncontr2=AliITSVertexerZ::FindSecondPeak(fZHisto,binmin,binmax,secPeakPos);
940 if(ncontr2>=fMinTrackletsForPilup){
943 AliESDVertex secondVert(secPeakPos,0.1,ncontr2);
944 fVertArray = new AliESDVertex[2];
945 fVertArray[0]=(*fCurrentVertex);
946 fVertArray[1]=secondVert;
953 //________________________________________________________
954 void AliITSVertexer3D::PrintStatus() const {
955 // Print current status
956 printf("========= First step selections =====================\n");
957 printf("Cut on diamond (Z) %f\n",fZCutDiamond);
958 printf("Loose cut on Delta Phi %f\n",fCoarseDiffPhiCut);
959 printf("Loose cut on tracklet DCA to Z axis %f\n",fCoarseMaxRCut);
960 printf("Cut on DCA - tracklet to tracklet and to vertex %f\n",fDCAcut);
961 printf("========= Second step selections ====================\n");
962 printf("Cut on tracklet-to-first-vertex Z distance %f\n",fMaxZCut);
963 printf("Max Phi difference: %f\n",fDiffPhiMax);
964 printf("Cut on tracklet DCA to beam axis %f\n",fMaxRCut);
965 printf("Cut on tracklet DCA to beam axis (algo2) %f\n",fMaxRCut2);
966 printf("========= Pileup selections =========================\n");
967 printf("Pileup algo: %d\n",fPileupAlgo);
968 printf("Min DCA to 1st vertex for pileup (algo 0 and 1): %f\n",fDCAforPileup);
969 printf("Cut on distance between pair-vertices (algo 2): %f\n",fCutOnPairs);
970 printf("=======================================================\n");