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():
56 // Default constructor
64 //-----------------------------------------------------------------------------
65 AliVertexerTracks::AliVertexerTracks(Double_t xStart, Double_t yStart):
80 // Standard constructor
82 SetVtxStart(xStart,yStart);
88 //-----------------------------------------------------------------------------
89 AliVertexerTracks::~AliVertexerTracks() {
91 // The objects pointed by the following pointer are not owned
92 // by this class and are not deleted
96 //---------------------------------------------------------------------------
97 void AliVertexerTracks::SetSkipTracks(Int_t n,Int_t *skipped) {
99 // Mark the tracks not ot be used in the vertex finding
102 fTrksToSkip = new Int_t[n];
103 for(Int_t i=0;i<n;i++) fTrksToSkip[i] = skipped[i];
106 //----------------------------------------------------------------------------
107 Int_t AliVertexerTracks::PrepareTracks(TTree &trkTree, Int_t OptImpParCut) {
109 // Propagate tracks to initial vertex position and store them in a TObjArray
111 Double_t maxd0rphi = 3.;
113 Double_t sigmaCurr[3],sigmaVtx[3],posVtx[3];
114 Double_t covVtx[6],covVtxXY;
115 Float_t d0z0[2],covd0z0[3];
116 Double_t momentum[3],postrk[3];
117 Double_t pt,sigma,sigmavtxphi,phitrk;
118 Double_t field=GetField();
120 AliESDVertex *initVertex = new AliESDVertex(fNominalPos,fNominalSigma);
122 Int_t nEntries = (Int_t)trkTree.GetEntries();
123 if(!fTrkArray.IsEmpty()) fTrkArray.Clear();
124 fTrkArray.Expand(nEntries);
127 printf(" PrepareTracks()\n");
130 for(Int_t i=0; i<nEntries; i++) {
131 AliESDtrack *track = new AliESDtrack;
132 trkTree.SetBranchAddress("tracks",&track);
137 // propagate track to vertex
138 if(OptImpParCut==1) { // OptImpParCut==1
139 track->RelateToVertex(initVertex,field,100.);
140 initVertex->GetXYZ(posVtx);
141 initVertex->GetSigmaXYZ(sigmaVtx);
143 } else { // OptImpParCut==2
144 fCurrentVertex->GetSigmaXYZ(sigmaCurr);
145 if((sigmaCurr[0]+sigmaCurr[1])<(fNominalSigma[0]+fNominalSigma[1])) {
146 track->RelateToVertex(fCurrentVertex,field,100.);
147 fCurrentVertex->GetXYZ(posVtx);
148 fCurrentVertex->GetSigmaXYZ(sigmaVtx);
149 fCurrentVertex->GetCovMatrix(covVtx);
150 covVtxXY = covVtx[1];
152 track->RelateToVertex(initVertex,field,100.);
153 initVertex->GetXYZ(posVtx);
154 initVertex->GetSigmaXYZ(sigmaVtx);
159 // select tracks with d0rphi < maxd0rphi
160 track->GetImpactParameters(d0z0,covd0z0);
162 track->GetXYZ(postrk);
163 track->GetPxPyPz(momentum);
164 pt = TMath::Sqrt(momentum[0]*momentum[0]+momentum[1]*momentum[1]);
166 phitrk = TMath::ATan2(postrk[1]-posVtx[1],postrk[0]-posVtx[0]);
167 sigmavtxphi = TMath::Sqrt(sigmaVtx[0]*sigmaVtx[0]*
168 TMath::Cos(phitrk)*TMath::Cos(phitrk)+
169 sigmaVtx[1]*sigmaVtx[1]*
170 TMath::Sin(phitrk)*TMath::Sin(phitrk)+
172 TMath::Cos(phitrk)*TMath::Sin(phitrk));
173 sigma = TMath::Sqrt(Sigmad0rphi(pt)*Sigmad0rphi(pt)+
174 sigmavtxphi*sigmavtxphi);
175 maxd0rphi = fNSigma*sigma;
176 if(OptImpParCut==1) maxd0rphi *= 5.;
178 if(fDebug) printf("trk %d; lab %d; |d0| = %f; cut = %f\n",i,track->GetLabel(),TMath::Abs(d0z0[0]),maxd0rphi);
179 if(TMath::Abs(d0z0[0]) > maxd0rphi) {
180 if(fDebug) printf(" rejected\n");
181 delete track; continue;
184 fTrkArray.AddLast(track);
192 //----------------------------------------------------------------------------
193 Double_t AliVertexerTracks::Sigmad0rphi(Double_t pt) const {
195 // Impact parameter resolution in rphi [cm] vs pt [GeV/c]
197 Double_t a_meas = 11.6;
198 Double_t a_scatt = 65.8;
201 Double_t sigma = a_meas*a_meas+a_scatt*a_scatt/TMath::Power(pt,b);
202 sigma = 0.0001*TMath::Sqrt(sigma);
206 //----------------------------------------------------------------------------
207 AliVertex* AliVertexerTracks::VertexForSelectedTracks(TTree *trkTree) {
209 // Return vertex from tracks in trkTree
212 // get tracks and propagate them to initial vertex position
213 Int_t nTrks = PrepareTracks(*trkTree,1);
214 if(nTrks < fMinTracks) {
215 printf("TooFewTracks\n");
216 Double_t vtx[3]={0,0,0};
218 fVert.SetDispersion(999);
219 fVert.SetNContributors(-5);
223 // Set initial vertex position from ESD
224 if(fAlgo==1) StrLinVertexFinderMinDist(1);
225 if(fAlgo==2) StrLinVertexFinderMinDist(0);
226 if(fAlgo==3) HelixVertexFinder();
227 if(fAlgo==4) VertexFinder(1);
228 if(fAlgo==5) VertexFinder(0);
231 //----------------------------------------------------------------------------
232 AliESDVertex* AliVertexerTracks::FindPrimaryVertex(const AliESD *esdEvent)
235 // Primary vertex for current ESD event
237 // 1st with 5*fNSigma*sigma(pt) cut w.r.t. to initial vertex;
238 // 2nd with fNSigma*sigma(pt) cut w.r.t. to vertex found in 1st iteration)
239 // All ESD tracks with inside the beam pipe are then propagated to found vertex
243 Int_t entr = (Int_t)esdEvent->GetNumberOfTracks();
244 TTree *trkTree = new TTree("TreeT","tracks");
245 AliESDtrack *esdTrack = 0;
246 trkTree->Branch("tracks","AliESDtrack",&esdTrack);
249 for(Int_t i=0; i<entr; i++) {
250 // check tracks to skip
252 for(Int_t j=0; j<fNTrksToSkip; j++) {
253 if(i==fTrksToSkip[j]) {
254 if(fDebug) printf("skipping track: %d\n",i);
258 if(skipThis) continue;
259 AliESDtrack *et = esdEvent->GetTrack(i);
260 esdTrack = new AliESDtrack(*et);
261 if(!esdTrack->GetStatus()&AliESDtrack::kITSin) continue;
262 if(!esdTrack->GetStatus()&AliESDtrack::kITSrefit) continue;
263 Int_t nclus=esdTrack->GetNcls(0); // check number of clusters in ITS
264 if(nclus<fMinITSClusters) continue;
271 // propagate tracks to initVertex
272 // preselect them (reject for |d0|>5*fNSigma*sigma w.r.t. initVertex)
274 nTrks = PrepareTracks(*trkTree,1);
275 if(fDebug) printf(" tracks prepared - step 1: %d\n",nTrks);
276 if(nTrks < fMinTracks) {
277 printf("TooFewTracks\n");
278 TooFewTracks(esdEvent);
279 if(fDebug) fCurrentVertex->PrintStatus();
280 return fCurrentVertex;
285 case 1: StrLinVertexFinderMinDist(1); break;
286 case 2: StrLinVertexFinderMinDist(0); break;
287 case 3: HelixVertexFinder(); break;
288 case 4: VertexFinder(1); break;
289 case 5: VertexFinder(0); break;
290 default: printf("Wrong algorithm\n"); break;
292 if(fDebug) printf(" vertex finding completed\n");
296 if(fDebug) printf(" vertex fit completed\n");
299 // propagate tracks to best between initVertex and fCurrentVertex
300 // preselect tracks (reject for |d0|>fNSigma*sigma w.r.t. best
301 // between initVertex and fCurrentVertex)
302 nTrks = PrepareTracks(*trkTree,2);
304 if(fDebug) printf(" tracks prepared - step 2: %d\n",nTrks);
305 if(nTrks < fMinTracks) {
306 printf("TooFewTracks\n");
307 TooFewTracks(esdEvent);
308 if(fDebug) fCurrentVertex->PrintStatus();
309 return fCurrentVertex;
314 case 1: StrLinVertexFinderMinDist(1); break;
315 case 2: StrLinVertexFinderMinDist(0); break;
316 case 3: HelixVertexFinder(); break;
317 case 4: VertexFinder(1); break;
318 case 5: VertexFinder(0); break;
319 default: printf("Wrong algorithm\n"); break;
321 if(fDebug) printf(" vertex finding completed\n");
325 if(fDebug) printf(" vertex fit completed\n");
330 // take true pos from ESD
332 esdEvent->GetVertex()->GetTruePos(tp);
333 fCurrentVertex->SetTruePos(tp);
334 if(fNominalSigma[0]>1.) {
335 fCurrentVertex->SetTitle("VertexerTracksNoConstraint");
337 fCurrentVertex->SetTitle("VertexerTracksWithConstraint");
340 if(fDebug) fCurrentVertex->PrintStatus();
341 if(fTrksToSkip) delete [] fTrksToSkip;
343 // propagate tracks to found vertex
344 if(TMath::Sqrt(fCurrentVertex->GetXv()*fCurrentVertex->GetXv()+fCurrentVertex->GetYv()*fCurrentVertex->GetYv())<3.) {
345 for(Int_t ii=0; ii<entr; ii++) {
346 AliESDtrack *et = esdEvent->GetTrack(ii);
347 if(!et->GetStatus()&AliESDtrack::kITSin) continue;
348 if(et->GetX()>3.) continue;
349 et->RelateToVertex(fCurrentVertex,GetField(),100.);
352 AliWarning("Found vertex outside beam pipe!");
355 return fCurrentVertex;
357 //----------------------------------------------------------------------------
358 AliESDVertex* AliVertexerTracks::FindPrimaryVertexOld(const AliESD *esdEvent)
361 // Primary vertex for current ESD event (one iteration with |d0rphi|<3cm)
365 Int_t entr = (Int_t)esdEvent->GetNumberOfTracks();
366 TTree *trkTree = new TTree("TreeT","tracks");
367 AliESDtrack *esdTrack = 0;
368 trkTree->Branch("tracks","AliESDtrack",&esdTrack);
370 for(Int_t i=0; i<entr; i++) {
371 AliESDtrack *et = esdEvent->GetTrack(i);
372 esdTrack = new AliESDtrack(*et);
373 if(!esdTrack->GetStatus()&AliESDtrack::kITSin) continue;
374 if(!esdTrack->GetStatus()&AliESDtrack::kITSrefit) continue;
375 Int_t nclus=esdTrack->GetNcls(0); // check number of clusters in ITS
376 if(nclus<fMinITSClusters) continue;
382 // preselect tracks and propagate them to initial vertex position
383 Int_t nTrks = PrepareTracks(*trkTree,1);
385 if(fDebug) printf(" tracks prepared: %d\n",nTrks);
386 if(nTrks < fMinTracks) {
387 printf("TooFewTracks\n");
388 fCurrentVertex = new AliESDVertex(0.,0.,-1);
389 return fCurrentVertex;
394 case 1: StrLinVertexFinderMinDist(1); break;
395 case 2: StrLinVertexFinderMinDist(0); break;
396 case 3: HelixVertexFinder(); break;
397 case 4: VertexFinder(1); break;
398 case 5: VertexFinder(0); break;
399 default: printf("Wrong algorithm\n"); break;
404 if(fDebug) printf(" vertex fit completed\n");
408 esdEvent->GetVertex()->GetTruePos(tp);
409 fCurrentVertex->SetTruePos(tp);
410 fCurrentVertex->SetTitle("VertexerTracks");
413 return fCurrentVertex;
415 //----------------------------------------------------------------------------
416 AliVertex* AliVertexerTracks::VertexForSelectedTracks(TObjArray *trkArray) {
418 // Return vertex from array of tracks
421 // get tracks and propagate them to initial vertex position
422 Int_t nTrks = trkArray->GetEntriesFast();
423 if(nTrks < fMinTracks) {
424 printf("TooFewTracks\n");
425 Double_t vtx[3]={0,0,0};
427 fVert.SetDispersion(999);
428 fVert.SetNContributors(-5);
431 TTree *trkTree = new TTree("TreeT","tracks");
432 AliESDtrack *esdTrack = 0;
433 trkTree->Branch("tracks","AliESDtrack",&esdTrack);
434 for(Int_t i=0; i<nTrks; i++){
435 esdTrack = (AliESDtrack*)trkArray->At(i);
439 AliVertex *vtx = VertexForSelectedTracks(trkTree);
444 //---------------------------------------------------------------------------
445 void AliVertexerTracks::VertexFinder(Int_t optUseWeights) {
447 // Get estimate of vertex position in (x,y) from tracks DCA
451 for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
452 Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
453 Double_t aver[3]={0.,0.,0.};
454 Double_t aversq[3]={0.,0.,0.};
455 Double_t sigmasq[3]={0.,0.,0.};
460 Double_t pos[3],dir[3];
461 Double_t alpha,mindist;
462 Double_t field=GetField();
464 for(Int_t i=0; i<nacc; i++){
465 track1 = (AliESDtrack*)fTrkArray.At(i);
466 alpha=track1->GetAlpha();
467 mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
468 track1->GetXYZAt(mindist,field,pos);
469 track1->GetPxPyPzAt(mindist,field,dir);
470 AliStrLine *line1 = new AliStrLine(pos,dir);
472 // AliStrLine *line1 = new AliStrLine();
473 // track1->ApproximateHelixWithLine(mindist,field,line1);
475 for(Int_t j=i+1; j<nacc; j++){
476 track2 = (AliESDtrack*)fTrkArray.At(j);
477 alpha=track2->GetAlpha();
478 mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
479 track2->GetXYZAt(mindist,field,pos);
480 track2->GetPxPyPzAt(mindist,field,dir);
481 AliStrLine *line2 = new AliStrLine(pos,dir);
482 // AliStrLine *line2 = new AliStrLine();
483 // track2->ApproximateHelixWithLine(mindist,field,line2);
484 Double_t distCA=line2->GetDCA(line1);
485 if(fDCAcut<=0 || (fDCAcut>0&&distCA<fDCAcut)){
486 Double_t pnt1[3],pnt2[3],crosspoint[3];
488 if(optUseWeights<=0){
489 Int_t retcode = line2->Cross(line1,crosspoint);
492 for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
493 for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]);
497 Int_t retcode = line1->CrossPoints(line2,pnt1,pnt2);
499 Double_t alpha, cs, sn;
500 alpha=track1->GetAlpha();
501 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
502 Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2();
503 alpha=track2->GetAlpha();
504 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
505 Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2();
506 Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2();
507 Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
508 Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
509 Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
510 crosspoint[0]=wx1*pnt1[0] + wx2*pnt2[0];
511 crosspoint[1]=wy1*pnt1[1] + wy2*pnt2[1];
512 crosspoint[2]=wz1*pnt1[2] + wz2*pnt2[2];
515 for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
516 for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]);
525 for(Int_t jj=0;jj<3;jj++){
526 initPos[jj] = aver[jj]/ncombi;
528 sigmasq[jj]=aversq[jj]-initPos[jj]*initPos[jj];
531 sigma=TMath::Sqrt(TMath::Abs(sigma));
534 Warning("VertexFinder","Finder did not succed");
537 fVert.SetXYZ(initPos);
538 fVert.SetDispersion(sigma);
539 fVert.SetNContributors(ncombi);
541 //---------------------------------------------------------------------------
542 void AliVertexerTracks::HelixVertexFinder() {
544 // Get estimate of vertex position in (x,y) from tracks DCA
549 for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
550 Double_t field=GetField();
552 Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
554 Double_t aver[3]={0.,0.,0.};
555 Double_t averquad[3]={0.,0.,0.};
556 Double_t sigmaquad[3]={0.,0.,0.};
563 Double_t alpha, cs, sn;
564 Double_t crosspoint[3];
565 for(Int_t i=0; i<nacc; i++){
566 track1 = (AliESDtrack*)fTrkArray.At(i);
569 for(Int_t j=i+1; j<nacc; j++){
570 track2 = (AliESDtrack*)fTrkArray.At(j);
572 distCA=track2->PropagateToDCA(track1,field);
574 if(fDCAcut<=0 ||(fDCAcut>0&&distCA<fDCAcut)){
575 track1->GetExternalParameters(x,par);
576 alpha=track1->GetAlpha();
577 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
578 Double_t x1=x*cs - par[0]*sn;
579 Double_t y1=x*sn + par[0]*cs;
581 Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2();
582 track2->GetExternalParameters(x,par);
583 alpha=track2->GetAlpha();
584 cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
585 Double_t x2=x*cs - par[0]*sn;
586 Double_t y2=x*sn + par[0]*cs;
588 Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2();
589 Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2();
590 Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
591 Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
592 Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
593 crosspoint[0]=wx1*x1 + wx2*x2;
594 crosspoint[1]=wy1*y1 + wy2*y2;
595 crosspoint[2]=wz1*z1 + wz2*z2;
598 for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
599 for(Int_t jj=0;jj<3;jj++)averquad[jj]+=(crosspoint[jj]*crosspoint[jj]);
605 for(Int_t jj=0;jj<3;jj++){
606 initPos[jj] = aver[jj]/ncombi;
607 averquad[jj]/=ncombi;
608 sigmaquad[jj]=averquad[jj]-initPos[jj]*initPos[jj];
609 sigma+=sigmaquad[jj];
611 sigma=TMath::Sqrt(TMath::Abs(sigma));
614 Warning("HelixVertexFinder","Finder did not succed");
617 fVert.SetXYZ(initPos);
618 fVert.SetDispersion(sigma);
619 fVert.SetNContributors(ncombi);
621 //---------------------------------------------------------------------------
622 void AliVertexerTracks::StrLinVertexFinderMinDist(Int_t optUseWeights){
624 // Calculate the point at minimum distance to prepared tracks
629 for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
630 const Int_t knacc = (Int_t)fTrkArray.GetEntriesFast();
631 Double_t field=GetField();
634 Double_t (*vectP0)[3]=new Double_t [knacc][3];
635 Double_t (*vectP1)[3]=new Double_t [knacc][3];
638 Double_t dsum[3]={0,0,0};
639 for(Int_t i=0;i<3;i++)
640 for(Int_t j=0;j<3;j++)sum[i][j]=0;
641 for(Int_t i=0; i<knacc; i++){
642 track1 = (AliESDtrack*)fTrkArray.At(i);
643 Double_t alpha=track1->GetAlpha();
644 Double_t mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
645 Double_t pos[3],dir[3];
646 track1->GetXYZAt(mindist,field,pos);
647 track1->GetPxPyPzAt(mindist,field,dir);
648 AliStrLine *line1 = new AliStrLine(pos,dir);
649 // AliStrLine *line1 = new AliStrLine();
650 // track1->ApproximateHelixWithLine(mindist,field,line1);
652 Double_t p0[3],cd[3];
655 Double_t p1[3]={p0[0]+cd[0],p0[1]+cd[1],p0[2]+cd[2]};
665 if(optUseWeights==0)GetStrLinDerivMatrix(p0,p1,matr,dknow);
666 if(optUseWeights==1){
668 sigmasq[0]=track1->GetSigmaY2();
669 sigmasq[1]=track1->GetSigmaY2();
670 sigmasq[2]=track1->GetSigmaZ2();
671 GetStrLinDerivMatrix(p0,p1,sigmasq,matr,dknow);
674 for(Int_t iii=0;iii<3;iii++){
675 dsum[iii]+=dknow[iii];
676 for(Int_t lj=0;lj<3;lj++) sum[iii][lj]+=matr[iii][lj];
682 Double_t det=GetDeterminant3X3(sum);
685 for(Int_t zz=0;zz<3;zz++){
686 for(Int_t ww=0;ww<3;ww++){
687 for(Int_t kk=0;kk<3;kk++) vett[ww][kk]=sum[ww][kk];
689 for(Int_t kk=0;kk<3;kk++) vett[kk][zz]=dsum[kk];
690 initPos[zz]=GetDeterminant3X3(vett)/det;
694 for(Int_t i=0; i<knacc; i++){
695 Double_t p0[3]={0,0,0},p1[3]={0,0,0};
696 for(Int_t ii=0;ii<3;ii++){
697 p0[ii]=vectP0[i][ii];
698 p1[ii]=vectP1[i][ii];
700 sigma+=GetStrLinMinDist(p0,p1,initPos);
703 sigma=TMath::Sqrt(sigma);
705 Warning("StrLinVertexFinderMinDist","Finder did not succed");
710 fVert.SetXYZ(initPos);
711 fVert.SetDispersion(sigma);
712 fVert.SetNContributors(knacc);
714 //_______________________________________________________________________
715 Double_t AliVertexerTracks::GetDeterminant3X3(Double_t matr[][3]){
717 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];
720 //____________________________________________________________________________
721 void AliVertexerTracks::GetStrLinDerivMatrix(Double_t *p0,Double_t *p1,Double_t (*m)[3],Double_t *d){
724 Double_t x12=p0[0]-p1[0];
725 Double_t y12=p0[1]-p1[1];
726 Double_t z12=p0[2]-p1[2];
727 Double_t kk=x12*x12+y12*y12+z12*z12;
728 m[0][0]=2-2/kk*x12*x12;
729 m[0][1]=-2/kk*x12*y12;
730 m[0][2]=-2/kk*x12*z12;
731 m[1][0]=-2/kk*x12*y12;
732 m[1][1]=2-2/kk*y12*y12;
733 m[1][2]=-2/kk*y12*z12;
734 m[2][0]=-2/kk*x12*z12;
736 m[2][2]=2-2/kk*z12*z12;
737 d[0]=2*p0[0]-2/kk*p0[0]*x12*x12-2/kk*p0[2]*x12*z12-2/kk*p0[1]*x12*y12;
738 d[1]=2*p0[1]-2/kk*p0[1]*y12*y12-2/kk*p0[0]*x12*y12-2/kk*p0[2]*z12*y12;
739 d[2]=2*p0[2]-2/kk*p0[2]*z12*z12-2/kk*p0[0]*x12*z12-2/kk*p0[1]*z12*y12;
742 //____________________________________________________________________________
743 void AliVertexerTracks::GetStrLinDerivMatrix(Double_t *p0,Double_t *p1,Double_t *sigmasq,Double_t (*m)[3],Double_t *d){
745 Double_t x12=p1[0]-p0[0];
746 Double_t y12=p1[1]-p0[1];
747 Double_t z12=p1[2]-p0[2];
749 Double_t den= x12*x12*sigmasq[1]*sigmasq[2]+y12*y12*sigmasq[0]*sigmasq[2]+z12*z12*sigmasq[0]*sigmasq[1];
751 Double_t kk= 2*(x12*x12/sigmasq[0]+y12*y12/sigmasq[1]+z12*z12/sigmasq[2]);
754 cc[0]=-x12/sigmasq[0];
755 cc[1]=-y12/sigmasq[1];
756 cc[2]=-z12/sigmasq[2];
758 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;
760 Double_t ss= -p0[0]*cc[0]-p0[1]*cc[1]-p0[2]*cc[2];
763 aa[0]=x12*sigmasq[1]*sigmasq[2]/den;
764 aa[1]=y12*sigmasq[0]*sigmasq[2]/den;
765 aa[2]=z12*sigmasq[0]*sigmasq[1]/den;
767 m[0][0]=aa[0]*(aa[0]*kk+2*cc[0])+2*cc[0]*aa[0]+2/sigmasq[0];
768 m[0][1]=aa[1]*(aa[0]*kk+2*cc[0])+2*cc[1]*aa[0];
769 m[0][2]=aa[2]*(aa[0]*kk+2*cc[0])+2*cc[2]*aa[0];
771 m[1][0]=aa[0]*(aa[1]*kk+2*cc[1])+2*cc[0]*aa[1];
772 m[1][1]=aa[1]*(aa[1]*kk+2*cc[1])+2*cc[1]*aa[1]+2/sigmasq[1];
773 m[1][2]=aa[2]*(aa[1]*kk+2*cc[1])+2*cc[2]*aa[1];
775 m[2][0]=aa[0]*(aa[2]*kk+2*cc[2])+2*cc[0]*aa[2];
776 m[2][1]=aa[1]*(aa[2]*kk+2*cc[2])+2*cc[1]*aa[2];
777 m[2][2]=aa[2]*(aa[2]*kk+2*cc[2])+2*cc[2]*aa[2]+2/sigmasq[2];
779 d[0]=-ww*(aa[0]*kk+2*cc[0])-2*ss*aa[0]+2*p0[0]/sigmasq[0];
780 d[1]=-ww*(aa[1]*kk+2*cc[1])-2*ss*aa[1]+2*p0[1]/sigmasq[1];
781 d[2]=-ww*(aa[2]*kk+2*cc[2])-2*ss*aa[2]+2*p0[2]/sigmasq[2];
784 //_____________________________________________________________________________
785 Double_t AliVertexerTracks::GetStrLinMinDist(Double_t *p0,Double_t *p1,Double_t *x0){
787 Double_t x12=p0[0]-p1[0];
788 Double_t y12=p0[1]-p1[1];
789 Double_t z12=p0[2]-p1[2];
790 Double_t x10=p0[0]-x0[0];
791 Double_t y10=p0[1]-x0[1];
792 Double_t z10=p0[2]-x0[2];
793 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);
795 //---------------------------------------------------------------------------
796 void AliVertexerTracks::VertexFitter(Bool_t useNominalVtx) {
798 // The optimal estimate of the vertex position is given by a "weighted
799 // average of tracks positions"
800 // Original method: CMS Note 97/0051
803 fVert.GetXYZ(initPos);
805 printf(" VertexFitter(): start\n");
806 printf(" Number of tracks in array: %d\n",(Int_t)fTrkArray.GetEntriesFast());
807 printf(" Minimum # tracks required in fit: %d\n",fMinTracks);
808 printf("Vertex position after finder: %f,%f,%f\n",initPos[0],initPos[1],initPos[2]);
809 if(useNominalVtx) printf(" This vertex will be used in fit: (%f+-%f,%f+-%f)\n",fNominalPos[0],fNominalSigma[0],fNominalPos[1],fNominalSigma[1]);
815 rv(0,0) = initPos[0];
816 rv(1,0) = initPos[1];
818 Double_t xlStart,alpha;
820 Double_t cosRot,sinRot;
824 Int_t arrEntries = (Int_t)fTrkArray.GetEntries();
829 // 1st - estimate of vtx using all tracks
830 // 2nd - estimate of global chi2
831 for(step=0; step<2; step++) {
832 if(fDebug) printf(" step = %d\n",step);
836 TMatrixD sumWiri(3,1);
840 for(j=0; j<3; j++) sumWi(j,i) = 0.;
845 sumWiri(i,0) += fNominalPos[i]/fNominalSigma[i]/fNominalSigma[i];
846 sumWi(i,i) += 1./fNominalSigma[i]/fNominalSigma[i];
852 for(k=0; k<arrEntries; k++) {
853 // get track from track array
854 t = (AliESDtrack*)fTrkArray.At(k);
855 alpha = t->GetAlpha();
856 xlStart = initPos[0]*TMath::Cos(alpha)+initPos[1]*TMath::Sin(alpha);
857 t->AliExternalTrackParam::PropagateTo(xlStart,AliTracker::GetBz()); // to vtxSeed
859 if(alpha<0.) rotAngle += 2.*TMath::Pi();
860 cosRot = TMath::Cos(rotAngle);
861 sinRot = TMath::Sin(rotAngle);
863 // vector of track global coordinates
865 ri(0,0) = t->GetX()*cosRot-t->GetY()*sinRot;
866 ri(1,0) = t->GetX()*sinRot+t->GetY()*cosRot;
869 // matrix to go from global (x,y,z) to local (y,z);
878 // covariance matrix of local (y,z) - inverted
880 t->GetExternalCovariance(cc);
886 // weights matrix: wWi = qQiT * uUiInv * qQi
887 if(uUi.Determinant() <= 0.) continue;
888 TMatrixD uUiInv(TMatrixD::kInverted,uUi);
889 TMatrixD uUiInvQi(uUiInv,TMatrixD::kMult,qQi);
890 TMatrixD wWi(qQi,TMatrixD::kTransposeMult,uUiInvQi);
893 TMatrixD deltar = rv; deltar -= ri;
894 TMatrixD wWideltar(wWi,TMatrixD::kMult,deltar);
895 chi2i = deltar(0,0)*wWideltar(0,0)+
896 deltar(1,0)*wWideltar(1,0)+
897 deltar(2,0)*wWideltar(2,0);
903 TMatrixD wWiri(wWi,TMatrixD::kMult,ri);
909 } // end loop on tracks
911 if(nUsedTrks < fMinTracks) {
916 Double_t determinant = sumWi.Determinant();
917 //cerr<<" determinant: "<<determinant<<endl;
918 if(determinant < 100.) {
919 printf("det(V) = 0\n");
924 // inverted of weights matrix
925 TMatrixD invsumWi(TMatrixD::kInverted,sumWi);
928 // position of primary vertex
931 } // end loop on the 2 steps
936 printf("TooFewTracks\n");
937 fCurrentVertex = new AliESDVertex(0.,0.,-1);
941 Double_t position[3];
942 position[0] = rv(0,0);
943 position[1] = rv(1,0);
944 position[2] = rv(2,0);
945 Double_t covmatrix[6];
946 covmatrix[0] = vV(0,0);
947 covmatrix[1] = vV(0,1);
948 covmatrix[2] = vV(1,1);
949 covmatrix[3] = vV(0,2);
950 covmatrix[4] = vV(1,2);
951 covmatrix[5] = vV(2,2);
953 // store data in the vertex object
954 fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nUsedTrks);
957 printf(" VertexFitter(): finish\n");
958 printf(" rv = ( %f , %f , %f )\n\n",rv(0,0),rv(1,0),rv(2,0));
959 fCurrentVertex->PrintStatus();
964 //---------------------------------------------------------------------------
965 void AliVertexerTracks::TooFewTracks(const AliESD* esdEvent) {
967 // When the number of tracks is < fMinTracks
970 // deal with vertices not found
971 Double_t pos[3],err[3];
973 pos[0] = fNominalPos[0];
974 err[0] = fNominalSigma[0];
975 pos[1] = fNominalPos[1];
976 err[1] = fNominalSigma[1];
977 pos[2] = esdEvent->GetVertex()->GetZv();
978 err[2] = esdEvent->GetVertex()->GetZRes();
979 if(err[0]>1. && esdEvent->GetVertex()->GetNContributors()<=0)
980 ncontr = -1; // (x,y,z) = (0,0,0)
981 if(err[0]>1. && esdEvent->GetVertex()->GetNContributors()>0)
982 ncontr = -2; // (x,y,z) = (0,0,z_fromSPD)
983 if(err[0]<1. && esdEvent->GetVertex()->GetNContributors()<=0)
984 ncontr = -3; // (x,y,z) = (x_mean,y_mean,0)
985 if(err[0]<1. && esdEvent->GetVertex()->GetNContributors()>0)
986 ncontr = -4; // (x,y,z) = (x_mean,y_mean,z_fromSPD)
988 fCurrentVertex = new AliESDVertex(pos,err);
989 fCurrentVertex->SetNContributors(ncontr);
992 esdEvent->GetVertex()->GetTruePos(tp);
993 fCurrentVertex->SetTruePos(tp);
994 fCurrentVertex->SetTitle("VertexerTracksWithConstraint");
995 if(ncontr==-1||ncontr==-2)
996 fCurrentVertex->SetTitle("VertexerTracksNoConstraint");