1 /**************************************************************************
2 * Copyright(c) 1998-2003, 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 //-----------------------------------------------------------------
17 // Implementation of the vertexer from ESD tracks
19 // Origin: AliITSVertexerTracks
21 // andrea.dainese@pd.infn.it
23 // massimo.masera@to.infn.it
24 // Moved to STEER and adapted to ESD tracks:
25 // F.Prino, Torino, prino@to.infn.it
26 //-----------------------------------------------------------------
28 //---- Root headers --------
31 //---- AliRoot headers -----
32 #include "AliStrLine.h"
33 #include "AliVertexerTracks.h"
35 #include "AliESDtrack.h"
37 ClassImp(AliVertexerTracks)
40 //----------------------------------------------------------------------------
41 AliVertexerTracks::AliVertexerTracks():
61 // Default constructor
70 //-----------------------------------------------------------------------------
71 AliVertexerTracks::AliVertexerTracks(Double_t xStart, Double_t yStart):
91 // Standard constructor
93 SetVtxStart(xStart,yStart);
100 //-----------------------------------------------------------------------------
101 AliVertexerTracks::~AliVertexerTracks()
103 // Default Destructor
104 // The objects pointed by the following pointer are not owned
105 // by this class and are not deleted
107 delete[] fTrksToSkip;
109 //----------------------------------------------------------------------------
110 AliESDVertex* AliVertexerTracks::FindPrimaryVertex(const AliESD *esdEvent)
113 // Primary vertex for current ESD event
115 // 1st with 5*fNSigma*sigma cut w.r.t. to initial vertex
116 // + cut on sqrt(d0d0+z0z0) if fConstraint=kFALSE
117 // 2nd with fNSigma*sigma cut w.r.t. to vertex found in 1st iteration)
118 // All ESD tracks with inside the beam pipe are then propagated to found vertex
122 // read tracks from ESD
123 Int_t nTrksTot = (Int_t)esdEvent->GetNumberOfTracks();
124 if (nTrksTot<=0) return fCurrentVertex;
126 TTree *trkTree = new TTree("TreeT","tracks");
127 AliESDtrack *esdTrack = 0;
128 trkTree->Branch("tracks","AliESDtrack",&esdTrack);
131 for(Int_t i=0; i<nTrksTot; i++) {
132 AliESDtrack *et = esdEvent->GetTrack(i);
133 esdTrack = new AliESDtrack(*et);
134 // check tracks to skip
136 for(Int_t j=0; j<fNTrksToSkip; j++) {
137 if(et->GetID()==fTrksToSkip[j]) {
138 if(fDebug) printf("skipping track: %d\n",i);
142 if(skipThis) {delete esdTrack;continue;}
144 if(!(esdTrack->GetStatus()&AliESDtrack::kITSin)) {delete esdTrack;continue;}
145 if(fITSrefit && !(esdTrack->GetStatus()&AliESDtrack::kITSrefit)) {delete esdTrack;continue;}
146 Int_t nclus=esdTrack->GetNcls(0); // check number of clusters in ITS
147 if(nclus<fMinITSClusters) {delete esdTrack;continue;}
153 // If fConstraint=kFALSE
154 // run VertexFinder(1) to get rough estimate of initVertex (x,y)
156 // fill fTrkArray, for VertexFinder()
157 if(!fTrkArray.IsEmpty()) fTrkArray.Delete();
158 PrepareTracks(*trkTree,0);
159 Double_t cutsave = fDCAcut; fDCAcut = 0.1; // 1 mm
160 VertexFinder(1); // using weights, cutting dca < fDCAcut
163 if(fVert.GetNContributors()>0) {
164 fVert.GetXYZ(fNominalPos);
165 fNominalPos[0] = fVert.GetXv();
166 fNominalPos[1] = fVert.GetYv();
167 fNominalPos[2] = fVert.GetZv();
168 if(fDebug) printf("No mean vertex: VertexFinder gives (%f, %f, %f)\n",fNominalPos[0],fNominalPos[1],fNominalPos[2]);
173 if(fDebug) printf("No mean vertex and VertexFinder failed\n");
180 // propagate tracks to fNominalPos vertex
182 // if(constraint) reject for |d0|>5*fNSigma*sigma w.r.t. fNominal... vertex
183 // else reject for |d0|\oplus|z0| > 5 mm w.r.t. fNominal... vertex
185 // propagate tracks to best between initVertex and fCurrentVertex
186 // preselect tracks (reject for |d0|>fNSigma*sigma w.r.t. best
187 // between initVertex and fCurrentVertex)
188 for(Int_t iter=0; iter<2; iter++) {
189 if(fOnlyFitter && iter==0) continue;
190 Int_t nTrksPrep = PrepareTracks(*trkTree,iter+1);
191 if(fDebug) printf(" tracks prepared - iteration %d: %d\n",iter+1,nTrksPrep);
192 if(nTrksPrep < fMinTracks) {
193 if(fDebug) printf("TooFewTracks\n");
194 TooFewTracks(esdEvent);
195 if(fDebug) fCurrentVertex->PrintStatus();
198 return fCurrentVertex;
204 if(fDebug) printf("Just one track\n");
205 OneTrackVertFinder();
208 case 1: StrLinVertexFinderMinDist(1); break;
209 case 2: StrLinVertexFinderMinDist(0); break;
210 case 3: HelixVertexFinder(); break;
211 case 4: VertexFinder(1); break;
212 case 5: VertexFinder(0); break;
213 default: printf("Wrong algorithm\n"); break;
216 if(fDebug) printf(" Vertex finding completed\n");
220 VertexFitter(fConstraint);
221 if(fDebug) printf(" Vertex fit completed\n");
222 if(iter==0) fTrkArray.Delete();
223 } // end loop on the two iterations
226 // take true pos from SPD vertex in ESD and write it in tracks' vertex
228 esdEvent->GetVertex()->GetTruePos(tp);
229 fCurrentVertex->SetTruePos(tp);
233 fCurrentVertex->SetTitle("VertexerTracksWithConstraintOnlyFitter");
235 fCurrentVertex->SetTitle("VertexerTracksWithConstraint");
238 fCurrentVertex->SetTitle("VertexerTracksNoConstraint");
242 // set indices of used tracks
243 UShort_t *indices = 0;
244 AliESDtrack *ett = 0;
245 if(fCurrentVertex->GetNContributors()>0) {
246 indices = new UShort_t[fCurrentVertex->GetNContributors()];
247 for(Int_t jj=0;jj<(Int_t)fTrkArray.GetEntriesFast();jj++) {
248 ett = (AliESDtrack*)fTrkArray.At(jj);
249 indices[jj] = (UShort_t)ett->GetID();
251 fCurrentVertex->SetIndices(fCurrentVertex->GetNContributors(),indices);
259 if(fTrksToSkip) delete [] fTrksToSkip;
262 if(fDebug) fCurrentVertex->PrintStatus();
263 if(fDebug) fCurrentVertex->PrintIndices();
266 return fCurrentVertex;
268 //------------------------------------------------------------------------
269 Double_t AliVertexerTracks::GetDeterminant3X3(Double_t matr[][3])
272 Double_t det=matr[0][0]*matr[1][1]*matr[2][2]-matr[0][0]*matr[1][2]*matr[2][1]-matr[0][1]*matr[1][0]*matr[2][2]+matr[0][1]*matr[1][2]*matr[2][0]+matr[0][2]*matr[1][0]*matr[2][1]-matr[0][2]*matr[1][1]*matr[2][0];
275 //-------------------------------------------------------------------------
276 void AliVertexerTracks::GetStrLinDerivMatrix(Double_t *p0,Double_t *p1,Double_t (*m)[3],Double_t *d)
279 Double_t x12=p0[0]-p1[0];
280 Double_t y12=p0[1]-p1[1];
281 Double_t z12=p0[2]-p1[2];
282 Double_t kk=x12*x12+y12*y12+z12*z12;
283 m[0][0]=2-2/kk*x12*x12;
284 m[0][1]=-2/kk*x12*y12;
285 m[0][2]=-2/kk*x12*z12;
286 m[1][0]=-2/kk*x12*y12;
287 m[1][1]=2-2/kk*y12*y12;
288 m[1][2]=-2/kk*y12*z12;
289 m[2][0]=-2/kk*x12*z12;
291 m[2][2]=2-2/kk*z12*z12;
292 d[0]=2*p0[0]-2/kk*p0[0]*x12*x12-2/kk*p0[2]*x12*z12-2/kk*p0[1]*x12*y12;
293 d[1]=2*p0[1]-2/kk*p0[1]*y12*y12-2/kk*p0[0]*x12*y12-2/kk*p0[2]*z12*y12;
294 d[2]=2*p0[2]-2/kk*p0[2]*z12*z12-2/kk*p0[0]*x12*z12-2/kk*p0[1]*z12*y12;
297 //--------------------------------------------------------------------------
298 void AliVertexerTracks::GetStrLinDerivMatrix(Double_t *p0,Double_t *p1,Double_t *sigmasq,Double_t (*m)[3],Double_t *d)
301 Double_t x12=p1[0]-p0[0];
302 Double_t y12=p1[1]-p0[1];
303 Double_t z12=p1[2]-p0[2];
305 Double_t den= x12*x12*sigmasq[1]*sigmasq[2]+y12*y12*sigmasq[0]*sigmasq[2]+z12*z12*sigmasq[0]*sigmasq[1];
307 Double_t kk= 2*(x12*x12/sigmasq[0]+y12*y12/sigmasq[1]+z12*z12/sigmasq[2]);
310 cc[0]=-x12/sigmasq[0];
311 cc[1]=-y12/sigmasq[1];
312 cc[2]=-z12/sigmasq[2];
314 Double_t ww=(-p0[0]*x12*sigmasq[1]*sigmasq[2]-p0[1]*y12*sigmasq[0]*sigmasq[2]-p0[2]*z12*sigmasq[0]*sigmasq[1])/den;
316 Double_t ss= -p0[0]*cc[0]-p0[1]*cc[1]-p0[2]*cc[2];
319 aa[0]=x12*sigmasq[1]*sigmasq[2]/den;
320 aa[1]=y12*sigmasq[0]*sigmasq[2]/den;
321 aa[2]=z12*sigmasq[0]*sigmasq[1]/den;
323 m[0][0]=aa[0]*(aa[0]*kk+2*cc[0])+2*cc[0]*aa[0]+2/sigmasq[0];
324 m[0][1]=aa[1]*(aa[0]*kk+2*cc[0])+2*cc[1]*aa[0];
325 m[0][2]=aa[2]*(aa[0]*kk+2*cc[0])+2*cc[2]*aa[0];
327 m[1][0]=aa[0]*(aa[1]*kk+2*cc[1])+2*cc[0]*aa[1];
328 m[1][1]=aa[1]*(aa[1]*kk+2*cc[1])+2*cc[1]*aa[1]+2/sigmasq[1];
329 m[1][2]=aa[2]*(aa[1]*kk+2*cc[1])+2*cc[2]*aa[1];
331 m[2][0]=aa[0]*(aa[2]*kk+2*cc[2])+2*cc[0]*aa[2];
332 m[2][1]=aa[1]*(aa[2]*kk+2*cc[2])+2*cc[1]*aa[2];
333 m[2][2]=aa[2]*(aa[2]*kk+2*cc[2])+2*cc[2]*aa[2]+2/sigmasq[2];
335 d[0]=-ww*(aa[0]*kk+2*cc[0])-2*ss*aa[0]+2*p0[0]/sigmasq[0];
336 d[1]=-ww*(aa[1]*kk+2*cc[1])-2*ss*aa[1]+2*p0[1]/sigmasq[1];
337 d[2]=-ww*(aa[2]*kk+2*cc[2])-2*ss*aa[2]+2*p0[2]/sigmasq[2];
340 //--------------------------------------------------------------------------
341 Double_t AliVertexerTracks::GetStrLinMinDist(Double_t *p0,Double_t *p1,Double_t *x0)
344 Double_t x12=p0[0]-p1[0];
345 Double_t y12=p0[1]-p1[1];
346 Double_t z12=p0[2]-p1[2];
347 Double_t x10=p0[0]-x0[0];
348 Double_t y10=p0[1]-x0[1];
349 Double_t z10=p0[2]-x0[2];
350 return ((x10*x10+y10*y10+z10*z10)*(x12*x12+y12*y12+z12*z12)-(x10*x12+y10*y12+z10*z12)*(x10*x12+y10*y12+z10*z12))/(x12*x12+y12*y12+z12*z12);
352 //---------------------------------------------------------------------------
353 void AliVertexerTracks::OneTrackVertFinder()
355 // find vertex for events with 1 track, using DCA to nominal beam axis
356 if(fDebug) printf("Number of prepared tracks =%d - Call OneTrackVertFinder",fTrkArray.GetEntries());
358 track1 = (AliESDtrack*)fTrkArray.At(0);
359 Double_t field=GetField();
360 Double_t alpha=track1->GetAlpha();
361 Double_t mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
362 Double_t pos[3],dir[3];
363 track1->GetXYZAt(mindist,field,pos);
364 track1->GetPxPyPzAt(mindist,field,dir);
365 AliStrLine *line1 = new AliStrLine(pos,dir);
366 Double_t p1[3]={fNominalPos[0],fNominalPos[1],0.};
367 Double_t p2[3]={fNominalPos[0],fNominalPos[1],10.};
368 AliStrLine *zeta=new AliStrLine(p1,p2,kTRUE);
369 Double_t crosspoint[3]={0.,0.,0.};
372 Int_t retcode = zeta->Cross(line1,crosspoint);
374 sigma=line1->GetDistFromPoint(crosspoint);
379 fVert.SetXYZ(crosspoint);
380 fVert.SetDispersion(sigma);
381 fVert.SetNContributors(nContrib);
383 //---------------------------------------------------------------------------
384 void AliVertexerTracks::HelixVertexFinder()
386 // Get estimate of vertex position in (x,y) from tracks DCA
391 for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
392 Double_t field=GetField();
394 Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
396 Double_t aver[3]={0.,0.,0.};
397 Double_t averquad[3]={0.,0.,0.};
398 Double_t sigmaquad[3]={0.,0.,0.};
405 Double_t alpha, cs, sn;
406 Double_t crosspoint[3];
407 for(Int_t i=0; i<nacc; i++){
408 track1 = (AliESDtrack*)fTrkArray.At(i);
411 for(Int_t j=i+1; j<nacc; j++){
412 track2 = (AliESDtrack*)fTrkArray.At(j);
414 distCA=track2->PropagateToDCA(track1,field);
415 if(fDCAcut<=0 ||(fDCAcut>0&&distCA<fDCAcut)){
416 track1->GetExternalParameters(x,par);
417 alpha=track1->GetAlpha();
418 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
419 Double_t x1=x*cs - par[0]*sn;
420 Double_t y1=x*sn + par[0]*cs;
422 Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2();
423 track2->GetExternalParameters(x,par);
424 alpha=track2->GetAlpha();
425 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
426 Double_t x2=x*cs - par[0]*sn;
427 Double_t y2=x*sn + par[0]*cs;
429 Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2();
430 Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2();
431 Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
432 Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
433 Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
434 crosspoint[0]=wx1*x1 + wx2*x2;
435 crosspoint[1]=wy1*y1 + wy2*y2;
436 crosspoint[2]=wz1*z1 + wz2*z2;
439 for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
440 for(Int_t jj=0;jj<3;jj++)averquad[jj]+=(crosspoint[jj]*crosspoint[jj]);
446 for(Int_t jj=0;jj<3;jj++){
447 initPos[jj] = aver[jj]/ncombi;
448 averquad[jj]/=ncombi;
449 sigmaquad[jj]=averquad[jj]-initPos[jj]*initPos[jj];
450 sigma+=sigmaquad[jj];
452 sigma=TMath::Sqrt(TMath::Abs(sigma));
455 Warning("HelixVertexFinder","Finder did not succed");
458 fVert.SetXYZ(initPos);
459 fVert.SetDispersion(sigma);
460 fVert.SetNContributors(ncombi);
462 //----------------------------------------------------------------------------
463 Int_t AliVertexerTracks::PrepareTracks(TTree &trkTree,Int_t optImpParCut)
466 // Propagate tracks to initial vertex position and store them in a TObjArray
469 Double_t maxd0z0 = fMaxd0z0; // default is 5 mm
471 Double_t sigmaCurr[3];
472 Double_t normdistx,normdisty;
473 Float_t d0z0[2],covd0z0[3];
475 Double_t field=GetField();
477 AliESDVertex *initVertex = new AliESDVertex(fNominalPos,fNominalCov,1,1);
479 Int_t nEntries = (Int_t)trkTree.GetEntries();
480 if(!fTrkArray.IsEmpty()) fTrkArray.Delete();
483 printf(" PrepareTracks()\n");
486 for(Int_t i=0; i<nEntries; i++) {
487 AliESDtrack *track = new AliESDtrack;
488 trkTree.SetBranchAddress("tracks",&track);
491 // propagate track to vertex
492 if(optImpParCut<=1 || fOnlyFitter) { // optImpParCut==1 or 0
493 track->RelateToVertex(initVertex,field,100.);
494 } else { // optImpParCut==2
495 fCurrentVertex->GetSigmaXYZ(sigmaCurr);
496 normdistx = TMath::Abs(fCurrentVertex->GetXv()-fNominalPos[0])/TMath::Sqrt(sigmaCurr[0]*sigmaCurr[0]+fNominalCov[0]);
497 normdisty = TMath::Abs(fCurrentVertex->GetYv()-fNominalPos[1])/TMath::Sqrt(sigmaCurr[1]*sigmaCurr[1]+fNominalCov[1]);
498 if(normdistx < 3. && normdisty < 3. &&
499 (sigmaCurr[0]+sigmaCurr[1])<(TMath::Sqrt(fNominalCov[0])+TMath::Sqrt(fNominalCov[1]))) {
500 track->RelateToVertex(fCurrentVertex,field,100.);
502 track->RelateToVertex(initVertex,field,100.);
506 track->GetImpactParameters(d0z0,covd0z0);
507 sigma = TMath::Sqrt(covd0z0[0]);
508 maxd0rphi = fNSigma*sigma;
509 if(optImpParCut==1) maxd0rphi *= 5.;
513 if(fDebug) printf("trk %d; lab %d; |d0| = %f; d0 cut = %f; |z0| = %f; |d0|oplus|z0| = %f; d0z0 cut = %f\n",i,track->GetLabel(),TMath::Abs(d0z0[0]),maxd0rphi,TMath::Abs(d0z0[1]),TMath::Sqrt(d0z0[0]*d0z0[0]+d0z0[1]*d0z0[1]),maxd0z0);
515 // during iteration 1, if fConstraint=kFALSE,
516 // select tracks with d0oplusz0 < maxd0z0
517 if(optImpParCut==1 && !fConstraint && nEntries>=3 &&
518 fVert.GetNContributors()>0) {
519 if(TMath::Sqrt(d0z0[0]*d0z0[0]+d0z0[1]*d0z0[1]) > maxd0z0) {
520 if(fDebug) printf(" rejected\n");
521 delete track; continue;
525 // select tracks with d0rphi < maxd0rphi
526 if(optImpParCut>0 && TMath::Abs(d0z0[0]) > maxd0rphi) {
527 if(fDebug) printf(" rejected\n");
528 delete track; continue;
531 fTrkArray.AddLast(track);
539 //---------------------------------------------------------------------------
540 AliESDVertex* AliVertexerTracks::RemoveTracksFromVertex(AliESDVertex *inVtx,
545 // Removes tracks in trksTree from fit of inVtx
548 if(!strstr(inVtx->GetTitle(),"VertexerTracksWithConstraintOnlyFitter"))
549 printf("WARNING: result of tracks' removal will be only approximately correct\n");
552 rv(0,0) = inVtx->GetXv();
553 rv(1,0) = inVtx->GetYv();
554 rv(2,0) = inVtx->GetZv();
557 inVtx->GetCovMatrix(cov);
559 vV(0,1) = cov[1]; vV(1,0) = cov[1];
561 vV(0,2) = cov[3]; vV(2,0) = cov[3];
562 vV(1,2) = cov[4]; vV(2,1) = cov[4];
565 TMatrixD sumWi(TMatrixD::kInverted,vV);
566 TMatrixD sumWiri(sumWi,TMatrixD::kMult,rv);
568 Int_t nUsedTrks = inVtx->GetNContributors();
569 Double_t chi2 = inVtx->GetChi2();
571 AliESDtrack *track = 0;
572 trksTree->SetBranchAddress("tracks",&track);
573 Int_t ntrks = trksTree->GetEntries();
574 for(Int_t i=0;i<ntrks;i++) {
575 trksTree->GetEvent(i);
576 if(!inVtx->UsesTrack(track->GetID())) {
577 printf("track %d was not used in vertex fit\n",track->GetID());
580 Double_t alpha = track->GetAlpha();
581 Double_t xl = diamondxy[0]*TMath::Cos(alpha)+diamondxy[1]*TMath::Sin(alpha);
582 track->AliExternalTrackParam::PropagateTo(xl,AliTracker::GetBz());
583 // vector of track global coordinates
585 // covariance matrix of ri
588 // get space point from track
589 if(!TrackToPoint(track,ri,wWi)) continue;
591 TMatrixD wWiri(wWi,TMatrixD::kMult,ri);
597 TMatrixD deltar = rv; deltar -= ri;
598 TMatrixD wWideltar(wWi,TMatrixD::kMult,deltar);
599 Double_t chi2i = deltar(0,0)*wWideltar(0,0)+
600 deltar(1,0)*wWideltar(1,0)+
601 deltar(2,0)*wWideltar(2,0);
602 // remove from total chi2
607 printf("Trying to remove too many tracks!\n");
615 // new inverted of weights matrix
616 TMatrixD invsumWi(TMatrixD::kInverted,sumWi);
618 // new position of primary vertex
619 rvnew.Mult(vVnew,sumWiri);
621 Double_t position[3];
622 position[0] = rvnew(0,0);
623 position[1] = rvnew(1,0);
624 position[2] = rvnew(2,0);
632 // store data in the vertex object
633 AliESDVertex *outVtx = new AliESDVertex(position,cov,chi2,nUsedTrks);
634 outVtx->SetTitle(inVtx->GetTitle());
636 inVtx->GetTruePos(tp);
637 outVtx->SetTruePos(tp);
638 UShort_t *inindices = inVtx->GetIndices();
639 UShort_t *outindices = new UShort_t[outVtx->GetNContributors()];
642 for(Int_t k=0; k<inVtx->GetNIndices(); k++) {
644 for(Int_t l=0; l<ntrks; l++) {
645 trksTree->GetEvent(l);
646 if(inindices[k]==track->GetID()) copyindex=kFALSE;
649 outindices[j] = inindices[k]; j++;
652 outVtx->SetIndices(outVtx->GetNContributors(),outindices);
653 delete [] outindices;
656 printf("Vertex before removing tracks:\n");
657 inVtx->PrintStatus();
658 inVtx->PrintIndices();
659 printf("Vertex after removing tracks:\n");
660 outVtx->PrintStatus();
661 outVtx->PrintIndices();
666 //---------------------------------------------------------------------------
667 void AliVertexerTracks::SetSkipTracks(Int_t n,Int_t *skipped)
670 // Mark the tracks not to be used in the vertex reconstruction.
671 // Tracks are identified by AliESDtrack::GetID()
673 fNTrksToSkip = n; fTrksToSkip = new Int_t[n];
674 for(Int_t i=0;i<n;i++) fTrksToSkip[i] = skipped[i];
677 //---------------------------------------------------------------------------
678 void AliVertexerTracks::SetVtxStart(AliESDVertex *vtx)
681 // Set initial vertex knowledge
683 vtx->GetXYZ(fNominalPos);
684 vtx->GetCovMatrix(fNominalCov);
688 //---------------------------------------------------------------------------
689 void AliVertexerTracks::StrLinVertexFinderMinDist(Int_t optUseWeights)
691 // Calculate the point at minimum distance to prepared tracks
696 for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
697 const Int_t knacc = (Int_t)fTrkArray.GetEntriesFast();
698 Double_t field=GetField();
701 Double_t (*vectP0)[3]=new Double_t [knacc][3];
702 Double_t (*vectP1)[3]=new Double_t [knacc][3];
705 Double_t dsum[3]={0,0,0};
706 for(Int_t i=0;i<3;i++)
707 for(Int_t j=0;j<3;j++)sum[i][j]=0;
708 for(Int_t i=0; i<knacc; i++){
709 track1 = (AliESDtrack*)fTrkArray.At(i);
710 Double_t alpha=track1->GetAlpha();
711 Double_t mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
712 Double_t pos[3],dir[3];
713 track1->GetXYZAt(mindist,field,pos);
714 track1->GetPxPyPzAt(mindist,field,dir);
715 AliStrLine *line1 = new AliStrLine(pos,dir);
716 // AliStrLine *line1 = new AliStrLine();
717 // track1->ApproximateHelixWithLine(mindist,field,line1);
719 Double_t p0[3],cd[3];
722 Double_t p1[3]={p0[0]+cd[0],p0[1]+cd[1],p0[2]+cd[2]};
732 if(optUseWeights==0)GetStrLinDerivMatrix(p0,p1,matr,dknow);
733 if(optUseWeights==1){
735 sigmasq[0]=track1->GetSigmaY2();
736 sigmasq[1]=track1->GetSigmaY2();
737 sigmasq[2]=track1->GetSigmaZ2();
738 GetStrLinDerivMatrix(p0,p1,sigmasq,matr,dknow);
741 for(Int_t iii=0;iii<3;iii++){
742 dsum[iii]+=dknow[iii];
743 for(Int_t lj=0;lj<3;lj++) sum[iii][lj]+=matr[iii][lj];
749 Double_t det=GetDeterminant3X3(sum);
752 for(Int_t zz=0;zz<3;zz++){
753 for(Int_t ww=0;ww<3;ww++){
754 for(Int_t kk=0;kk<3;kk++) vett[ww][kk]=sum[ww][kk];
756 for(Int_t kk=0;kk<3;kk++) vett[kk][zz]=dsum[kk];
757 initPos[zz]=GetDeterminant3X3(vett)/det;
761 for(Int_t i=0; i<knacc; i++){
762 Double_t p0[3]={0,0,0},p1[3]={0,0,0};
763 for(Int_t ii=0;ii<3;ii++){
764 p0[ii]=vectP0[i][ii];
765 p1[ii]=vectP1[i][ii];
767 sigma+=GetStrLinMinDist(p0,p1,initPos);
770 sigma=TMath::Sqrt(sigma);
772 Warning("StrLinVertexFinderMinDist","Finder did not succed");
777 fVert.SetXYZ(initPos);
778 fVert.SetDispersion(sigma);
779 fVert.SetNContributors(knacc);
781 //---------------------------------------------------------------------------
782 Bool_t AliVertexerTracks::TrackToPoint(AliESDtrack *t,
783 TMatrixD &ri,TMatrixD &wWi) const
786 // Extract from the AliESDtrack the global coordinates ri and covariance matrix
787 // wWi of the space point that it represents (to be used in VertexFitter())
791 Double_t rotAngle = t->GetAlpha();
792 if(rotAngle<0.) rotAngle += 2.*TMath::Pi();
793 Double_t cosRot = TMath::Cos(rotAngle);
794 Double_t sinRot = TMath::Sin(rotAngle);
796 ri(0,0) = t->GetX()*cosRot-t->GetY()*sinRot;
797 ri(1,0) = t->GetX()*sinRot+t->GetY()*cosRot;
800 // matrix to go from global (x,y,z) to local (y,z);
809 // covariance matrix of local (y,z) - inverted
812 t->GetExternalCovariance(cc);
817 if(uUi.Determinant() <= 0.) return kFALSE;
818 TMatrixD uUiInv(TMatrixD::kInverted,uUi);
820 // weights matrix: wWi = qQiT * uUiInv * qQi
821 TMatrixD uUiInvQi(uUiInv,TMatrixD::kMult,qQi);
822 TMatrixD m(qQi,TMatrixD::kTransposeMult,uUiInvQi);
827 //---------------------------------------------------------------------------
828 void AliVertexerTracks::TooFewTracks(const AliESD* esdEvent)
831 // When the number of tracks is < fMinTracks
834 // deal with vertices not found
835 Double_t pos[3],err[3];
837 pos[0] = fNominalPos[0];
838 err[0] = TMath::Sqrt(fNominalCov[0]);
839 pos[1] = fNominalPos[1];
840 err[1] = TMath::Sqrt(fNominalCov[2]);
841 pos[2] = esdEvent->GetVertex()->GetZv();
842 err[2] = esdEvent->GetVertex()->GetZRes();
843 if(err[0]>1. && esdEvent->GetVertex()->GetNContributors()<=0)
844 ncontr = -1; // (x,y,z) = (0,0,0)
845 if(err[0]>1. && esdEvent->GetVertex()->GetNContributors()>0)
846 ncontr = -2; // (x,y,z) = (0,0,z_fromSPD)
847 if(err[0]<1. && esdEvent->GetVertex()->GetNContributors()<=0)
848 ncontr = -3; // (x,y,z) = (x_mean,y_mean,0)
849 if(err[0]<1. && esdEvent->GetVertex()->GetNContributors()>0)
850 ncontr = -4; // (x,y,z) = (x_mean,y_mean,z_fromSPD)
852 fCurrentVertex = new AliESDVertex(pos,err);
853 fCurrentVertex->SetNContributors(ncontr);
856 esdEvent->GetVertex()->GetTruePos(tp);
857 fCurrentVertex->SetTruePos(tp);
859 fCurrentVertex->SetTitle("VertexerTracksWithConstraint");
861 fCurrentVertex->SetTitle("VertexerTracksNoConstraint");
866 //---------------------------------------------------------------------------
867 void AliVertexerTracks::VertexFinder(Int_t optUseWeights)
870 // Get estimate of vertex position in (x,y) from tracks DCA
874 for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
875 Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
876 Double_t aver[3]={0.,0.,0.};
877 Double_t aversq[3]={0.,0.,0.};
878 Double_t sigmasq[3]={0.,0.,0.};
883 Double_t pos[3],dir[3];
884 Double_t alpha,mindist;
885 Double_t field=GetField();
887 for(Int_t i=0; i<nacc; i++){
888 track1 = (AliESDtrack*)fTrkArray.At(i);
889 alpha=track1->GetAlpha();
890 mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
891 track1->GetXYZAt(mindist,field,pos);
892 track1->GetPxPyPzAt(mindist,field,dir);
893 AliStrLine *line1 = new AliStrLine(pos,dir);
895 // AliStrLine *line1 = new AliStrLine();
896 // track1->ApproximateHelixWithLine(mindist,field,line1);
898 for(Int_t j=i+1; j<nacc; j++){
899 track2 = (AliESDtrack*)fTrkArray.At(j);
900 alpha=track2->GetAlpha();
901 mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
902 track2->GetXYZAt(mindist,field,pos);
903 track2->GetPxPyPzAt(mindist,field,dir);
904 AliStrLine *line2 = new AliStrLine(pos,dir);
905 // AliStrLine *line2 = new AliStrLine();
906 // track2->ApproximateHelixWithLine(mindist,field,line2);
907 Double_t distCA=line2->GetDCA(line1);
908 //printf("%d %d %f\n",i,j,distCA);
909 if(fDCAcut<=0 || (fDCAcut>0&&distCA<fDCAcut)){
910 Double_t pnt1[3],pnt2[3],crosspoint[3];
912 if(optUseWeights<=0){
913 Int_t retcode = line2->Cross(line1,crosspoint);
916 for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
917 for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]);
921 Int_t retcode = line1->CrossPoints(line2,pnt1,pnt2);
924 alpha=track1->GetAlpha();
925 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
926 Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2();
927 alpha=track2->GetAlpha();
928 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
929 Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2();
930 Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2();
931 Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
932 Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
933 Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
934 crosspoint[0]=wx1*pnt1[0] + wx2*pnt2[0];
935 crosspoint[1]=wy1*pnt1[1] + wy2*pnt2[1];
936 crosspoint[2]=wz1*pnt1[2] + wz2*pnt2[2];
939 for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
940 for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]);
949 for(Int_t jj=0;jj<3;jj++){
950 initPos[jj] = aver[jj]/ncombi;
951 //printf("%f\n",initPos[jj]);
953 sigmasq[jj]=aversq[jj]-initPos[jj]*initPos[jj];
956 sigma=TMath::Sqrt(TMath::Abs(sigma));
959 Warning("VertexFinder","Finder did not succed");
962 fVert.SetXYZ(initPos);
963 fVert.SetDispersion(sigma);
964 fVert.SetNContributors(ncombi);
966 //---------------------------------------------------------------------------
967 void AliVertexerTracks::VertexFitter(Bool_t useConstraint)
970 // The optimal estimate of the vertex position is given by a "weighted
971 // average of tracks positions"
972 // Original method: V. Karimaki, CMS Note 97/0051
975 fVert.GetXYZ(initPos);
976 Int_t arrEntries = (Int_t)fTrkArray.GetEntries();
977 if(arrEntries==1) useConstraint=kTRUE;
979 printf(" VertexFitter(): start\n");
980 printf(" Number of tracks in array: %d\n",(Int_t)fTrkArray.GetEntriesFast());
981 printf(" Minimum # tracks required in fit: %d\n",fMinTracks);
982 printf("Vertex position after finder: %f,%f,%f\n",initPos[0],initPos[1],initPos[2]);
983 if(useConstraint) printf(" This vertex will be used in fit: (%f+-%f,%f+-%f)\n",fNominalPos[0],TMath::Sqrt(fNominalCov[0]),fNominalPos[1],TMath::Sqrt(fNominalCov[2]));
989 rv(0,0) = initPos[0];
990 rv(1,0) = initPos[1];
992 Double_t xlStart,alpha;
994 Double_t chi2,chi2i,chi2b;
998 // initial vertex covariance matrix
1000 vVb(0,0) = fNominalCov[0];
1001 vVb(0,1) = fNominalCov[1];
1003 vVb(1,0) = fNominalCov[1];
1004 vVb(1,1) = fNominalCov[2];
1008 vVb(2,2) = fNominalCov[5];
1009 TMatrixD vVbInv(TMatrixD::kInverted,vVb);
1011 rb(0,0) = fNominalPos[0];
1012 rb(1,0) = fNominalPos[1];
1013 rb(2,0) = fNominalPos[2];
1014 TMatrixD vVbInvrb(vVbInv,TMatrixD::kMult,rb);
1018 // 1st - estimate of vtx using all tracks
1019 // 2nd - estimate of global chi2
1020 for(step=0; step<2; step++) {
1021 if(fDebug) printf(" step = %d\n",step);
1025 TMatrixD sumWiri(3,1);
1026 TMatrixD sumWi(3,3);
1027 for(i=0; i<3; i++) {
1029 for(j=0; j<3; j++) sumWi(j,i) = 0.;
1032 // mean vertex constraint
1035 sumWiri(i,0) += vVbInvrb(i,0);
1036 for(k=0;k<3;k++) sumWi(i,k) += vVbInv(i,k);
1039 TMatrixD deltar = rv; deltar -= rb;
1040 TMatrixD vVbInvdeltar(vVbInv,TMatrixD::kMult,deltar);
1041 chi2b = deltar(0,0)*vVbInvdeltar(0,0)+
1042 deltar(1,0)*vVbInvdeltar(1,0)+
1043 deltar(2,0)*vVbInvdeltar(2,0);
1049 for(k=0; k<arrEntries; k++) {
1050 // get track from track array
1051 t = (AliESDtrack*)fTrkArray.At(k);
1052 alpha = t->GetAlpha();
1053 xlStart = initPos[0]*TMath::Cos(alpha)+initPos[1]*TMath::Sin(alpha);
1054 // to vtxSeed (from finder)
1055 t->AliExternalTrackParam::PropagateTo(xlStart,AliTracker::GetBz());
1058 // vector of track global coordinates
1060 // covariance matrix of ri
1063 // get space point from track
1064 if(!TrackToPoint(t,ri,wWi)) continue;
1066 TMatrixD wWiri(wWi,TMatrixD::kMult,ri);
1069 TMatrixD deltar = rv; deltar -= ri;
1070 TMatrixD wWideltar(wWi,TMatrixD::kMult,deltar);
1071 chi2i = deltar(0,0)*wWideltar(0,0)+
1072 deltar(1,0)*wWideltar(1,0)+
1073 deltar(2,0)*wWideltar(2,0);
1075 // add to total chi2
1082 } // end loop on tracks
1084 if(nUsedTrks < fMinTracks) {
1089 Double_t determinant = sumWi.Determinant();
1090 //cerr<<" determinant: "<<determinant<<endl;
1091 if(determinant < 100.) {
1092 printf("det(V) = 0\n");
1097 // inverted of weights matrix
1098 TMatrixD invsumWi(TMatrixD::kInverted,sumWi);
1101 // position of primary vertex
1102 rv.Mult(vV,sumWiri);
1104 } // end loop on the 2 steps
1109 if(fDebug) printf("TooFewTracks\n");
1110 fCurrentVertex = new AliESDVertex(0.,0.,-1);
1114 Double_t position[3];
1115 position[0] = rv(0,0);
1116 position[1] = rv(1,0);
1117 position[2] = rv(2,0);
1118 Double_t covmatrix[6];
1119 covmatrix[0] = vV(0,0);
1120 covmatrix[1] = vV(0,1);
1121 covmatrix[2] = vV(1,1);
1122 covmatrix[3] = vV(0,2);
1123 covmatrix[4] = vV(1,2);
1124 covmatrix[5] = vV(2,2);
1126 // store data in the vertex object
1127 fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nUsedTrks);
1130 printf(" VertexFitter(): finish\n");
1131 printf(" rv = ( %f , %f , %f )\n\n",rv(0,0),rv(1,0),rv(2,0));
1132 fCurrentVertex->PrintStatus();
1137 //----------------------------------------------------------------------------
1138 AliESDVertex* AliVertexerTracks::VertexForSelectedTracks(TTree *trkTree,Bool_t optUseFitter,Bool_t optPropagate)
1141 // Return vertex from tracks in trkTree
1146 // get tracks and propagate them to initial vertex position
1147 Int_t nTrksPrep = PrepareTracks(*trkTree,0);
1148 if(nTrksPrep < TMath::Max(2,fMinTracks) ) {
1149 if(fDebug) printf("TooFewTracks\n");
1150 fCurrentVertex = new AliESDVertex(0.,0.,-1);
1151 return fCurrentVertex;
1155 case 1: StrLinVertexFinderMinDist(1); break;
1156 case 2: StrLinVertexFinderMinDist(0); break;
1157 case 3: HelixVertexFinder(); break;
1158 case 4: VertexFinder(1); break;
1159 case 5: VertexFinder(0); break;
1160 default: printf("Wrong algorithm\n"); break;
1163 if(fDebug) printf(" Vertex finding completed\n");
1167 VertexFitter(fConstraint);
1168 if(fDebug) printf(" Vertex fit completed\n");
1170 Double_t position[3]={fVert.GetXv(),fVert.GetYv(),fVert.GetZv()};
1171 Double_t covmatrix[6]={0.,0.,0.,0.,0.,0.};
1172 Double_t chi2=99999.;
1173 Int_t nUsedTrks=fVert.GetNContributors();
1174 fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nUsedTrks);
1176 fCurrentVertex->SetDispersion(fVert.GetDispersion());
1179 // set indices of used tracks and propagate track to found vertex
1180 UShort_t *indices = 0;
1181 AliESDtrack *eta = 0;
1182 if(fCurrentVertex->GetNContributors()>0) {
1183 indices = new UShort_t[fCurrentVertex->GetNContributors()];
1184 for(Int_t jj=0;jj<(Int_t)fTrkArray.GetEntriesFast();jj++) {
1185 eta = (AliESDtrack*)fTrkArray.At(jj);
1186 indices[jj] = (UShort_t)eta->GetID();
1187 if(optPropagate&&optUseFitter){
1188 if(TMath::Sqrt(fCurrentVertex->GetXv()*fCurrentVertex->GetXv()+fCurrentVertex->GetYv()*fCurrentVertex->GetYv())<3.) {
1189 eta->RelateToVertex(fCurrentVertex,GetField(),100.);
1190 if(fDebug) printf("Track %d propagated to found vertex\n",jj);
1192 AliWarning("Found vertex outside beam pipe!");
1196 fCurrentVertex->SetIndices(fCurrentVertex->GetNContributors(),indices);
1201 return fCurrentVertex;
1203 //----------------------------------------------------------------------------
1204 AliESDVertex* AliVertexerTracks::VertexForSelectedTracks(TObjArray *trkArray,Bool_t optUseFitter, Bool_t optPropagate)
1207 // Return vertex from array of tracks
1210 // get tracks and propagate them to initial vertex position
1211 Int_t nTrks = trkArray->GetEntriesFast();
1212 if(nTrks < TMath::Max(2,fMinTracks) ) {
1213 if(fDebug) printf("TooFewTracks\n");
1214 fCurrentVertex = new AliESDVertex(0.,0.,-1);
1215 return fCurrentVertex;
1217 TTree *trkTree = new TTree("TreeT","tracks");
1218 AliESDtrack *esdTrack = 0;
1219 trkTree->Branch("tracks","AliESDtrack",&esdTrack);
1220 for(Int_t i=0; i<nTrks; i++){
1221 esdTrack = (AliESDtrack*)trkArray->At(i);
1225 AliESDVertex *vtx = VertexForSelectedTracks(trkTree,optUseFitter,optPropagate);
1229 //--------------------------------------------------------------------------