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"
34 #include "AliESDEvent.h"
35 #include "AliESDtrack.h"
37 ClassImp(AliVertexerTracks)
40 //----------------------------------------------------------------------------
41 AliVertexerTracks::AliVertexerTracks():
63 // Default constructor
72 //----------------------------------------------------------------------------
73 AliVertexerTracks::AliVertexerTracks(Double_t fieldkG):
95 // Standard constructor
105 //-----------------------------------------------------------------------------
106 AliVertexerTracks::~AliVertexerTracks()
108 // Default Destructor
109 // The objects pointed by the following pointer are not owned
110 // by this class and are not deleted
112 if(fTrksToSkip) { delete [] fTrksToSkip; fTrksToSkip=NULL; }
114 //----------------------------------------------------------------------------
115 AliESDVertex* AliVertexerTracks::FindPrimaryVertex(const AliESDEvent *esdEvent)
118 // Primary vertex for current ESD event
120 // 1st with 5*fNSigma*sigma cut w.r.t. to initial vertex
121 // + cut on sqrt(d0d0+z0z0) if fConstraint=kFALSE
122 // 2nd with fNSigma*sigma cut w.r.t. to vertex found in 1st iteration)
123 // All ESD tracks with inside the beam pipe are then propagated to found vertex
127 // accept 1-track case only if constraint is available
128 if(!fConstraint && fMinTracks==1) fMinTracks=2;
130 // read tracks from ESD
131 Int_t nTrksTot = (Int_t)esdEvent->GetNumberOfTracks();
133 if(fDebug) printf("TooFewTracks\n");
134 TooFewTracks(esdEvent);
135 return fCurrentVertex;
138 TTree *trkTree = new TTree("TreeT","tracks");
139 AliESDtrack *esdTrack = 0;
140 trkTree->Branch("tracks","AliESDtrack",&esdTrack);
143 for(Int_t i=0; i<nTrksTot; i++) {
144 AliESDtrack *et = esdEvent->GetTrack(i);
145 esdTrack = new AliESDtrack(*et);
146 // check tracks to skip
148 for(Int_t j=0; j<fNTrksToSkip; j++) {
149 if(et->GetID()==fTrksToSkip[j]) {
150 if(fDebug) printf("skipping track: %d\n",i);
154 if(skipThis) {delete esdTrack;continue;}
156 if(!(esdTrack->GetStatus()&AliESDtrack::kITSin)) {delete esdTrack;continue;}
157 if(fITSrefit && !(esdTrack->GetStatus()&AliESDtrack::kITSrefit)) {delete esdTrack;continue;}
158 Int_t nclus=esdTrack->GetNcls(0); // check number of clusters in ITS
159 if(nclus<fMinITSClusters) {delete esdTrack;continue;}
165 // If fConstraint=kFALSE
166 // run VertexFinder(1) to get rough estimate of initVertex (x,y)
168 // fill fTrkArray, for VertexFinder()
169 if(!fTrkArray.IsEmpty()) fTrkArray.Delete();
170 PrepareTracks(*trkTree,0);
171 Double_t cutsave = fDCAcut; fDCAcut = 0.1; // 1 mm
172 VertexFinder(1); // using weights, cutting dca < fDCAcut
175 if(fVert.GetNContributors()>0) {
176 fVert.GetXYZ(fNominalPos);
177 fNominalPos[0] = fVert.GetXv();
178 fNominalPos[1] = fVert.GetYv();
179 fNominalPos[2] = fVert.GetZv();
180 if(fDebug) printf("No mean vertex: VertexFinder gives (%f, %f, %f)\n",fNominalPos[0],fNominalPos[1],fNominalPos[2]);
185 if(fDebug) printf("No mean vertex and VertexFinder failed\n");
192 // propagate tracks to fNominalPos vertex
194 // if(constraint) reject for |d0|>5*fNSigma*sigma w.r.t. fNominal... vertex
195 // else reject for |d0|\oplus|z0| > 5 mm w.r.t. fNominal... vertex
197 // propagate tracks to best between initVertex and fCurrentVertex
198 // preselect tracks (reject for |d0|>fNSigma*sigma w.r.t. best
199 // between initVertex and fCurrentVertex)
200 for(Int_t iter=0; iter<2; iter++) {
201 if(fOnlyFitter && iter==0) continue;
202 Int_t nTrksPrep = PrepareTracks(*trkTree,iter+1);
203 if(fDebug) printf(" tracks prepared - iteration %d: %d\n",iter+1,nTrksPrep);
204 if(nTrksPrep < fMinTracks) {
205 if(fDebug) printf("TooFewTracks\n");
206 TooFewTracks(esdEvent);
207 if(fDebug) fCurrentVertex->PrintStatus();
210 return fCurrentVertex;
216 if(fDebug) printf("Just one track\n");
217 OneTrackVertFinder();
220 case 1: StrLinVertexFinderMinDist(1); break;
221 case 2: StrLinVertexFinderMinDist(0); break;
222 case 3: HelixVertexFinder(); break;
223 case 4: VertexFinder(1); break;
224 case 5: VertexFinder(0); break;
225 default: printf("Wrong algorithm\n"); break;
228 if(fDebug) printf(" Vertex finding completed\n");
232 VertexFitter(fConstraint);
233 if(fDebug) printf(" Vertex fit completed\n");
234 if(iter==0) fTrkArray.Delete();
235 } // end loop on the two iterations
240 fCurrentVertex->SetTitle("VertexerTracksWithConstraintOnlyFitter");
242 fCurrentVertex->SetTitle("VertexerTracksWithConstraint");
245 fCurrentVertex->SetTitle("VertexerTracksNoConstraint");
249 // set indices of used tracks
250 UShort_t *indices = 0;
251 AliESDtrack *ett = 0;
252 if(fCurrentVertex->GetNContributors()>0) {
253 indices = new UShort_t[fCurrentVertex->GetNContributors()];
254 for(Int_t jj=0;jj<(Int_t)fTrkArray.GetEntriesFast();jj++) {
255 ett = (AliESDtrack*)fTrkArray.At(jj);
256 indices[jj] = (UShort_t)ett->GetID();
258 fCurrentVertex->SetIndices(fCurrentVertex->GetNContributors(),indices);
266 if(fTrksToSkip) { delete [] fTrksToSkip; fTrksToSkip=NULL; }
269 if(fDebug) fCurrentVertex->PrintStatus();
270 if(fDebug) fCurrentVertex->PrintIndices();
273 return fCurrentVertex;
275 //------------------------------------------------------------------------
276 Double_t AliVertexerTracks::GetDeterminant3X3(Double_t matr[][3])
279 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];
282 //-------------------------------------------------------------------------
283 void AliVertexerTracks::GetStrLinDerivMatrix(Double_t *p0,Double_t *p1,Double_t (*m)[3],Double_t *d)
286 Double_t x12=p0[0]-p1[0];
287 Double_t y12=p0[1]-p1[1];
288 Double_t z12=p0[2]-p1[2];
289 Double_t kk=x12*x12+y12*y12+z12*z12;
290 m[0][0]=2-2/kk*x12*x12;
291 m[0][1]=-2/kk*x12*y12;
292 m[0][2]=-2/kk*x12*z12;
293 m[1][0]=-2/kk*x12*y12;
294 m[1][1]=2-2/kk*y12*y12;
295 m[1][2]=-2/kk*y12*z12;
296 m[2][0]=-2/kk*x12*z12;
298 m[2][2]=2-2/kk*z12*z12;
299 d[0]=2*p0[0]-2/kk*p0[0]*x12*x12-2/kk*p0[2]*x12*z12-2/kk*p0[1]*x12*y12;
300 d[1]=2*p0[1]-2/kk*p0[1]*y12*y12-2/kk*p0[0]*x12*y12-2/kk*p0[2]*z12*y12;
301 d[2]=2*p0[2]-2/kk*p0[2]*z12*z12-2/kk*p0[0]*x12*z12-2/kk*p0[1]*z12*y12;
304 //--------------------------------------------------------------------------
305 void AliVertexerTracks::GetStrLinDerivMatrix(Double_t *p0,Double_t *p1,Double_t *sigmasq,Double_t (*m)[3],Double_t *d)
308 Double_t x12=p1[0]-p0[0];
309 Double_t y12=p1[1]-p0[1];
310 Double_t z12=p1[2]-p0[2];
312 Double_t den= x12*x12*sigmasq[1]*sigmasq[2]+y12*y12*sigmasq[0]*sigmasq[2]+z12*z12*sigmasq[0]*sigmasq[1];
314 Double_t kk= 2*(x12*x12/sigmasq[0]+y12*y12/sigmasq[1]+z12*z12/sigmasq[2]);
317 cc[0]=-x12/sigmasq[0];
318 cc[1]=-y12/sigmasq[1];
319 cc[2]=-z12/sigmasq[2];
321 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;
323 Double_t ss= -p0[0]*cc[0]-p0[1]*cc[1]-p0[2]*cc[2];
326 aa[0]=x12*sigmasq[1]*sigmasq[2]/den;
327 aa[1]=y12*sigmasq[0]*sigmasq[2]/den;
328 aa[2]=z12*sigmasq[0]*sigmasq[1]/den;
330 m[0][0]=aa[0]*(aa[0]*kk+2*cc[0])+2*cc[0]*aa[0]+2/sigmasq[0];
331 m[0][1]=aa[1]*(aa[0]*kk+2*cc[0])+2*cc[1]*aa[0];
332 m[0][2]=aa[2]*(aa[0]*kk+2*cc[0])+2*cc[2]*aa[0];
334 m[1][0]=aa[0]*(aa[1]*kk+2*cc[1])+2*cc[0]*aa[1];
335 m[1][1]=aa[1]*(aa[1]*kk+2*cc[1])+2*cc[1]*aa[1]+2/sigmasq[1];
336 m[1][2]=aa[2]*(aa[1]*kk+2*cc[1])+2*cc[2]*aa[1];
338 m[2][0]=aa[0]*(aa[2]*kk+2*cc[2])+2*cc[0]*aa[2];
339 m[2][1]=aa[1]*(aa[2]*kk+2*cc[2])+2*cc[1]*aa[2];
340 m[2][2]=aa[2]*(aa[2]*kk+2*cc[2])+2*cc[2]*aa[2]+2/sigmasq[2];
342 d[0]=-ww*(aa[0]*kk+2*cc[0])-2*ss*aa[0]+2*p0[0]/sigmasq[0];
343 d[1]=-ww*(aa[1]*kk+2*cc[1])-2*ss*aa[1]+2*p0[1]/sigmasq[1];
344 d[2]=-ww*(aa[2]*kk+2*cc[2])-2*ss*aa[2]+2*p0[2]/sigmasq[2];
347 //--------------------------------------------------------------------------
348 Double_t AliVertexerTracks::GetStrLinMinDist(Double_t *p0,Double_t *p1,Double_t *x0)
351 Double_t x12=p0[0]-p1[0];
352 Double_t y12=p0[1]-p1[1];
353 Double_t z12=p0[2]-p1[2];
354 Double_t x10=p0[0]-x0[0];
355 Double_t y10=p0[1]-x0[1];
356 Double_t z10=p0[2]-x0[2];
357 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);
359 //---------------------------------------------------------------------------
360 void AliVertexerTracks::OneTrackVertFinder()
362 // find vertex for events with 1 track, using DCA to nominal beam axis
363 if(fDebug) printf("Number of prepared tracks =%d - Call OneTrackVertFinder",fTrkArray.GetEntries());
365 track1 = (AliESDtrack*)fTrkArray.At(0);
366 Double_t field=GetFieldkG();
367 Double_t alpha=track1->GetAlpha();
368 Double_t mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
369 Double_t pos[3],dir[3];
370 track1->GetXYZAt(mindist,field,pos);
371 track1->GetPxPyPzAt(mindist,field,dir);
372 AliStrLine *line1 = new AliStrLine(pos,dir);
373 Double_t p1[3]={fNominalPos[0],fNominalPos[1],0.};
374 Double_t p2[3]={fNominalPos[0],fNominalPos[1],10.};
375 AliStrLine *zeta=new AliStrLine(p1,p2,kTRUE);
376 Double_t crosspoint[3]={0.,0.,0.};
379 Int_t retcode = zeta->Cross(line1,crosspoint);
381 sigma=line1->GetDistFromPoint(crosspoint);
386 fVert.SetXYZ(crosspoint);
387 fVert.SetDispersion(sigma);
388 fVert.SetNContributors(nContrib);
390 //---------------------------------------------------------------------------
391 void AliVertexerTracks::HelixVertexFinder()
393 // Get estimate of vertex position in (x,y) from tracks DCA
398 for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
399 Double_t field=GetFieldkG();
401 Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
403 Double_t aver[3]={0.,0.,0.};
404 Double_t averquad[3]={0.,0.,0.};
405 Double_t sigmaquad[3]={0.,0.,0.};
412 Double_t alpha, cs, sn;
413 Double_t crosspoint[3];
414 for(Int_t i=0; i<nacc; i++){
415 track1 = (AliESDtrack*)fTrkArray.At(i);
418 for(Int_t j=i+1; j<nacc; j++){
419 track2 = (AliESDtrack*)fTrkArray.At(j);
421 distCA=track2->PropagateToDCA(track1,field);
422 if(fDCAcut<=0 ||(fDCAcut>0&&distCA<fDCAcut)){
423 track1->GetExternalParameters(x,par);
424 alpha=track1->GetAlpha();
425 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
426 Double_t x1=x*cs - par[0]*sn;
427 Double_t y1=x*sn + par[0]*cs;
429 Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2();
430 track2->GetExternalParameters(x,par);
431 alpha=track2->GetAlpha();
432 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
433 Double_t x2=x*cs - par[0]*sn;
434 Double_t y2=x*sn + par[0]*cs;
436 Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2();
437 Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2();
438 Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
439 Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
440 Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
441 crosspoint[0]=wx1*x1 + wx2*x2;
442 crosspoint[1]=wy1*y1 + wy2*y2;
443 crosspoint[2]=wz1*z1 + wz2*z2;
446 for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
447 for(Int_t jj=0;jj<3;jj++)averquad[jj]+=(crosspoint[jj]*crosspoint[jj]);
453 for(Int_t jj=0;jj<3;jj++){
454 initPos[jj] = aver[jj]/ncombi;
455 averquad[jj]/=ncombi;
456 sigmaquad[jj]=averquad[jj]-initPos[jj]*initPos[jj];
457 sigma+=sigmaquad[jj];
459 sigma=TMath::Sqrt(TMath::Abs(sigma));
462 Warning("HelixVertexFinder","Finder did not succed");
465 fVert.SetXYZ(initPos);
466 fVert.SetDispersion(sigma);
467 fVert.SetNContributors(ncombi);
469 //----------------------------------------------------------------------------
470 Int_t AliVertexerTracks::PrepareTracks(TTree &trkTree,Int_t optImpParCut)
473 // Propagate tracks to initial vertex position and store them in a TObjArray
476 Double_t maxd0z0 = fMaxd0z0; // default is 5 mm
478 Double_t sigmaCurr[3];
479 Double_t normdistx,normdisty;
480 Float_t d0z0[2],covd0z0[3];
482 Double_t field=GetFieldkG();
484 AliESDVertex *initVertex = new AliESDVertex(fNominalPos,fNominalCov,1,1);
486 Int_t nEntries = (Int_t)trkTree.GetEntries();
487 if(!fTrkArray.IsEmpty()) fTrkArray.Delete();
490 printf(" PrepareTracks()\n");
493 for(Int_t i=0; i<nEntries; i++) {
494 AliESDtrack *track = new AliESDtrack;
495 trkTree.SetBranchAddress("tracks",&track);
498 // propagate track to vertex
499 if(optImpParCut<=1 || fOnlyFitter) { // optImpParCut==1 or 0
500 track->RelateToVertex(initVertex,field,100.);
501 } else { // optImpParCut==2
502 fCurrentVertex->GetSigmaXYZ(sigmaCurr);
503 normdistx = TMath::Abs(fCurrentVertex->GetXv()-fNominalPos[0])/TMath::Sqrt(sigmaCurr[0]*sigmaCurr[0]+fNominalCov[0]);
504 normdisty = TMath::Abs(fCurrentVertex->GetYv()-fNominalPos[1])/TMath::Sqrt(sigmaCurr[1]*sigmaCurr[1]+fNominalCov[2]);
505 if(normdistx < 3. && normdisty < 3. &&
506 (sigmaCurr[0]+sigmaCurr[1])<(TMath::Sqrt(fNominalCov[0])+TMath::Sqrt(fNominalCov[2]))) {
507 track->RelateToVertex(fCurrentVertex,field,100.);
509 track->RelateToVertex(initVertex,field,100.);
513 track->GetImpactParameters(d0z0,covd0z0);
514 sigma = TMath::Sqrt(covd0z0[0]);
515 maxd0rphi = fNSigma*sigma;
516 if(optImpParCut==1) maxd0rphi *= 5.;
520 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);
522 // during iterations 1 and 2, if fConstraint=kFALSE,
523 // select tracks with d0oplusz0 < maxd0z0
524 if(optImpParCut>=1 && !fConstraint && nEntries>=3 &&
525 fVert.GetNContributors()>0) {
526 if(TMath::Sqrt(d0z0[0]*d0z0[0]+d0z0[1]*d0z0[1]) > maxd0z0) {
527 if(fDebug) printf(" rejected\n");
528 delete track; continue;
532 // select tracks with d0rphi < maxd0rphi
533 if(optImpParCut>0 && TMath::Abs(d0z0[0]) > maxd0rphi) {
534 if(fDebug) printf(" rejected\n");
535 delete track; continue;
538 fTrkArray.AddLast(track);
546 //---------------------------------------------------------------------------
547 AliESDVertex* AliVertexerTracks::RemoveTracksFromVertex(AliESDVertex *inVtx,
552 // Removes tracks in trksTree from fit of inVtx
555 if(!strstr(inVtx->GetTitle(),"VertexerTracksWithConstraint")) {
556 printf("ERROR: primary vertex has no constraint: cannot remove tracks\n");
559 if(!strstr(inVtx->GetTitle(),"VertexerTracksWithConstraintOnlyFitter"))
560 printf("WARNING: result of tracks' removal will be only approximately correct\n");
563 rv(0,0) = inVtx->GetXv();
564 rv(1,0) = inVtx->GetYv();
565 rv(2,0) = inVtx->GetZv();
568 inVtx->GetCovMatrix(cov);
570 vV(0,1) = cov[1]; vV(1,0) = cov[1];
572 vV(0,2) = cov[3]; vV(2,0) = cov[3];
573 vV(1,2) = cov[4]; vV(2,1) = cov[4];
576 TMatrixD sumWi(TMatrixD::kInverted,vV);
577 TMatrixD sumWiri(sumWi,TMatrixD::kMult,rv);
579 Int_t nUsedTrks = inVtx->GetNContributors();
580 Double_t chi2 = inVtx->GetChi2();
582 AliESDtrack *track = 0;
583 trksTree->SetBranchAddress("tracks",&track);
584 Int_t ntrks = trksTree->GetEntries();
585 for(Int_t i=0;i<ntrks;i++) {
586 trksTree->GetEvent(i);
587 if(!inVtx->UsesTrack(track->GetID())) {
588 printf("track %d was not used in vertex fit\n",track->GetID());
591 Double_t alpha = track->GetAlpha();
592 Double_t xl = diamondxy[0]*TMath::Cos(alpha)+diamondxy[1]*TMath::Sin(alpha);
593 track->AliExternalTrackParam::PropagateTo(xl,GetFieldkG());
594 // vector of track global coordinates
596 // covariance matrix of ri
599 // get space point from track
600 if(!TrackToPoint(track,ri,wWi)) continue;
602 TMatrixD wWiri(wWi,TMatrixD::kMult,ri);
608 TMatrixD deltar = rv; deltar -= ri;
609 TMatrixD wWideltar(wWi,TMatrixD::kMult,deltar);
610 Double_t chi2i = deltar(0,0)*wWideltar(0,0)+
611 deltar(1,0)*wWideltar(1,0)+
612 deltar(2,0)*wWideltar(2,0);
613 // remove from total chi2
618 printf("Trying to remove too many tracks!\n");
626 // new inverted of weights matrix
627 TMatrixD invsumWi(TMatrixD::kInverted,sumWi);
629 // new position of primary vertex
630 rvnew.Mult(vVnew,sumWiri);
632 Double_t position[3];
633 position[0] = rvnew(0,0);
634 position[1] = rvnew(1,0);
635 position[2] = rvnew(2,0);
643 // store data in the vertex object
644 AliESDVertex *outVtx = new AliESDVertex(position,cov,chi2,nUsedTrks);
645 outVtx->SetTitle(inVtx->GetTitle());
646 UShort_t *inindices = inVtx->GetIndices();
647 UShort_t *outindices = new UShort_t[outVtx->GetNContributors()];
650 for(Int_t k=0; k<inVtx->GetNIndices(); k++) {
652 for(Int_t l=0; l<ntrks; l++) {
653 trksTree->GetEvent(l);
654 if(inindices[k]==track->GetID()) copyindex=kFALSE;
657 outindices[j] = inindices[k]; j++;
660 outVtx->SetIndices(outVtx->GetNContributors(),outindices);
661 delete [] outindices;
664 printf("Vertex before removing tracks:\n");
665 inVtx->PrintStatus();
666 inVtx->PrintIndices();
667 printf("Vertex after removing tracks:\n");
668 outVtx->PrintStatus();
669 outVtx->PrintIndices();
674 //---------------------------------------------------------------------------
675 void AliVertexerTracks::SetSkipTracks(Int_t n,Int_t *skipped)
678 // Mark the tracks not to be used in the vertex reconstruction.
679 // Tracks are identified by AliESDtrack::GetID()
681 fNTrksToSkip = n; fTrksToSkip = new Int_t[n];
682 for(Int_t i=0;i<n;i++) fTrksToSkip[i] = skipped[i];
685 //---------------------------------------------------------------------------
686 void AliVertexerTracks::SetVtxStart(AliESDVertex *vtx)
689 // Set initial vertex knowledge
691 vtx->GetXYZ(fNominalPos);
692 vtx->GetCovMatrix(fNominalCov);
696 //---------------------------------------------------------------------------
697 void AliVertexerTracks::StrLinVertexFinderMinDist(Int_t optUseWeights)
700 Double_t field=GetFieldkG();
701 const Int_t knacc = (Int_t)fTrkArray.GetEntriesFast();
702 TClonesArray *linarray = new TClonesArray("AliStrLine",1000);
703 TClonesArray &lines = *linarray;
704 for(Int_t i=0; i<knacc; i++){
705 track1 = (AliESDtrack*)fTrkArray.At(i);
706 Double_t alpha=track1->GetAlpha();
707 Double_t mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
708 Double_t pos[3],dir[3],sigmasq[3];
709 track1->GetXYZAt(mindist,field,pos);
710 track1->GetPxPyPzAt(mindist,field,dir);
711 sigmasq[0]=TMath::Sin(alpha)*TMath::Sin(alpha)*track1->GetSigmaY2();
712 sigmasq[1]=TMath::Cos(alpha)*TMath::Cos(alpha)*track1->GetSigmaY2();
713 sigmasq[2]=track1->GetSigmaZ2();
716 if(!TrackToPoint(track1,ri,wWi)) continue;
719 for(Int_t ia=0;ia<3;ia++){
720 for(Int_t ib=0;ib<3;ib++){
721 wmat[iel]=wWi(ia,ib);
725 new(lines[i]) AliStrLine(pos,sigmasq,wmat,dir);
727 fVert=TrackletVertexFinder(linarray,optUseWeights);
731 //---------------------------------------------------------------------------
732 AliESDVertex AliVertexerTracks::TrackletVertexFinder(TClonesArray *lines, Int_t optUseWeights)
734 // Calculate the point at minimum distance to prepared tracks
736 const Int_t knacc = (Int_t)lines->GetEntriesFast();
737 Double_t initPos[3]={0.,0.,0.};
739 Double_t (*vectP0)[3]=new Double_t [knacc][3];
740 Double_t (*vectP1)[3]=new Double_t [knacc][3];
743 Double_t dsum[3]={0,0,0};
745 for(Int_t i=0;i<3;i++){
746 for(Int_t j=0;j<3;j++){
752 for(Int_t i=0; i<knacc; i++){
753 AliStrLine* line1 = (AliStrLine*)lines->At(i);
754 Double_t p0[3],cd[3],sigmasq[3];
758 line1->GetSigma2P0(sigmasq);
759 line1->GetWMatrix(wmat);
762 for(Int_t ia=0;ia<3;ia++){
763 for(Int_t ib=0;ib<3;ib++){
764 wWi(ia,ib)=wmat[iel];
771 Double_t p1[3]={p0[0]+cd[0],p0[1]+cd[1],p0[2]+cd[2]};
781 if(optUseWeights==0)GetStrLinDerivMatrix(p0,p1,matr,dknow);
782 else GetStrLinDerivMatrix(p0,p1,sigmasq,matr,dknow);
785 for(Int_t iii=0;iii<3;iii++){
786 dsum[iii]+=dknow[iii];
787 for(Int_t lj=0;lj<3;lj++) sum[iii][lj]+=matr[iii][lj];
791 TMatrixD invsumWi(TMatrixD::kInverted,sumWi);
792 Double_t covmatrix[6];
793 covmatrix[0] = invsumWi(0,0);
794 covmatrix[1] = invsumWi(0,1);
795 covmatrix[2] = invsumWi(1,1);
796 covmatrix[3] = invsumWi(0,2);
797 covmatrix[4] = invsumWi(1,2);
798 covmatrix[5] = invsumWi(2,2);
801 Double_t det=GetDeterminant3X3(sum);
805 for(Int_t zz=0;zz<3;zz++){
806 for(Int_t ww=0;ww<3;ww++){
807 for(Int_t kk=0;kk<3;kk++) vett[ww][kk]=sum[ww][kk];
809 for(Int_t kk=0;kk<3;kk++) vett[kk][zz]=dsum[kk];
810 initPos[zz]=GetDeterminant3X3(vett)/det;
814 for(Int_t i=0; i<knacc; i++){
815 Double_t p0[3]={0,0,0},p1[3]={0,0,0};
816 for(Int_t ii=0;ii<3;ii++){
817 p0[ii]=vectP0[i][ii];
818 p1[ii]=vectP1[i][ii];
820 sigma+=GetStrLinMinDist(p0,p1,initPos);
823 sigma=TMath::Sqrt(sigma);
827 AliESDVertex theVert(initPos,covmatrix,99999.,knacc);
828 theVert.SetDispersion(sigma);
833 //---------------------------------------------------------------------------
834 Bool_t AliVertexerTracks::TrackToPoint(AliESDtrack *t,
835 TMatrixD &ri,TMatrixD &wWi,
836 Bool_t uUi3by3) const
839 // Extract from the AliESDtrack the global coordinates ri and covariance matrix
840 // wWi of the space point that it represents (to be used in VertexFitter())
844 Double_t rotAngle = t->GetAlpha();
845 if(rotAngle<0.) rotAngle += 2.*TMath::Pi();
846 Double_t cosRot = TMath::Cos(rotAngle);
847 Double_t sinRot = TMath::Sin(rotAngle);
849 ri(0,0) = t->GetX()*cosRot-t->GetY()*sinRot;
850 ri(1,0) = t->GetX()*sinRot+t->GetY()*cosRot;
856 // matrix to go from global (x,y,z) to local (y,z);
865 // covariance matrix of local (y,z) - inverted
867 t->GetExternalCovariance(cc);
874 //printf(" %f %f\n",uUi(0,0),uUi(0,1));
875 //printf(" %f %f\n",uUi(1,0),uUi(1,1));
877 if(uUi.Determinant() <= 0.) return kFALSE;
878 TMatrixD uUiInv(TMatrixD::kInverted,uUi);
880 // weights matrix: wWi = qQiT * uUiInv * qQi
881 TMatrixD uUiInvQi(uUiInv,TMatrixD::kMult,qQi);
882 TMatrixD m(qQi,TMatrixD::kTransposeMult,uUiInvQi);
885 if(fVert.GetNContributors()<1) AliFatal("Vertex from finder is empty");
886 // matrix to go from global (x,y,z) to local (x,y,z);
898 // covariance of fVert along the track
899 Double_t p[3],pt,ptot;
901 pt = TMath::Sqrt(p[0]*p[0]+p[1]*p[1]);
902 ptot = TMath::Sqrt(pt*pt+p[2]*p[2]);
903 Double_t cphi = p[0]/pt; //cos(phi)=px/pt
904 Double_t sphi = p[1]/pt; //sin(phi)=py/pt
905 Double_t clambda = pt/ptot; //cos(lambda)=pt/ptot
906 Double_t slambda = p[2]/ptot; //sin(lambda)=pz/ptot
907 Double_t covfVert[6];
908 fVert.GetCovMatrix(covfVert);
909 Double_t covfVertalongt =
910 covfVert[0]*cphi*cphi*clambda*clambda
911 +covfVert[1]*2.*cphi*sphi*clambda*clambda
912 +covfVert[3]*2.*cphi*clambda*slambda
913 +covfVert[2]*sphi*sphi*clambda*clambda
914 +covfVert[4]*2.*sphi*clambda*slambda
915 +covfVert[5]*slambda*slambda;
917 t->GetExternalCovariance(cc);
918 // covariance matrix of local (x,y,z) - inverted
920 uUi(0,0) = covfVertalongt * fnSigmaForUi00 * fnSigmaForUi00;
921 if(fDebug) printf("=====> sqrtUi00 cm %f\n",TMath::Sqrt(uUi(0,0)));
931 //printf(" %f %f\n",uUi(0,0),uUi(0,1));
932 //printf(" %f %f\n",uUi(1,0),uUi(1,1));
934 if(uUi.Determinant() <= 0.) return kFALSE;
935 TMatrixD uUiInv(TMatrixD::kInverted,uUi);
937 // weights matrix: wWi = qQiT * uUiInv * qQi
938 TMatrixD uUiInvQi(uUiInv,TMatrixD::kMult,qQi);
939 TMatrixD m(qQi,TMatrixD::kTransposeMult,uUiInvQi);
946 //---------------------------------------------------------------------------
947 void AliVertexerTracks::TooFewTracks(const AliESDEvent* esdEvent)
950 // When the number of tracks is < fMinTracks
953 // deal with vertices not found
954 Double_t pos[3],err[3];
956 pos[0] = fNominalPos[0];
957 err[0] = TMath::Sqrt(fNominalCov[0]);
958 pos[1] = fNominalPos[1];
959 err[1] = TMath::Sqrt(fNominalCov[2]);
960 pos[2] = esdEvent->GetVertex()->GetZv();
961 err[2] = esdEvent->GetVertex()->GetZRes();
962 if(err[0]>1. && esdEvent->GetVertex()->GetNContributors()<=0)
963 ncontr = -1; // (x,y,z) = (0,0,0)
964 if(err[0]>1. && esdEvent->GetVertex()->GetNContributors()>0)
965 ncontr = -2; // (x,y,z) = (0,0,z_fromSPD)
966 if(err[0]<1. && esdEvent->GetVertex()->GetNContributors()<=0)
967 ncontr = -3; // (x,y,z) = (x_mean,y_mean,0)
968 if(err[0]<1. && esdEvent->GetVertex()->GetNContributors()>0)
969 ncontr = -4; // (x,y,z) = (x_mean,y_mean,z_fromSPD)
971 fCurrentVertex = new AliESDVertex(pos,err);
972 fCurrentVertex->SetNContributors(ncontr);
975 fCurrentVertex->SetTitle("VertexerTracksWithConstraint");
977 fCurrentVertex->SetTitle("VertexerTracksNoConstraint");
982 //---------------------------------------------------------------------------
983 void AliVertexerTracks::VertexFinder(Int_t optUseWeights)
986 // Get estimate of vertex position in (x,y) from tracks DCA
990 for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
991 Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
992 Double_t aver[3]={0.,0.,0.};
993 Double_t aversq[3]={0.,0.,0.};
994 Double_t sigmasq[3]={0.,0.,0.};
999 Double_t pos[3],dir[3];
1000 Double_t alpha,mindist;
1001 Double_t field=GetFieldkG();
1003 for(Int_t i=0; i<nacc; i++){
1004 track1 = (AliESDtrack*)fTrkArray.At(i);
1005 alpha=track1->GetAlpha();
1006 mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
1007 track1->GetXYZAt(mindist,field,pos);
1008 track1->GetPxPyPzAt(mindist,field,dir);
1009 AliStrLine *line1 = new AliStrLine(pos,dir);
1011 // AliStrLine *line1 = new AliStrLine();
1012 // track1->ApproximateHelixWithLine(mindist,field,line1);
1014 for(Int_t j=i+1; j<nacc; j++){
1015 track2 = (AliESDtrack*)fTrkArray.At(j);
1016 alpha=track2->GetAlpha();
1017 mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
1018 track2->GetXYZAt(mindist,field,pos);
1019 track2->GetPxPyPzAt(mindist,field,dir);
1020 AliStrLine *line2 = new AliStrLine(pos,dir);
1021 // AliStrLine *line2 = new AliStrLine();
1022 // track2->ApproximateHelixWithLine(mindist,field,line2);
1023 Double_t distCA=line2->GetDCA(line1);
1024 //printf("%d %d %f\n",i,j,distCA);
1025 if(fDCAcut<=0 || (fDCAcut>0&&distCA<fDCAcut)){
1026 Double_t pnt1[3],pnt2[3],crosspoint[3];
1028 if(optUseWeights<=0){
1029 Int_t retcode = line2->Cross(line1,crosspoint);
1032 for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
1033 for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]);
1036 if(optUseWeights>0){
1037 Int_t retcode = line1->CrossPoints(line2,pnt1,pnt2);
1040 alpha=track1->GetAlpha();
1041 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
1042 Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2();
1043 alpha=track2->GetAlpha();
1044 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
1045 Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2();
1046 Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2();
1047 Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
1048 Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
1049 Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
1050 crosspoint[0]=wx1*pnt1[0] + wx2*pnt2[0];
1051 crosspoint[1]=wy1*pnt1[1] + wy2*pnt2[1];
1052 crosspoint[2]=wz1*pnt1[2] + wz2*pnt2[2];
1055 for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
1056 for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]);
1065 for(Int_t jj=0;jj<3;jj++){
1066 initPos[jj] = aver[jj]/ncombi;
1067 //printf("%f\n",initPos[jj]);
1069 sigmasq[jj]=aversq[jj]-initPos[jj]*initPos[jj];
1072 sigma=TMath::Sqrt(TMath::Abs(sigma));
1075 Warning("VertexFinder","Finder did not succed");
1078 fVert.SetXYZ(initPos);
1079 fVert.SetDispersion(sigma);
1080 fVert.SetNContributors(ncombi);
1082 //---------------------------------------------------------------------------
1083 void AliVertexerTracks::VertexFitter(Bool_t useConstraint)
1086 // The optimal estimate of the vertex position is given by a "weighted
1087 // average of tracks positions".
1088 // Original method: V. Karimaki, CMS Note 97/0051
1090 Double_t initPos[3];
1092 fVert.GetXYZ(initPos);
1094 initPos[0]=fNominalPos[0];
1095 initPos[1]=fNominalPos[1];
1096 initPos[2]=fNominalPos[2];
1099 Int_t nTrks = (Int_t)fTrkArray.GetEntries();
1100 if(nTrks==1) useConstraint=kTRUE;
1102 printf("--- VertexFitter(): start\n");
1103 printf(" Number of tracks in array: %d\n",nTrks);
1104 printf(" Minimum # tracks required in fit: %d\n",fMinTracks);
1105 printf(" Vertex position after finder: %f,%f,%f\n",initPos[0],initPos[1],initPos[2]);
1106 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]));
1109 // special treatment for few-tracks fits (e.g. secondary vertices)
1110 Bool_t uUi3by3 = kFALSE; if(nTrks<5 && !useConstraint) uUi3by3 = kTRUE;
1115 rv(0,0) = initPos[0];
1116 rv(1,0) = initPos[1];
1118 Double_t xlStart,alpha;
1120 Double_t chi2,chi2i,chi2b;
1124 // initial vertex covariance matrix
1126 vVb(0,0) = fNominalCov[0];
1127 vVb(0,1) = fNominalCov[1];
1129 vVb(1,0) = fNominalCov[1];
1130 vVb(1,1) = fNominalCov[2];
1134 vVb(2,2) = fNominalCov[5];
1135 TMatrixD vVbInv(TMatrixD::kInverted,vVb);
1137 rb(0,0) = fNominalPos[0];
1138 rb(1,0) = fNominalPos[1];
1139 rb(2,0) = fNominalPos[2];
1140 TMatrixD vVbInvrb(vVbInv,TMatrixD::kMult,rb);
1144 // 1st - estimate of vtx using all tracks
1145 // 2nd - estimate of global chi2
1146 for(step=0; step<2; step++) {
1147 if(fDebug) printf(" step = %d\n",step);
1151 if(step==1) { initPos[0]=rv(0,0); initPos[0]=rv(1,0); }
1153 TMatrixD sumWiri(3,1);
1154 TMatrixD sumWi(3,3);
1155 for(i=0; i<3; i++) {
1157 for(j=0; j<3; j++) sumWi(j,i) = 0.;
1160 // mean vertex constraint
1163 sumWiri(i,0) += vVbInvrb(i,0);
1164 for(k=0;k<3;k++) sumWi(i,k) += vVbInv(i,k);
1167 TMatrixD deltar = rv; deltar -= rb;
1168 TMatrixD vVbInvdeltar(vVbInv,TMatrixD::kMult,deltar);
1169 chi2b = deltar(0,0)*vVbInvdeltar(0,0)+
1170 deltar(1,0)*vVbInvdeltar(1,0)+
1171 deltar(2,0)*vVbInvdeltar(2,0);
1176 for(k=0; k<nTrks; k++) {
1177 // get track from track array
1178 t = (AliESDtrack*)fTrkArray.At(k);
1179 alpha = t->GetAlpha();
1180 xlStart = initPos[0]*TMath::Cos(alpha)+initPos[1]*TMath::Sin(alpha);
1181 // to vtxSeed (from finder)
1182 t->AliExternalTrackParam::PropagateTo(xlStart,GetFieldkG());
1184 // vector of track global coordinates
1186 // covariance matrix of ri
1189 // get space point from track
1190 if(!TrackToPoint(t,ri,wWi,uUi3by3)) continue;
1191 TMatrixD wWiri(wWi,TMatrixD::kMult,ri);
1194 TMatrixD deltar = rv; deltar -= ri;
1195 TMatrixD wWideltar(wWi,TMatrixD::kMult,deltar);
1196 chi2i = deltar(0,0)*wWideltar(0,0)+
1197 deltar(1,0)*wWideltar(1,0)+
1198 deltar(2,0)*wWideltar(2,0);
1200 // add to total chi2
1207 } // end loop on tracks
1209 if(nUsedTrks < fMinTracks) {
1214 Double_t determinant = sumWi.Determinant();
1215 //cerr<<" determinant: "<<determinant<<endl;
1216 if(determinant < 100.) {
1217 printf("det(V) = 0\n");
1223 // inverted of weights matrix
1224 TMatrixD invsumWi(TMatrixD::kInverted,sumWi);
1226 // position of primary vertex
1227 rv.Mult(vV,sumWiri);
1229 } // end loop on the 2 steps
1232 if(fDebug) printf("TooFewTracks\n");
1233 fCurrentVertex = new AliESDVertex(0.,0.,-1);
1237 Double_t position[3];
1238 position[0] = rv(0,0);
1239 position[1] = rv(1,0);
1240 position[2] = rv(2,0);
1241 Double_t covmatrix[6];
1242 covmatrix[0] = vV(0,0);
1243 covmatrix[1] = vV(0,1);
1244 covmatrix[2] = vV(1,1);
1245 covmatrix[3] = vV(0,2);
1246 covmatrix[4] = vV(1,2);
1247 covmatrix[5] = vV(2,2);
1249 // for correct chi2/ndf, count constraint as additional "track"
1250 if(fConstraint) nUsedTrks++;
1251 // store data in the vertex object
1252 fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nUsedTrks);
1255 printf(" Vertex after fit:\n");
1256 fCurrentVertex->PrintStatus();
1257 printf("--- VertexFitter(): finish\n");
1262 //----------------------------------------------------------------------------
1263 AliESDVertex* AliVertexerTracks::VertexForSelectedTracks(TTree *trkTree,Bool_t optUseFitter,Bool_t optPropagate)
1266 // Return vertex from tracks in trkTree
1271 // get tracks and propagate them to initial vertex position
1272 Int_t nTrksPrep = PrepareTracks(*trkTree,0);
1273 if(nTrksPrep < TMath::Max(2,fMinTracks) ) {
1274 if(fDebug) printf("TooFewTracks\n");
1275 fCurrentVertex = new AliESDVertex(0.,0.,-1);
1276 return fCurrentVertex;
1280 case 1: StrLinVertexFinderMinDist(1); break;
1281 case 2: StrLinVertexFinderMinDist(0); break;
1282 case 3: HelixVertexFinder(); break;
1283 case 4: VertexFinder(1); break;
1284 case 5: VertexFinder(0); break;
1285 default: printf("Wrong algorithm\n"); break;
1288 if(fDebug) printf(" Vertex finding completed\n");
1292 //SetVtxStart(&fVert);
1293 VertexFitter(fConstraint);
1294 if(fDebug) printf(" Vertex fit completed\n");
1296 Double_t position[3]={fVert.GetXv(),fVert.GetYv(),fVert.GetZv()};
1297 Double_t covmatrix[6];
1298 fVert.GetCovMatrix(covmatrix);
1299 Double_t chi2=99999.;
1300 Int_t nUsedTrks=fVert.GetNContributors();
1301 fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nUsedTrks);
1303 fCurrentVertex->SetDispersion(fVert.GetDispersion());
1306 // set indices of used tracks and propagate track to found vertex
1307 UShort_t *indices = 0;
1308 AliESDtrack *eta = 0;
1309 if(fCurrentVertex->GetNContributors()>0) {
1310 indices = new UShort_t[fCurrentVertex->GetNContributors()];
1311 for(Int_t jj=0;jj<(Int_t)fTrkArray.GetEntriesFast();jj++) {
1312 eta = (AliESDtrack*)fTrkArray.At(jj);
1313 indices[jj] = (UShort_t)eta->GetID();
1314 if(optPropagate&&optUseFitter){
1315 if(TMath::Sqrt(fCurrentVertex->GetXv()*fCurrentVertex->GetXv()+fCurrentVertex->GetYv()*fCurrentVertex->GetYv())<3.) {
1316 eta->RelateToVertex(fCurrentVertex,GetFieldkG(),100.);
1317 if(fDebug) printf("Track %d propagated to found vertex\n",jj);
1319 AliWarning("Found vertex outside beam pipe!");
1323 fCurrentVertex->SetIndices(fCurrentVertex->GetNContributors(),indices);
1328 return fCurrentVertex;
1330 //----------------------------------------------------------------------------
1331 AliESDVertex* AliVertexerTracks::VertexForSelectedTracks(TObjArray *trkArray,Bool_t optUseFitter, Bool_t optPropagate)
1334 // Return vertex from array of tracks
1337 // get tracks and propagate them to initial vertex position
1338 Int_t nTrks = trkArray->GetEntriesFast();
1339 if(nTrks < TMath::Max(2,fMinTracks) ) {
1340 if(fDebug) printf("TooFewTracks\n");
1341 fCurrentVertex = new AliESDVertex(0.,0.,-1);
1342 return fCurrentVertex;
1344 TTree *trkTree = new TTree("TreeT","tracks");
1345 AliESDtrack *esdTrack = 0;
1346 trkTree->Branch("tracks","AliESDtrack",&esdTrack);
1347 for(Int_t i=0; i<nTrks; i++){
1348 esdTrack = (AliESDtrack*)trkArray->At(i);
1352 AliESDVertex *vtx = VertexForSelectedTracks(trkTree,optUseFitter,optPropagate);
1356 //--------------------------------------------------------------------------