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();
125 if(fDebug) printf("TooFewTracks\n");
126 TooFewTracks(esdEvent);
127 return fCurrentVertex;
130 TTree *trkTree = new TTree("TreeT","tracks");
131 AliESDtrack *esdTrack = 0;
132 trkTree->Branch("tracks","AliESDtrack",&esdTrack);
135 for(Int_t i=0; i<nTrksTot; i++) {
136 AliESDtrack *et = esdEvent->GetTrack(i);
137 esdTrack = new AliESDtrack(*et);
138 // check tracks to skip
140 for(Int_t j=0; j<fNTrksToSkip; j++) {
141 if(et->GetID()==fTrksToSkip[j]) {
142 if(fDebug) printf("skipping track: %d\n",i);
146 if(skipThis) {delete esdTrack;continue;}
148 if(!(esdTrack->GetStatus()&AliESDtrack::kITSin)) {delete esdTrack;continue;}
149 if(fITSrefit && !(esdTrack->GetStatus()&AliESDtrack::kITSrefit)) {delete esdTrack;continue;}
150 Int_t nclus=esdTrack->GetNcls(0); // check number of clusters in ITS
151 if(nclus<fMinITSClusters) {delete esdTrack;continue;}
157 // If fConstraint=kFALSE
158 // run VertexFinder(1) to get rough estimate of initVertex (x,y)
160 // fill fTrkArray, for VertexFinder()
161 if(!fTrkArray.IsEmpty()) fTrkArray.Delete();
162 PrepareTracks(*trkTree,0);
163 Double_t cutsave = fDCAcut; fDCAcut = 0.1; // 1 mm
164 VertexFinder(1); // using weights, cutting dca < fDCAcut
167 if(fVert.GetNContributors()>0) {
168 fVert.GetXYZ(fNominalPos);
169 fNominalPos[0] = fVert.GetXv();
170 fNominalPos[1] = fVert.GetYv();
171 fNominalPos[2] = fVert.GetZv();
172 if(fDebug) printf("No mean vertex: VertexFinder gives (%f, %f, %f)\n",fNominalPos[0],fNominalPos[1],fNominalPos[2]);
177 if(fDebug) printf("No mean vertex and VertexFinder failed\n");
184 // propagate tracks to fNominalPos vertex
186 // if(constraint) reject for |d0|>5*fNSigma*sigma w.r.t. fNominal... vertex
187 // else reject for |d0|\oplus|z0| > 5 mm w.r.t. fNominal... vertex
189 // propagate tracks to best between initVertex and fCurrentVertex
190 // preselect tracks (reject for |d0|>fNSigma*sigma w.r.t. best
191 // between initVertex and fCurrentVertex)
192 for(Int_t iter=0; iter<2; iter++) {
193 if(fOnlyFitter && iter==0) continue;
194 Int_t nTrksPrep = PrepareTracks(*trkTree,iter+1);
195 if(fDebug) printf(" tracks prepared - iteration %d: %d\n",iter+1,nTrksPrep);
196 if(nTrksPrep < fMinTracks) {
197 if(fDebug) printf("TooFewTracks\n");
198 TooFewTracks(esdEvent);
199 if(fDebug) fCurrentVertex->PrintStatus();
202 return fCurrentVertex;
208 if(fDebug) printf("Just one track\n");
209 OneTrackVertFinder();
212 case 1: StrLinVertexFinderMinDist(1); break;
213 case 2: StrLinVertexFinderMinDist(0); break;
214 case 3: HelixVertexFinder(); break;
215 case 4: VertexFinder(1); break;
216 case 5: VertexFinder(0); break;
217 default: printf("Wrong algorithm\n"); break;
220 if(fDebug) printf(" Vertex finding completed\n");
224 VertexFitter(fConstraint);
225 if(fDebug) printf(" Vertex fit completed\n");
226 if(iter==0) fTrkArray.Delete();
227 } // end loop on the two iterations
230 // take true pos from SPD vertex in ESD and write it in tracks' vertex
232 esdEvent->GetVertex()->GetTruePos(tp);
233 fCurrentVertex->SetTruePos(tp);
237 fCurrentVertex->SetTitle("VertexerTracksWithConstraintOnlyFitter");
239 fCurrentVertex->SetTitle("VertexerTracksWithConstraint");
242 fCurrentVertex->SetTitle("VertexerTracksNoConstraint");
246 // set indices of used tracks
247 UShort_t *indices = 0;
248 AliESDtrack *ett = 0;
249 if(fCurrentVertex->GetNContributors()>0) {
250 indices = new UShort_t[fCurrentVertex->GetNContributors()];
251 for(Int_t jj=0;jj<(Int_t)fTrkArray.GetEntriesFast();jj++) {
252 ett = (AliESDtrack*)fTrkArray.At(jj);
253 indices[jj] = (UShort_t)ett->GetID();
255 fCurrentVertex->SetIndices(fCurrentVertex->GetNContributors(),indices);
263 if(fTrksToSkip) delete [] fTrksToSkip;
266 if(fDebug) fCurrentVertex->PrintStatus();
267 if(fDebug) fCurrentVertex->PrintIndices();
270 return fCurrentVertex;
272 //------------------------------------------------------------------------
273 Double_t AliVertexerTracks::GetDeterminant3X3(Double_t matr[][3])
276 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];
279 //-------------------------------------------------------------------------
280 void AliVertexerTracks::GetStrLinDerivMatrix(Double_t *p0,Double_t *p1,Double_t (*m)[3],Double_t *d)
283 Double_t x12=p0[0]-p1[0];
284 Double_t y12=p0[1]-p1[1];
285 Double_t z12=p0[2]-p1[2];
286 Double_t kk=x12*x12+y12*y12+z12*z12;
287 m[0][0]=2-2/kk*x12*x12;
288 m[0][1]=-2/kk*x12*y12;
289 m[0][2]=-2/kk*x12*z12;
290 m[1][0]=-2/kk*x12*y12;
291 m[1][1]=2-2/kk*y12*y12;
292 m[1][2]=-2/kk*y12*z12;
293 m[2][0]=-2/kk*x12*z12;
295 m[2][2]=2-2/kk*z12*z12;
296 d[0]=2*p0[0]-2/kk*p0[0]*x12*x12-2/kk*p0[2]*x12*z12-2/kk*p0[1]*x12*y12;
297 d[1]=2*p0[1]-2/kk*p0[1]*y12*y12-2/kk*p0[0]*x12*y12-2/kk*p0[2]*z12*y12;
298 d[2]=2*p0[2]-2/kk*p0[2]*z12*z12-2/kk*p0[0]*x12*z12-2/kk*p0[1]*z12*y12;
301 //--------------------------------------------------------------------------
302 void AliVertexerTracks::GetStrLinDerivMatrix(Double_t *p0,Double_t *p1,Double_t *sigmasq,Double_t (*m)[3],Double_t *d)
305 Double_t x12=p1[0]-p0[0];
306 Double_t y12=p1[1]-p0[1];
307 Double_t z12=p1[2]-p0[2];
309 Double_t den= x12*x12*sigmasq[1]*sigmasq[2]+y12*y12*sigmasq[0]*sigmasq[2]+z12*z12*sigmasq[0]*sigmasq[1];
311 Double_t kk= 2*(x12*x12/sigmasq[0]+y12*y12/sigmasq[1]+z12*z12/sigmasq[2]);
314 cc[0]=-x12/sigmasq[0];
315 cc[1]=-y12/sigmasq[1];
316 cc[2]=-z12/sigmasq[2];
318 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;
320 Double_t ss= -p0[0]*cc[0]-p0[1]*cc[1]-p0[2]*cc[2];
323 aa[0]=x12*sigmasq[1]*sigmasq[2]/den;
324 aa[1]=y12*sigmasq[0]*sigmasq[2]/den;
325 aa[2]=z12*sigmasq[0]*sigmasq[1]/den;
327 m[0][0]=aa[0]*(aa[0]*kk+2*cc[0])+2*cc[0]*aa[0]+2/sigmasq[0];
328 m[0][1]=aa[1]*(aa[0]*kk+2*cc[0])+2*cc[1]*aa[0];
329 m[0][2]=aa[2]*(aa[0]*kk+2*cc[0])+2*cc[2]*aa[0];
331 m[1][0]=aa[0]*(aa[1]*kk+2*cc[1])+2*cc[0]*aa[1];
332 m[1][1]=aa[1]*(aa[1]*kk+2*cc[1])+2*cc[1]*aa[1]+2/sigmasq[1];
333 m[1][2]=aa[2]*(aa[1]*kk+2*cc[1])+2*cc[2]*aa[1];
335 m[2][0]=aa[0]*(aa[2]*kk+2*cc[2])+2*cc[0]*aa[2];
336 m[2][1]=aa[1]*(aa[2]*kk+2*cc[2])+2*cc[1]*aa[2];
337 m[2][2]=aa[2]*(aa[2]*kk+2*cc[2])+2*cc[2]*aa[2]+2/sigmasq[2];
339 d[0]=-ww*(aa[0]*kk+2*cc[0])-2*ss*aa[0]+2*p0[0]/sigmasq[0];
340 d[1]=-ww*(aa[1]*kk+2*cc[1])-2*ss*aa[1]+2*p0[1]/sigmasq[1];
341 d[2]=-ww*(aa[2]*kk+2*cc[2])-2*ss*aa[2]+2*p0[2]/sigmasq[2];
344 //--------------------------------------------------------------------------
345 Double_t AliVertexerTracks::GetStrLinMinDist(Double_t *p0,Double_t *p1,Double_t *x0)
348 Double_t x12=p0[0]-p1[0];
349 Double_t y12=p0[1]-p1[1];
350 Double_t z12=p0[2]-p1[2];
351 Double_t x10=p0[0]-x0[0];
352 Double_t y10=p0[1]-x0[1];
353 Double_t z10=p0[2]-x0[2];
354 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);
356 //---------------------------------------------------------------------------
357 void AliVertexerTracks::OneTrackVertFinder()
359 // find vertex for events with 1 track, using DCA to nominal beam axis
360 if(fDebug) printf("Number of prepared tracks =%d - Call OneTrackVertFinder",fTrkArray.GetEntries());
362 track1 = (AliESDtrack*)fTrkArray.At(0);
363 Double_t field=GetField();
364 Double_t alpha=track1->GetAlpha();
365 Double_t mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
366 Double_t pos[3],dir[3];
367 track1->GetXYZAt(mindist,field,pos);
368 track1->GetPxPyPzAt(mindist,field,dir);
369 AliStrLine *line1 = new AliStrLine(pos,dir);
370 Double_t p1[3]={fNominalPos[0],fNominalPos[1],0.};
371 Double_t p2[3]={fNominalPos[0],fNominalPos[1],10.};
372 AliStrLine *zeta=new AliStrLine(p1,p2,kTRUE);
373 Double_t crosspoint[3]={0.,0.,0.};
376 Int_t retcode = zeta->Cross(line1,crosspoint);
378 sigma=line1->GetDistFromPoint(crosspoint);
383 fVert.SetXYZ(crosspoint);
384 fVert.SetDispersion(sigma);
385 fVert.SetNContributors(nContrib);
387 //---------------------------------------------------------------------------
388 void AliVertexerTracks::HelixVertexFinder()
390 // Get estimate of vertex position in (x,y) from tracks DCA
395 for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
396 Double_t field=GetField();
398 Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
400 Double_t aver[3]={0.,0.,0.};
401 Double_t averquad[3]={0.,0.,0.};
402 Double_t sigmaquad[3]={0.,0.,0.};
409 Double_t alpha, cs, sn;
410 Double_t crosspoint[3];
411 for(Int_t i=0; i<nacc; i++){
412 track1 = (AliESDtrack*)fTrkArray.At(i);
415 for(Int_t j=i+1; j<nacc; j++){
416 track2 = (AliESDtrack*)fTrkArray.At(j);
418 distCA=track2->PropagateToDCA(track1,field);
419 if(fDCAcut<=0 ||(fDCAcut>0&&distCA<fDCAcut)){
420 track1->GetExternalParameters(x,par);
421 alpha=track1->GetAlpha();
422 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
423 Double_t x1=x*cs - par[0]*sn;
424 Double_t y1=x*sn + par[0]*cs;
426 Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2();
427 track2->GetExternalParameters(x,par);
428 alpha=track2->GetAlpha();
429 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
430 Double_t x2=x*cs - par[0]*sn;
431 Double_t y2=x*sn + par[0]*cs;
433 Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2();
434 Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2();
435 Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
436 Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
437 Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
438 crosspoint[0]=wx1*x1 + wx2*x2;
439 crosspoint[1]=wy1*y1 + wy2*y2;
440 crosspoint[2]=wz1*z1 + wz2*z2;
443 for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
444 for(Int_t jj=0;jj<3;jj++)averquad[jj]+=(crosspoint[jj]*crosspoint[jj]);
450 for(Int_t jj=0;jj<3;jj++){
451 initPos[jj] = aver[jj]/ncombi;
452 averquad[jj]/=ncombi;
453 sigmaquad[jj]=averquad[jj]-initPos[jj]*initPos[jj];
454 sigma+=sigmaquad[jj];
456 sigma=TMath::Sqrt(TMath::Abs(sigma));
459 Warning("HelixVertexFinder","Finder did not succed");
462 fVert.SetXYZ(initPos);
463 fVert.SetDispersion(sigma);
464 fVert.SetNContributors(ncombi);
466 //----------------------------------------------------------------------------
467 Int_t AliVertexerTracks::PrepareTracks(TTree &trkTree,Int_t optImpParCut)
470 // Propagate tracks to initial vertex position and store them in a TObjArray
473 Double_t maxd0z0 = fMaxd0z0; // default is 5 mm
475 Double_t sigmaCurr[3];
476 Double_t normdistx,normdisty;
477 Float_t d0z0[2],covd0z0[3];
479 Double_t field=GetField();
481 AliESDVertex *initVertex = new AliESDVertex(fNominalPos,fNominalCov,1,1);
483 Int_t nEntries = (Int_t)trkTree.GetEntries();
484 if(!fTrkArray.IsEmpty()) fTrkArray.Delete();
487 printf(" PrepareTracks()\n");
490 for(Int_t i=0; i<nEntries; i++) {
491 AliESDtrack *track = new AliESDtrack;
492 trkTree.SetBranchAddress("tracks",&track);
495 // propagate track to vertex
496 if(optImpParCut<=1 || fOnlyFitter) { // optImpParCut==1 or 0
497 track->RelateToVertex(initVertex,field,100.);
498 } else { // optImpParCut==2
499 fCurrentVertex->GetSigmaXYZ(sigmaCurr);
500 normdistx = TMath::Abs(fCurrentVertex->GetXv()-fNominalPos[0])/TMath::Sqrt(sigmaCurr[0]*sigmaCurr[0]+fNominalCov[0]);
501 normdisty = TMath::Abs(fCurrentVertex->GetYv()-fNominalPos[1])/TMath::Sqrt(sigmaCurr[1]*sigmaCurr[1]+fNominalCov[2]);
502 if(normdistx < 3. && normdisty < 3. &&
503 (sigmaCurr[0]+sigmaCurr[1])<(TMath::Sqrt(fNominalCov[0])+TMath::Sqrt(fNominalCov[2]))) {
504 track->RelateToVertex(fCurrentVertex,field,100.);
506 track->RelateToVertex(initVertex,field,100.);
510 track->GetImpactParameters(d0z0,covd0z0);
511 sigma = TMath::Sqrt(covd0z0[0]);
512 maxd0rphi = fNSigma*sigma;
513 if(optImpParCut==1) maxd0rphi *= 5.;
517 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);
519 // during iterations 1 and 2, if fConstraint=kFALSE,
520 // select tracks with d0oplusz0 < maxd0z0
521 if(optImpParCut>=1 && !fConstraint && nEntries>=3 &&
522 fVert.GetNContributors()>0) {
523 if(TMath::Sqrt(d0z0[0]*d0z0[0]+d0z0[1]*d0z0[1]) > maxd0z0) {
524 if(fDebug) printf(" rejected\n");
525 delete track; continue;
529 // select tracks with d0rphi < maxd0rphi
530 if(optImpParCut>0 && TMath::Abs(d0z0[0]) > maxd0rphi) {
531 if(fDebug) printf(" rejected\n");
532 delete track; continue;
535 fTrkArray.AddLast(track);
543 //---------------------------------------------------------------------------
544 AliESDVertex* AliVertexerTracks::RemoveTracksFromVertex(AliESDVertex *inVtx,
549 // Removes tracks in trksTree from fit of inVtx
552 if(!strstr(inVtx->GetTitle(),"VertexerTracksWithConstraintOnlyFitter"))
553 printf("WARNING: result of tracks' removal will be only approximately correct\n");
556 rv(0,0) = inVtx->GetXv();
557 rv(1,0) = inVtx->GetYv();
558 rv(2,0) = inVtx->GetZv();
561 inVtx->GetCovMatrix(cov);
563 vV(0,1) = cov[1]; vV(1,0) = cov[1];
565 vV(0,2) = cov[3]; vV(2,0) = cov[3];
566 vV(1,2) = cov[4]; vV(2,1) = cov[4];
569 TMatrixD sumWi(TMatrixD::kInverted,vV);
570 TMatrixD sumWiri(sumWi,TMatrixD::kMult,rv);
572 Int_t nUsedTrks = inVtx->GetNContributors();
573 Double_t chi2 = inVtx->GetChi2();
575 AliESDtrack *track = 0;
576 trksTree->SetBranchAddress("tracks",&track);
577 Int_t ntrks = trksTree->GetEntries();
578 for(Int_t i=0;i<ntrks;i++) {
579 trksTree->GetEvent(i);
580 if(!inVtx->UsesTrack(track->GetID())) {
581 printf("track %d was not used in vertex fit\n",track->GetID());
584 Double_t alpha = track->GetAlpha();
585 Double_t xl = diamondxy[0]*TMath::Cos(alpha)+diamondxy[1]*TMath::Sin(alpha);
586 track->AliExternalTrackParam::PropagateTo(xl,AliTracker::GetBz());
587 // vector of track global coordinates
589 // covariance matrix of ri
592 // get space point from track
593 if(!TrackToPoint(track,ri,wWi)) continue;
595 TMatrixD wWiri(wWi,TMatrixD::kMult,ri);
601 TMatrixD deltar = rv; deltar -= ri;
602 TMatrixD wWideltar(wWi,TMatrixD::kMult,deltar);
603 Double_t chi2i = deltar(0,0)*wWideltar(0,0)+
604 deltar(1,0)*wWideltar(1,0)+
605 deltar(2,0)*wWideltar(2,0);
606 // remove from total chi2
611 printf("Trying to remove too many tracks!\n");
619 // new inverted of weights matrix
620 TMatrixD invsumWi(TMatrixD::kInverted,sumWi);
622 // new position of primary vertex
623 rvnew.Mult(vVnew,sumWiri);
625 Double_t position[3];
626 position[0] = rvnew(0,0);
627 position[1] = rvnew(1,0);
628 position[2] = rvnew(2,0);
636 // store data in the vertex object
637 AliESDVertex *outVtx = new AliESDVertex(position,cov,chi2,nUsedTrks);
638 outVtx->SetTitle(inVtx->GetTitle());
640 inVtx->GetTruePos(tp);
641 outVtx->SetTruePos(tp);
642 UShort_t *inindices = inVtx->GetIndices();
643 UShort_t *outindices = new UShort_t[outVtx->GetNContributors()];
646 for(Int_t k=0; k<inVtx->GetNIndices(); k++) {
648 for(Int_t l=0; l<ntrks; l++) {
649 trksTree->GetEvent(l);
650 if(inindices[k]==track->GetID()) copyindex=kFALSE;
653 outindices[j] = inindices[k]; j++;
656 outVtx->SetIndices(outVtx->GetNContributors(),outindices);
657 delete [] outindices;
660 printf("Vertex before removing tracks:\n");
661 inVtx->PrintStatus();
662 inVtx->PrintIndices();
663 printf("Vertex after removing tracks:\n");
664 outVtx->PrintStatus();
665 outVtx->PrintIndices();
670 //---------------------------------------------------------------------------
671 void AliVertexerTracks::SetSkipTracks(Int_t n,Int_t *skipped)
674 // Mark the tracks not to be used in the vertex reconstruction.
675 // Tracks are identified by AliESDtrack::GetID()
677 fNTrksToSkip = n; fTrksToSkip = new Int_t[n];
678 for(Int_t i=0;i<n;i++) fTrksToSkip[i] = skipped[i];
681 //---------------------------------------------------------------------------
682 void AliVertexerTracks::SetVtxStart(AliESDVertex *vtx)
685 // Set initial vertex knowledge
687 vtx->GetXYZ(fNominalPos);
688 vtx->GetCovMatrix(fNominalCov);
692 //---------------------------------------------------------------------------
693 void AliVertexerTracks::StrLinVertexFinderMinDist(Int_t optUseWeights)
695 // Calculate the point at minimum distance to prepared tracks
700 for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
701 const Int_t knacc = (Int_t)fTrkArray.GetEntriesFast();
702 Double_t field=GetField();
705 Double_t (*vectP0)[3]=new Double_t [knacc][3];
706 Double_t (*vectP1)[3]=new Double_t [knacc][3];
709 Double_t dsum[3]={0,0,0};
710 for(Int_t i=0;i<3;i++)
711 for(Int_t j=0;j<3;j++)sum[i][j]=0;
712 for(Int_t i=0; i<knacc; i++){
713 track1 = (AliESDtrack*)fTrkArray.At(i);
714 Double_t alpha=track1->GetAlpha();
715 Double_t mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
716 Double_t pos[3],dir[3];
717 track1->GetXYZAt(mindist,field,pos);
718 track1->GetPxPyPzAt(mindist,field,dir);
719 AliStrLine *line1 = new AliStrLine(pos,dir);
720 // AliStrLine *line1 = new AliStrLine();
721 // track1->ApproximateHelixWithLine(mindist,field,line1);
723 Double_t p0[3],cd[3];
726 Double_t p1[3]={p0[0]+cd[0],p0[1]+cd[1],p0[2]+cd[2]};
736 if(optUseWeights==0)GetStrLinDerivMatrix(p0,p1,matr,dknow);
737 if(optUseWeights==1){
739 sigmasq[0]=track1->GetSigmaY2();
740 sigmasq[1]=track1->GetSigmaY2();
741 sigmasq[2]=track1->GetSigmaZ2();
742 GetStrLinDerivMatrix(p0,p1,sigmasq,matr,dknow);
745 for(Int_t iii=0;iii<3;iii++){
746 dsum[iii]+=dknow[iii];
747 for(Int_t lj=0;lj<3;lj++) sum[iii][lj]+=matr[iii][lj];
753 Double_t det=GetDeterminant3X3(sum);
756 for(Int_t zz=0;zz<3;zz++){
757 for(Int_t ww=0;ww<3;ww++){
758 for(Int_t kk=0;kk<3;kk++) vett[ww][kk]=sum[ww][kk];
760 for(Int_t kk=0;kk<3;kk++) vett[kk][zz]=dsum[kk];
761 initPos[zz]=GetDeterminant3X3(vett)/det;
765 for(Int_t i=0; i<knacc; i++){
766 Double_t p0[3]={0,0,0},p1[3]={0,0,0};
767 for(Int_t ii=0;ii<3;ii++){
768 p0[ii]=vectP0[i][ii];
769 p1[ii]=vectP1[i][ii];
771 sigma+=GetStrLinMinDist(p0,p1,initPos);
774 sigma=TMath::Sqrt(sigma);
776 Warning("StrLinVertexFinderMinDist","Finder did not succed");
781 fVert.SetXYZ(initPos);
782 fVert.SetDispersion(sigma);
783 fVert.SetNContributors(knacc);
785 //---------------------------------------------------------------------------
786 Bool_t AliVertexerTracks::TrackToPoint(AliESDtrack *t,
787 TMatrixD &ri,TMatrixD &wWi) const
790 // Extract from the AliESDtrack the global coordinates ri and covariance matrix
791 // wWi of the space point that it represents (to be used in VertexFitter())
795 Double_t rotAngle = t->GetAlpha();
796 if(rotAngle<0.) rotAngle += 2.*TMath::Pi();
797 Double_t cosRot = TMath::Cos(rotAngle);
798 Double_t sinRot = TMath::Sin(rotAngle);
800 ri(0,0) = t->GetX()*cosRot-t->GetY()*sinRot;
801 ri(1,0) = t->GetX()*sinRot+t->GetY()*cosRot;
804 // matrix to go from global (x,y,z) to local (y,z);
813 // covariance matrix of local (y,z) - inverted
816 t->GetExternalCovariance(cc);
821 if(uUi.Determinant() <= 0.) return kFALSE;
822 TMatrixD uUiInv(TMatrixD::kInverted,uUi);
824 // weights matrix: wWi = qQiT * uUiInv * qQi
825 TMatrixD uUiInvQi(uUiInv,TMatrixD::kMult,qQi);
826 TMatrixD m(qQi,TMatrixD::kTransposeMult,uUiInvQi);
831 //---------------------------------------------------------------------------
832 void AliVertexerTracks::TooFewTracks(const AliESD* esdEvent)
835 // When the number of tracks is < fMinTracks
838 // deal with vertices not found
839 Double_t pos[3],err[3];
841 pos[0] = fNominalPos[0];
842 err[0] = TMath::Sqrt(fNominalCov[0]);
843 pos[1] = fNominalPos[1];
844 err[1] = TMath::Sqrt(fNominalCov[2]);
845 pos[2] = esdEvent->GetVertex()->GetZv();
846 err[2] = esdEvent->GetVertex()->GetZRes();
847 if(err[0]>1. && esdEvent->GetVertex()->GetNContributors()<=0)
848 ncontr = -1; // (x,y,z) = (0,0,0)
849 if(err[0]>1. && esdEvent->GetVertex()->GetNContributors()>0)
850 ncontr = -2; // (x,y,z) = (0,0,z_fromSPD)
851 if(err[0]<1. && esdEvent->GetVertex()->GetNContributors()<=0)
852 ncontr = -3; // (x,y,z) = (x_mean,y_mean,0)
853 if(err[0]<1. && esdEvent->GetVertex()->GetNContributors()>0)
854 ncontr = -4; // (x,y,z) = (x_mean,y_mean,z_fromSPD)
856 fCurrentVertex = new AliESDVertex(pos,err);
857 fCurrentVertex->SetNContributors(ncontr);
860 esdEvent->GetVertex()->GetTruePos(tp);
861 fCurrentVertex->SetTruePos(tp);
863 fCurrentVertex->SetTitle("VertexerTracksWithConstraint");
865 fCurrentVertex->SetTitle("VertexerTracksNoConstraint");
870 //---------------------------------------------------------------------------
871 void AliVertexerTracks::VertexFinder(Int_t optUseWeights)
874 // Get estimate of vertex position in (x,y) from tracks DCA
878 for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
879 Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
880 Double_t aver[3]={0.,0.,0.};
881 Double_t aversq[3]={0.,0.,0.};
882 Double_t sigmasq[3]={0.,0.,0.};
887 Double_t pos[3],dir[3];
888 Double_t alpha,mindist;
889 Double_t field=GetField();
891 for(Int_t i=0; i<nacc; i++){
892 track1 = (AliESDtrack*)fTrkArray.At(i);
893 alpha=track1->GetAlpha();
894 mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
895 track1->GetXYZAt(mindist,field,pos);
896 track1->GetPxPyPzAt(mindist,field,dir);
897 AliStrLine *line1 = new AliStrLine(pos,dir);
899 // AliStrLine *line1 = new AliStrLine();
900 // track1->ApproximateHelixWithLine(mindist,field,line1);
902 for(Int_t j=i+1; j<nacc; j++){
903 track2 = (AliESDtrack*)fTrkArray.At(j);
904 alpha=track2->GetAlpha();
905 mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
906 track2->GetXYZAt(mindist,field,pos);
907 track2->GetPxPyPzAt(mindist,field,dir);
908 AliStrLine *line2 = new AliStrLine(pos,dir);
909 // AliStrLine *line2 = new AliStrLine();
910 // track2->ApproximateHelixWithLine(mindist,field,line2);
911 Double_t distCA=line2->GetDCA(line1);
912 //printf("%d %d %f\n",i,j,distCA);
913 if(fDCAcut<=0 || (fDCAcut>0&&distCA<fDCAcut)){
914 Double_t pnt1[3],pnt2[3],crosspoint[3];
916 if(optUseWeights<=0){
917 Int_t retcode = line2->Cross(line1,crosspoint);
920 for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
921 for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]);
925 Int_t retcode = line1->CrossPoints(line2,pnt1,pnt2);
928 alpha=track1->GetAlpha();
929 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
930 Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2();
931 alpha=track2->GetAlpha();
932 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
933 Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2();
934 Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2();
935 Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
936 Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
937 Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
938 crosspoint[0]=wx1*pnt1[0] + wx2*pnt2[0];
939 crosspoint[1]=wy1*pnt1[1] + wy2*pnt2[1];
940 crosspoint[2]=wz1*pnt1[2] + wz2*pnt2[2];
943 for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
944 for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]);
953 for(Int_t jj=0;jj<3;jj++){
954 initPos[jj] = aver[jj]/ncombi;
955 //printf("%f\n",initPos[jj]);
957 sigmasq[jj]=aversq[jj]-initPos[jj]*initPos[jj];
960 sigma=TMath::Sqrt(TMath::Abs(sigma));
963 Warning("VertexFinder","Finder did not succed");
966 fVert.SetXYZ(initPos);
967 fVert.SetDispersion(sigma);
968 fVert.SetNContributors(ncombi);
970 //---------------------------------------------------------------------------
971 void AliVertexerTracks::VertexFitter(Bool_t useConstraint)
974 // The optimal estimate of the vertex position is given by a "weighted
975 // average of tracks positions"
976 // Original method: V. Karimaki, CMS Note 97/0051
979 fVert.GetXYZ(initPos);
980 Int_t arrEntries = (Int_t)fTrkArray.GetEntries();
981 if(arrEntries==1) useConstraint=kTRUE;
983 printf(" VertexFitter(): start\n");
984 printf(" Number of tracks in array: %d\n",(Int_t)fTrkArray.GetEntriesFast());
985 printf(" Minimum # tracks required in fit: %d\n",fMinTracks);
986 printf("Vertex position after finder: %f,%f,%f\n",initPos[0],initPos[1],initPos[2]);
987 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]));
993 rv(0,0) = initPos[0];
994 rv(1,0) = initPos[1];
996 Double_t xlStart,alpha;
998 Double_t chi2,chi2i,chi2b;
1002 // initial vertex covariance matrix
1004 vVb(0,0) = fNominalCov[0];
1005 vVb(0,1) = fNominalCov[1];
1007 vVb(1,0) = fNominalCov[1];
1008 vVb(1,1) = fNominalCov[2];
1012 vVb(2,2) = fNominalCov[5];
1013 TMatrixD vVbInv(TMatrixD::kInverted,vVb);
1015 rb(0,0) = fNominalPos[0];
1016 rb(1,0) = fNominalPos[1];
1017 rb(2,0) = fNominalPos[2];
1018 TMatrixD vVbInvrb(vVbInv,TMatrixD::kMult,rb);
1022 // 1st - estimate of vtx using all tracks
1023 // 2nd - estimate of global chi2
1024 for(step=0; step<2; step++) {
1025 if(fDebug) printf(" step = %d\n",step);
1029 TMatrixD sumWiri(3,1);
1030 TMatrixD sumWi(3,3);
1031 for(i=0; i<3; i++) {
1033 for(j=0; j<3; j++) sumWi(j,i) = 0.;
1036 // mean vertex constraint
1039 sumWiri(i,0) += vVbInvrb(i,0);
1040 for(k=0;k<3;k++) sumWi(i,k) += vVbInv(i,k);
1043 TMatrixD deltar = rv; deltar -= rb;
1044 TMatrixD vVbInvdeltar(vVbInv,TMatrixD::kMult,deltar);
1045 chi2b = deltar(0,0)*vVbInvdeltar(0,0)+
1046 deltar(1,0)*vVbInvdeltar(1,0)+
1047 deltar(2,0)*vVbInvdeltar(2,0);
1053 for(k=0; k<arrEntries; k++) {
1054 // get track from track array
1055 t = (AliESDtrack*)fTrkArray.At(k);
1056 alpha = t->GetAlpha();
1057 xlStart = initPos[0]*TMath::Cos(alpha)+initPos[1]*TMath::Sin(alpha);
1058 // to vtxSeed (from finder)
1059 t->AliExternalTrackParam::PropagateTo(xlStart,AliTracker::GetBz());
1062 // vector of track global coordinates
1064 // covariance matrix of ri
1067 // get space point from track
1068 if(!TrackToPoint(t,ri,wWi)) continue;
1070 TMatrixD wWiri(wWi,TMatrixD::kMult,ri);
1073 TMatrixD deltar = rv; deltar -= ri;
1074 TMatrixD wWideltar(wWi,TMatrixD::kMult,deltar);
1075 chi2i = deltar(0,0)*wWideltar(0,0)+
1076 deltar(1,0)*wWideltar(1,0)+
1077 deltar(2,0)*wWideltar(2,0);
1079 // add to total chi2
1086 } // end loop on tracks
1088 if(nUsedTrks < fMinTracks) {
1093 Double_t determinant = sumWi.Determinant();
1094 //cerr<<" determinant: "<<determinant<<endl;
1095 if(determinant < 100.) {
1096 printf("det(V) = 0\n");
1101 // inverted of weights matrix
1102 TMatrixD invsumWi(TMatrixD::kInverted,sumWi);
1105 // position of primary vertex
1106 rv.Mult(vV,sumWiri);
1108 } // end loop on the 2 steps
1113 if(fDebug) printf("TooFewTracks\n");
1114 fCurrentVertex = new AliESDVertex(0.,0.,-1);
1118 Double_t position[3];
1119 position[0] = rv(0,0);
1120 position[1] = rv(1,0);
1121 position[2] = rv(2,0);
1122 Double_t covmatrix[6];
1123 covmatrix[0] = vV(0,0);
1124 covmatrix[1] = vV(0,1);
1125 covmatrix[2] = vV(1,1);
1126 covmatrix[3] = vV(0,2);
1127 covmatrix[4] = vV(1,2);
1128 covmatrix[5] = vV(2,2);
1130 // store data in the vertex object
1131 fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nUsedTrks);
1134 printf(" VertexFitter(): finish\n");
1135 printf(" rv = ( %f , %f , %f )\n\n",rv(0,0),rv(1,0),rv(2,0));
1136 fCurrentVertex->PrintStatus();
1141 //----------------------------------------------------------------------------
1142 AliESDVertex* AliVertexerTracks::VertexForSelectedTracks(TTree *trkTree,Bool_t optUseFitter,Bool_t optPropagate)
1145 // Return vertex from tracks in trkTree
1150 // get tracks and propagate them to initial vertex position
1151 Int_t nTrksPrep = PrepareTracks(*trkTree,0);
1152 if(nTrksPrep < TMath::Max(2,fMinTracks) ) {
1153 if(fDebug) printf("TooFewTracks\n");
1154 fCurrentVertex = new AliESDVertex(0.,0.,-1);
1155 return fCurrentVertex;
1159 case 1: StrLinVertexFinderMinDist(1); break;
1160 case 2: StrLinVertexFinderMinDist(0); break;
1161 case 3: HelixVertexFinder(); break;
1162 case 4: VertexFinder(1); break;
1163 case 5: VertexFinder(0); break;
1164 default: printf("Wrong algorithm\n"); break;
1167 if(fDebug) printf(" Vertex finding completed\n");
1171 VertexFitter(fConstraint);
1172 if(fDebug) printf(" Vertex fit completed\n");
1174 Double_t position[3]={fVert.GetXv(),fVert.GetYv(),fVert.GetZv()};
1175 Double_t covmatrix[6]={0.,0.,0.,0.,0.,0.};
1176 Double_t chi2=99999.;
1177 Int_t nUsedTrks=fVert.GetNContributors();
1178 fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nUsedTrks);
1180 fCurrentVertex->SetDispersion(fVert.GetDispersion());
1183 // set indices of used tracks and propagate track to found vertex
1184 UShort_t *indices = 0;
1185 AliESDtrack *eta = 0;
1186 if(fCurrentVertex->GetNContributors()>0) {
1187 indices = new UShort_t[fCurrentVertex->GetNContributors()];
1188 for(Int_t jj=0;jj<(Int_t)fTrkArray.GetEntriesFast();jj++) {
1189 eta = (AliESDtrack*)fTrkArray.At(jj);
1190 indices[jj] = (UShort_t)eta->GetID();
1191 if(optPropagate&&optUseFitter){
1192 if(TMath::Sqrt(fCurrentVertex->GetXv()*fCurrentVertex->GetXv()+fCurrentVertex->GetYv()*fCurrentVertex->GetYv())<3.) {
1193 eta->RelateToVertex(fCurrentVertex,GetField(),100.);
1194 if(fDebug) printf("Track %d propagated to found vertex\n",jj);
1196 AliWarning("Found vertex outside beam pipe!");
1200 fCurrentVertex->SetIndices(fCurrentVertex->GetNContributors(),indices);
1205 return fCurrentVertex;
1207 //----------------------------------------------------------------------------
1208 AliESDVertex* AliVertexerTracks::VertexForSelectedTracks(TObjArray *trkArray,Bool_t optUseFitter, Bool_t optPropagate)
1211 // Return vertex from array of tracks
1214 // get tracks and propagate them to initial vertex position
1215 Int_t nTrks = trkArray->GetEntriesFast();
1216 if(nTrks < TMath::Max(2,fMinTracks) ) {
1217 if(fDebug) printf("TooFewTracks\n");
1218 fCurrentVertex = new AliESDVertex(0.,0.,-1);
1219 return fCurrentVertex;
1221 TTree *trkTree = new TTree("TreeT","tracks");
1222 AliESDtrack *esdTrack = 0;
1223 trkTree->Branch("tracks","AliESDtrack",&esdTrack);
1224 for(Int_t i=0; i<nTrks; i++){
1225 esdTrack = (AliESDtrack*)trkArray->At(i);
1229 AliESDVertex *vtx = VertexForSelectedTracks(trkTree,optUseFitter,optPropagate);
1233 //--------------------------------------------------------------------------