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():
45 fMaxChi2PerTrack(1e33),
54 // Default constructor
62 //-----------------------------------------------------------------------------
63 AliVertexerTracks::AliVertexerTracks(Double_t xStart, Double_t yStart):
67 fMaxChi2PerTrack(1e33),
76 // Standard constructor
78 SetVtxStart(xStart,yStart);
84 //-----------------------------------------------------------------------------
85 AliVertexerTracks::~AliVertexerTracks() {
87 // The objects poited by the following pointer are not owned
88 // by this class and are not deleted
92 //---------------------------------------------------------------------------
93 void AliVertexerTracks::SetSkipTracks(Int_t n,Int_t *skipped) {
95 // Mark the tracks not ot be used in the vertex finding
98 fTrksToSkip = new Int_t[n];
99 for(Int_t i=0;i<n;i++) fTrksToSkip[i] = skipped[i];
102 //---------------------------------------------------------------------------
103 void AliVertexerTracks::ComputeMaxChi2PerTrack(Int_t nTracks) {
105 // Max. contr. to the chi2 has been tuned as a function of multiplicity
107 if(nTracks < 7) { fMaxChi2PerTrack = 1.e6;
108 } else { fMaxChi2PerTrack = 100.; }
112 //----------------------------------------------------------------------------
113 Int_t AliVertexerTracks::PrepareTracks(TTree &trkTree, Int_t OptImpParCut) {
115 // Propagate tracks to initial vertex position and store them in a TObjArray
117 Double_t maxd0rphi = 3.;
120 Double_t sigmaCurr[3],sigmaVtx[3],posVtx[3];
121 Float_t d0z0[2],covd0z0[3];
122 Double_t momentum[3],postrk[3];
123 Double_t pt,sigma,sigmavtxphi,phitrk;
124 Double_t field=GetField();
126 AliESDVertex *initVertex = new AliESDVertex(fNominalPos,fNominalSigma);
128 Int_t nEntries = (Int_t)trkTree.GetEntries();
129 if(!fTrkArray.IsEmpty()) fTrkArray.Clear();
130 fTrkArray.Expand(nEntries);
133 printf(" PrepareTracks()\n");
136 for(Int_t i=0; i<nEntries; i++) {
137 // check tracks to skip
139 for(Int_t j=0; j<fNTrksToSkip; j++) {
140 if(i==fTrksToSkip[j]) {
141 if(fDebug) printf("skipping track: %d\n",i);
145 if(skipThis) continue;
147 AliESDtrack *track = new AliESDtrack;
148 trkTree.SetBranchAddress("tracks",&track);
153 // propagate track to vertex
154 if(OptImpParCut==1) {
155 track->RelateToVertex(initVertex,field,100.);
156 initVertex->GetXYZ(posVtx);
157 initVertex->GetSigmaXYZ(sigmaVtx);
159 fCurrentVertex->GetSigmaXYZ(sigmaCurr);
160 if((sigmaCurr[0]+sigmaCurr[1])<(fNominalSigma[0]+fNominalSigma[1])) {
161 track->RelateToVertex(fCurrentVertex,field,100.);
162 fCurrentVertex->GetXYZ(posVtx);
163 fCurrentVertex->GetSigmaXYZ(sigmaVtx);
165 track->RelateToVertex(initVertex,field,100.);
166 initVertex->GetXYZ(posVtx);
167 initVertex->GetSigmaXYZ(sigmaVtx);
171 // select tracks with d0rphi < maxd0rphi
172 track->GetImpactParameters(d0z0,covd0z0);
174 if(OptImpParCut==1) {
175 maxd0rphi = 3.; // cm
177 track->GetXYZ(postrk);
178 track->GetPxPyPz(momentum);
179 pt = TMath::Sqrt(momentum[0]*momentum[0]+momentum[1]*momentum[1]);
181 phitrk = TMath::ATan2(postrk[1]-posVtx[1],postrk[0]-posVtx[0]);
182 sigmavtxphi = TMath::Sqrt(sigmaVtx[0]*sigmaVtx[0]*
183 TMath::Cos(phitrk)*TMath::Cos(phitrk)+
184 sigmaVtx[1]*sigmaVtx[1]*
185 TMath::Sin(phitrk)*TMath::Sin(phitrk));
186 sigma = TMath::Sqrt(Sigmad0rphi(pt)*Sigmad0rphi(pt)+
187 sigmavtxphi*sigmavtxphi);
188 maxd0rphi = fNSigma*sigma;
190 if(TMath::Abs(d0z0[0]) > maxd0rphi) {
191 delete track; continue;
196 // propagate track to vtxSeed
197 alpha = track->GetAlpha();
198 xlStart = fNominalPos[0]*TMath::Cos(alpha)+fNominalPos[1]*TMath::Sin(alpha);
199 track->PropagateTo(xlStart,field); // to vtxSeed
200 d0rphi = TMath::Abs(track->GetD(fNominalPos[0],fNominalPos[1],field));
201 if(OptImpParCut==1 && d0rphi > maxd0rphi) { delete track; continue; }
202 if(OptImpParCut==2) {
203 // to be replaced by new Andrea's selection
204 if(d0rphi > maxd0rphi){
208 // end to be replaced by new Andrea's selection
212 fTrkArray.AddLast(track);
214 if(fDebug) printf(" :-) nTrks=%d, d0rphi=%f\n",nTrks,TMath::Abs(d0z0[0]));
217 if(fTrksToSkip) delete [] fTrksToSkip;
222 //----------------------------------------------------------------------------
223 Double_t AliVertexerTracks::Sigmad0rphi(Double_t pt) const {
225 // Impact parameter resolution in rphi [cm] vs pt [GeV/c]
227 Double_t a_meas = 11.6;
228 Double_t a_scatt = 65.8;
231 Double_t sigma = a_meas*a_meas+a_scatt*a_scatt/TMath::Power(pt,b);
232 sigma = 0.0001*TMath::Sqrt(sigma);
236 //----------------------------------------------------------------------------
237 AliVertex* AliVertexerTracks::VertexForSelectedTracks(TTree *trkTree) {
239 // Return vertex from tracks in trkTree
242 // get tracks and propagate them to initial vertex position
243 Int_t nTrks = PrepareTracks(*trkTree,1);
244 if(nTrks < fMinTracks) {
245 printf("TooFewTracks\n");
246 Double_t vtx[3]={0,0,0};
248 fVert.SetDispersion(999);
249 fVert.SetNContributors(-5);
253 // Set initial vertex position from ESD
254 if(fAlgo==1) StrLinVertexFinderMinDist(1);
255 if(fAlgo==2) StrLinVertexFinderMinDist(0);
256 if(fAlgo==3) HelixVertexFinder();
257 if(fAlgo==4) VertexFinder(1);
258 if(fAlgo==5) VertexFinder(0);
261 //----------------------------------------------------------------------------
262 AliESDVertex* AliVertexerTracks::FindPrimaryVertex(const AliESD *esdEvent)
265 // Primary vertex for current ESD event
266 // (Two iterations: 1st with |d0rphi|<3cm; 2nd with
267 // fNSigma*sigma(pt) cut w.r.t. to vertex found in 1st iteration)
271 Int_t entr = (Int_t)esdEvent->GetNumberOfTracks();
272 TTree *trkTree = new TTree("TreeT","tracks");
273 AliESDtrack *esdTrack = 0;
274 trkTree->Branch("tracks","AliESDtrack",&esdTrack);
276 for(Int_t i=0; i<entr; i++) {
277 AliESDtrack *et = esdEvent->GetTrack(i);
278 esdTrack = new AliESDtrack(*et);
279 if(!esdTrack->GetStatus()&AliESDtrack::kITSin) continue;
280 if(!esdTrack->GetStatus()&AliESDtrack::kITSrefit) continue;
281 Int_t nclus=esdTrack->GetNcls(0); // check number of clusters in ITS
282 if(nclus<fMinITSClusters) continue;
289 // propagate tracks to initVertex
290 // preselect them (reject only |d0|>3cm w.r.t. finitVertex)
292 nTrks = PrepareTracks(*trkTree,1);
293 if(fDebug) printf(" tracks prepared - step 1: %d\n",nTrks);
294 if(nTrks < fMinTracks) {
295 printf("TooFewTracks\n");
296 fCurrentVertex = new AliESDVertex(0.,0.,-1);
297 return fCurrentVertex;
302 case 1: StrLinVertexFinderMinDist(1); break;
303 case 2: StrLinVertexFinderMinDist(0); break;
304 case 3: HelixVertexFinder(); break;
305 case 4: VertexFinder(1); break;
306 case 5: VertexFinder(0); break;
307 default: printf("Wrong algorithm\n"); break;
309 if(fDebug) printf(" vertex finding completed\n");
312 ComputeMaxChi2PerTrack(nTrks);
314 if(fDebug) printf(" vertex fit completed\n");
317 // propagate tracks to best between initVertex and fCurrentVertex
318 // preselect tracks (reject for |d0|>fNSigma*sigma w.r.t. best
319 // between initVertex and fCurrentVertex)
320 nTrks = PrepareTracks(*trkTree,2);
322 if(fDebug) printf(" tracks prepared - step 2: %d\n",nTrks);
323 if(nTrks < fMinTracks) {
324 printf("TooFewTracks\n");
325 fCurrentVertex = new AliESDVertex(0.,0.,-1);
326 return fCurrentVertex;
331 case 1: StrLinVertexFinderMinDist(1); break;
332 case 2: StrLinVertexFinderMinDist(0); break;
333 case 3: HelixVertexFinder(); break;
334 case 4: VertexFinder(1); break;
335 case 5: VertexFinder(0); break;
336 default: printf("Wrong algorithm\n"); break;
338 if(fDebug) printf(" vertex finding completed\n");
341 ComputeMaxChi2PerTrack(nTrks);
343 if(fDebug) printf(" vertex fit completed\n");
346 esdEvent->GetVertex()->GetTruePos(tp);
347 fCurrentVertex->SetTruePos(tp);
348 fCurrentVertex->SetTitle("VertexerTracks");
351 return fCurrentVertex;
353 //----------------------------------------------------------------------------
354 AliESDVertex* AliVertexerTracks::FindPrimaryVertexOld(const AliESD *esdEvent)
357 // Primary vertex for current ESD event (one iteration with |d0rphi|<3cm)
361 Int_t entr = (Int_t)esdEvent->GetNumberOfTracks();
362 TTree *trkTree = new TTree("TreeT","tracks");
363 AliESDtrack *esdTrack = 0;
364 trkTree->Branch("tracks","AliESDtrack",&esdTrack);
366 for(Int_t i=0; i<entr; i++) {
367 AliESDtrack *et = esdEvent->GetTrack(i);
368 esdTrack = new AliESDtrack(*et);
369 if(!esdTrack->GetStatus()&AliESDtrack::kITSin) continue;
370 if(!esdTrack->GetStatus()&AliESDtrack::kITSrefit) continue;
371 Int_t nclus=esdTrack->GetNcls(0); // check number of clusters in ITS
372 if(nclus<fMinITSClusters) continue;
378 // preselect tracks and propagate them to initial vertex position
379 Int_t nTrks = PrepareTracks(*trkTree,1);
381 if(fDebug) printf(" tracks prepared: %d\n",nTrks);
382 if(nTrks < fMinTracks) {
383 printf("TooFewTracks\n");
384 fCurrentVertex = new AliESDVertex(0.,0.,-1);
385 return fCurrentVertex;
390 case 1: StrLinVertexFinderMinDist(1); break;
391 case 2: StrLinVertexFinderMinDist(0); break;
392 case 3: HelixVertexFinder(); break;
393 case 4: VertexFinder(1); break;
394 case 5: VertexFinder(0); break;
395 default: printf("Wrong algorithm\n"); break;
399 ComputeMaxChi2PerTrack(nTrks);
401 if(fDebug) printf(" vertex fit completed\n");
405 esdEvent->GetVertex()->GetTruePos(tp);
406 fCurrentVertex->SetTruePos(tp);
407 fCurrentVertex->SetTitle("VertexerTracks");
410 return fCurrentVertex;
412 //----------------------------------------------------------------------------
413 AliVertex* AliVertexerTracks::VertexForSelectedTracks(TObjArray *trkArray) {
415 // Return vertex from array of tracks
418 // get tracks and propagate them to initial vertex position
419 Int_t nTrks = trkArray->GetEntriesFast();
420 if(nTrks < fMinTracks) {
421 printf("TooFewTracks\n");
422 Double_t vtx[3]={0,0,0};
424 fVert.SetDispersion(999);
425 fVert.SetNContributors(-5);
428 TTree *trkTree = new TTree("TreeT","tracks");
429 AliESDtrack *esdTrack = 0;
430 trkTree->Branch("tracks","AliESDtrack",&esdTrack);
431 for(Int_t i=0; i<nTrks; i++){
432 esdTrack = (AliESDtrack*)trkArray->At(i);
436 AliVertex *vtx = VertexForSelectedTracks(trkTree);
441 //---------------------------------------------------------------------------
442 void AliVertexerTracks::VertexFinder(Int_t optUseWeights) {
444 // Get estimate of vertex position in (x,y) from tracks DCA
448 for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
449 Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
450 Double_t aver[3]={0.,0.,0.};
451 Double_t aversq[3]={0.,0.,0.};
452 Double_t sigmasq[3]={0.,0.,0.};
457 Double_t pos[3],dir[3];
458 Double_t alpha,mindist;
459 Double_t field=GetField();
461 for(Int_t i=0; i<nacc; i++){
462 track1 = (AliESDtrack*)fTrkArray.At(i);
463 alpha=track1->GetAlpha();
464 mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
465 track1->GetXYZAt(mindist,field,pos);
466 track1->GetPxPyPzAt(mindist,field,dir);
467 AliStrLine *line1 = new AliStrLine(pos,dir);
469 // AliStrLine *line1 = new AliStrLine();
470 // track1->ApproximateHelixWithLine(mindist,field,line1);
472 for(Int_t j=i+1; j<nacc; j++){
473 track2 = (AliESDtrack*)fTrkArray.At(j);
474 alpha=track2->GetAlpha();
475 mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
476 track2->GetXYZAt(mindist,field,pos);
477 track2->GetPxPyPzAt(mindist,field,dir);
478 AliStrLine *line2 = new AliStrLine(pos,dir);
479 // AliStrLine *line2 = new AliStrLine();
480 // track2->ApproximateHelixWithLine(mindist,field,line2);
481 Double_t distCA=line2->GetDCA(line1);
482 if(fDCAcut<=0 || (fDCAcut>0&&distCA<fDCAcut)){
483 Double_t pnt1[3],pnt2[3],crosspoint[3];
485 if(optUseWeights<=0){
486 Int_t retcode = line2->Cross(line1,crosspoint);
489 for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
490 for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]);
494 Int_t retcode = line1->CrossPoints(line2,pnt1,pnt2);
496 Double_t alpha, cs, sn;
497 alpha=track1->GetAlpha();
498 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
499 Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2();
500 alpha=track2->GetAlpha();
501 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
502 Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2();
503 Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2();
504 Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
505 Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
506 Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
507 crosspoint[0]=wx1*pnt1[0] + wx2*pnt2[0];
508 crosspoint[1]=wy1*pnt1[1] + wy2*pnt2[1];
509 crosspoint[2]=wz1*pnt1[2] + wz2*pnt2[2];
512 for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
513 for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]);
522 for(Int_t jj=0;jj<3;jj++){
523 initPos[jj] = aver[jj]/ncombi;
525 sigmasq[jj]=aversq[jj]-initPos[jj]*initPos[jj];
528 sigma=TMath::Sqrt(TMath::Abs(sigma));
531 Warning("VertexFinder","Finder did not succed");
534 fVert.SetXYZ(initPos);
535 fVert.SetDispersion(sigma);
536 fVert.SetNContributors(ncombi);
538 //---------------------------------------------------------------------------
539 void AliVertexerTracks::HelixVertexFinder() {
541 // Get estimate of vertex position in (x,y) from tracks DCA
546 for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
547 Double_t field=GetField();
549 Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
551 Double_t aver[3]={0.,0.,0.};
552 Double_t averquad[3]={0.,0.,0.};
553 Double_t sigmaquad[3]={0.,0.,0.};
560 Double_t alpha, cs, sn;
561 Double_t crosspoint[3];
562 for(Int_t i=0; i<nacc; i++){
563 track1 = (AliESDtrack*)fTrkArray.At(i);
566 for(Int_t j=i+1; j<nacc; j++){
567 track2 = (AliESDtrack*)fTrkArray.At(j);
569 distCA=track2->PropagateToDCA(track1,field);
571 if(fDCAcut<=0 ||(fDCAcut>0&&distCA<fDCAcut)){
572 track1->GetExternalParameters(x,par);
573 alpha=track1->GetAlpha();
574 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
575 Double_t x1=x*cs - par[0]*sn;
576 Double_t y1=x*sn + par[0]*cs;
578 Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2();
579 track2->GetExternalParameters(x,par);
580 alpha=track2->GetAlpha();
581 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
582 Double_t x2=x*cs - par[0]*sn;
583 Double_t y2=x*sn + par[0]*cs;
585 Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2();
586 Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2();
587 Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
588 Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
589 Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
590 crosspoint[0]=wx1*x1 + wx2*x2;
591 crosspoint[1]=wy1*y1 + wy2*y2;
592 crosspoint[2]=wz1*z1 + wz2*z2;
595 for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
596 for(Int_t jj=0;jj<3;jj++)averquad[jj]+=(crosspoint[jj]*crosspoint[jj]);
602 for(Int_t jj=0;jj<3;jj++){
603 initPos[jj] = aver[jj]/ncombi;
604 averquad[jj]/=ncombi;
605 sigmaquad[jj]=averquad[jj]-initPos[jj]*initPos[jj];
606 sigma+=sigmaquad[jj];
608 sigma=TMath::Sqrt(TMath::Abs(sigma));
611 Warning("HelixVertexFinder","Finder did not succed");
614 fVert.SetXYZ(initPos);
615 fVert.SetDispersion(sigma);
616 fVert.SetNContributors(ncombi);
618 //---------------------------------------------------------------------------
619 void AliVertexerTracks::StrLinVertexFinderMinDist(Int_t optUseWeights){
621 // Calculate the point at minimum distance to prepared tracks
626 for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
627 const Int_t knacc = (Int_t)fTrkArray.GetEntriesFast();
628 Double_t field=GetField();
631 Double_t (*vectP0)[3]=new Double_t [knacc][3];
632 Double_t (*vectP1)[3]=new Double_t [knacc][3];
635 Double_t dsum[3]={0,0,0};
636 for(Int_t i=0;i<3;i++)
637 for(Int_t j=0;j<3;j++)sum[i][j]=0;
638 for(Int_t i=0; i<knacc; i++){
639 track1 = (AliESDtrack*)fTrkArray.At(i);
640 Double_t alpha=track1->GetAlpha();
641 Double_t mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
642 Double_t pos[3],dir[3];
643 track1->GetXYZAt(mindist,field,pos);
644 track1->GetPxPyPzAt(mindist,field,dir);
645 AliStrLine *line1 = new AliStrLine(pos,dir);
646 // AliStrLine *line1 = new AliStrLine();
647 // track1->ApproximateHelixWithLine(mindist,field,line1);
649 Double_t p0[3],cd[3];
652 Double_t p1[3]={p0[0]+cd[0],p0[1]+cd[1],p0[2]+cd[2]};
662 if(optUseWeights==0)GetStrLinDerivMatrix(p0,p1,matr,dknow);
663 if(optUseWeights==1){
665 sigmasq[0]=track1->GetSigmaY2();
666 sigmasq[1]=track1->GetSigmaY2();
667 sigmasq[2]=track1->GetSigmaZ2();
668 GetStrLinDerivMatrix(p0,p1,sigmasq,matr,dknow);
671 for(Int_t iii=0;iii<3;iii++){
672 dsum[iii]+=dknow[iii];
673 for(Int_t lj=0;lj<3;lj++) sum[iii][lj]+=matr[iii][lj];
679 Double_t det=GetDeterminant3X3(sum);
682 for(Int_t zz=0;zz<3;zz++){
683 for(Int_t ww=0;ww<3;ww++){
684 for(Int_t kk=0;kk<3;kk++) vett[ww][kk]=sum[ww][kk];
686 for(Int_t kk=0;kk<3;kk++) vett[kk][zz]=dsum[kk];
687 initPos[zz]=GetDeterminant3X3(vett)/det;
691 for(Int_t i=0; i<knacc; i++){
692 Double_t p0[3]={0,0,0},p1[3]={0,0,0};
693 for(Int_t ii=0;ii<3;ii++){
694 p0[ii]=vectP0[i][ii];
695 p1[ii]=vectP1[i][ii];
697 sigma+=GetStrLinMinDist(p0,p1,initPos);
700 sigma=TMath::Sqrt(sigma);
702 Warning("StrLinVertexFinderMinDist","Finder did not succed");
707 fVert.SetXYZ(initPos);
708 fVert.SetDispersion(sigma);
709 fVert.SetNContributors(knacc);
711 //_______________________________________________________________________
712 Double_t AliVertexerTracks::GetDeterminant3X3(Double_t matr[][3]){
714 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];
717 //____________________________________________________________________________
718 void AliVertexerTracks::GetStrLinDerivMatrix(Double_t *p0,Double_t *p1,Double_t (*m)[3],Double_t *d){
721 Double_t x12=p0[0]-p1[0];
722 Double_t y12=p0[1]-p1[1];
723 Double_t z12=p0[2]-p1[2];
724 Double_t kk=x12*x12+y12*y12+z12*z12;
725 m[0][0]=2-2/kk*x12*x12;
726 m[0][1]=-2/kk*x12*y12;
727 m[0][2]=-2/kk*x12*z12;
728 m[1][0]=-2/kk*x12*y12;
729 m[1][1]=2-2/kk*y12*y12;
730 m[1][2]=-2/kk*y12*z12;
731 m[2][0]=-2/kk*x12*z12;
733 m[2][2]=2-2/kk*z12*z12;
734 d[0]=2*p0[0]-2/kk*p0[0]*x12*x12-2/kk*p0[2]*x12*z12-2/kk*p0[1]*x12*y12;
735 d[1]=2*p0[1]-2/kk*p0[1]*y12*y12-2/kk*p0[0]*x12*y12-2/kk*p0[2]*z12*y12;
736 d[2]=2*p0[2]-2/kk*p0[2]*z12*z12-2/kk*p0[0]*x12*z12-2/kk*p0[1]*z12*y12;
739 //____________________________________________________________________________
740 void AliVertexerTracks::GetStrLinDerivMatrix(Double_t *p0,Double_t *p1,Double_t *sigmasq,Double_t (*m)[3],Double_t *d){
742 Double_t x12=p1[0]-p0[0];
743 Double_t y12=p1[1]-p0[1];
744 Double_t z12=p1[2]-p0[2];
746 Double_t den= x12*x12*sigmasq[1]*sigmasq[2]+y12*y12*sigmasq[0]*sigmasq[2]+z12*z12*sigmasq[0]*sigmasq[1];
748 Double_t kk= 2*(x12*x12/sigmasq[0]+y12*y12/sigmasq[1]+z12*z12/sigmasq[2]);
751 cc[0]=-x12/sigmasq[0];
752 cc[1]=-y12/sigmasq[1];
753 cc[2]=-z12/sigmasq[2];
755 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;
757 Double_t ss= -p0[0]*cc[0]-p0[1]*cc[1]-p0[2]*cc[2];
760 aa[0]=x12*sigmasq[1]*sigmasq[2]/den;
761 aa[1]=y12*sigmasq[0]*sigmasq[2]/den;
762 aa[2]=z12*sigmasq[0]*sigmasq[1]/den;
764 m[0][0]=aa[0]*(aa[0]*kk+2*cc[0])+2*cc[0]*aa[0]+2/sigmasq[0];
765 m[0][1]=aa[1]*(aa[0]*kk+2*cc[0])+2*cc[1]*aa[0];
766 m[0][2]=aa[2]*(aa[0]*kk+2*cc[0])+2*cc[2]*aa[0];
768 m[1][0]=aa[0]*(aa[1]*kk+2*cc[1])+2*cc[0]*aa[1];
769 m[1][1]=aa[1]*(aa[1]*kk+2*cc[1])+2*cc[1]*aa[1]+2/sigmasq[1];
770 m[1][2]=aa[2]*(aa[1]*kk+2*cc[1])+2*cc[2]*aa[1];
772 m[2][0]=aa[0]*(aa[2]*kk+2*cc[2])+2*cc[0]*aa[2];
773 m[2][1]=aa[1]*(aa[2]*kk+2*cc[2])+2*cc[1]*aa[2];
774 m[2][2]=aa[2]*(aa[2]*kk+2*cc[2])+2*cc[2]*aa[2]+2/sigmasq[2];
776 d[0]=-ww*(aa[0]*kk+2*cc[0])-2*ss*aa[0]+2*p0[0]/sigmasq[0];
777 d[1]=-ww*(aa[1]*kk+2*cc[1])-2*ss*aa[1]+2*p0[1]/sigmasq[1];
778 d[2]=-ww*(aa[2]*kk+2*cc[2])-2*ss*aa[2]+2*p0[2]/sigmasq[2];
781 //_____________________________________________________________________________
782 Double_t AliVertexerTracks::GetStrLinMinDist(Double_t *p0,Double_t *p1,Double_t *x0){
784 Double_t x12=p0[0]-p1[0];
785 Double_t y12=p0[1]-p1[1];
786 Double_t z12=p0[2]-p1[2];
787 Double_t x10=p0[0]-x0[0];
788 Double_t y10=p0[1]-x0[1];
789 Double_t z10=p0[2]-x0[2];
790 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);
792 //---------------------------------------------------------------------------
793 void AliVertexerTracks::VertexFitter(Bool_t useNominalVtx) {
795 // The optimal estimate of the vertex position is given by a "weighted
796 // average of tracks positions"
797 // Original method: CMS Note 97/0051
800 fVert.GetXYZ(initPos);
802 printf(" VertexFitter(): start\n");
803 printf(" Number of tracks in array: %d\n",(Int_t)fTrkArray.GetEntriesFast());
804 printf(" Minimum # tracks required in fit: %d\n",fMinTracks);
805 printf("Vertex position after finder: %f,%f,%f\n",initPos[0],initPos[1],initPos[2]);
806 if(useNominalVtx) printf(" This vertex will be used in fit: (%f+-%f,%f+-%f) <- NOT IMPLEMENTED YET!n",fNominalPos[0],fNominalSigma[0],fNominalPos[1],fNominalSigma[1]);
812 rv(0,0) = initPos[0];
813 rv(1,0) = initPos[1];
815 Double_t xlStart,alpha;
817 Double_t cosRot,sinRot;
821 Int_t arrEntries = (Int_t)fTrkArray.GetEntries();
825 Int_t *skipTrack = new Int_t[arrEntries];
826 for(i=0; i<arrEntries; i++) skipTrack[i]=0;
829 // 1st - first estimate of vtx using all tracks
830 // 2nd - apply cut on chi2 max per track
831 // 3rd - estimate of global chi2
832 for(step=0; step<3; step++) {
833 if(fDebug) printf(" step = %d\n",step);
837 TMatrixD sumWiri(3,1);
841 for(j=0; j<3; j++) sumWi(j,i) = 0.;
846 sumWiri(i,0) += fNominalPos[i]/fNominalSigma[i]/fNominalSigma[i];
847 sumWi(i,i) += 1./fNominalSigma[i]/fNominalSigma[i];
853 for(k=0; k<arrEntries; k++) {
854 if(skipTrack[k]) continue;
855 // get track from track array
856 t = (AliESDtrack*)fTrkArray.At(k);
857 alpha = t->GetAlpha();
858 xlStart = initPos[0]*TMath::Cos(alpha)+initPos[1]*TMath::Sin(alpha);
859 t->PropagateTo(xlStart,AliTracker::GetBz()); // to vtxSeed
861 if(alpha<0.) rotAngle += 2.*TMath::Pi();
862 cosRot = TMath::Cos(rotAngle);
863 sinRot = TMath::Sin(rotAngle);
865 // vector of track global coordinates
867 ri(0,0) = t->GetX()*cosRot-t->GetY()*sinRot;
868 ri(1,0) = t->GetX()*sinRot+t->GetY()*cosRot;
871 // matrix to go from global (x,y,z) to local (y,z);
880 // covariance matrix of local (y,z) - inverted
882 t->GetExternalCovariance(cc);
888 // weights matrix: wWi = qQiT * uUiInv * qQi
889 if(uUi.Determinant() <= 0.) continue;
890 TMatrixD uUiInv(TMatrixD::kInverted,uUi);
891 TMatrixD uUiInvQi(uUiInv,TMatrixD::kMult,qQi);
892 TMatrixD wWi(qQi,TMatrixD::kTransposeMult,uUiInvQi);
895 TMatrixD deltar = rv; deltar -= ri;
896 TMatrixD wWideltar(wWi,TMatrixD::kMult,deltar);
897 chi2i = deltar(0,0)*wWideltar(0,0)+
898 deltar(1,0)*wWideltar(1,0)+
899 deltar(2,0)*wWideltar(2,0);
902 if(step==1 && chi2i > fMaxChi2PerTrack) {
910 TMatrixD wWiri(wWi,TMatrixD::kMult,ri);
916 } // end loop on tracks
918 if(nUsedTrks < fMinTracks) {
923 Double_t determinant = sumWi.Determinant();
924 //cerr<<" determinant: "<<determinant<<endl;
925 if(determinant < 100.) {
926 printf("det(V) = 0\n");
931 // inverted of weights matrix
932 TMatrixD invsumWi(TMatrixD::kInverted,sumWi);
935 // position of primary vertex
938 } // end loop on the 3 steps
944 printf("TooFewTracks\n");
945 fCurrentVertex = new AliESDVertex(0.,0.,-1);
949 Double_t position[3];
950 position[0] = rv(0,0);
951 position[1] = rv(1,0);
952 position[2] = rv(2,0);
953 Double_t covmatrix[6];
954 covmatrix[0] = vV(0,0);
955 covmatrix[1] = vV(0,1);
956 covmatrix[2] = vV(1,1);
957 covmatrix[3] = vV(0,2);
958 covmatrix[4] = vV(1,2);
959 covmatrix[5] = vV(2,2);
961 // store data in the vertex object
962 fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nUsedTrks);
965 printf(" VertexFitter(): finish\n");
966 printf(" rv = ( %f , %f , %f )\n\n",rv(0,0),rv(1,0),rv(2,0));
967 fCurrentVertex->PrintStatus();