/************************************************************************** * Copyright(c) 1998-2003, ALICE Experiment at CERN, All rights reserved. * * * * Author: The ALICE Off-line Project. * * Contributors are mentioned in the code where appropriate. * * * * Permission to use, copy, modify and distribute this software and its * * documentation strictly for non-commercial purposes is hereby granted * * without fee, provided that the above copyright notice appears in all * * copies and that both the copyright notice and this permission notice * * appear in the supporting documentation. The authors make no claims * * about the suitability of this software for any purpose. It is * * provided "as is" without express or implied warranty. * **************************************************************************/ //----------------------------------------------------------------- // Implementation of the vertexer from ESD tracks // // Origin: AliITSVertexerTracks // A.Dainese, Padova, // andrea.dainese@pd.infn.it // M.Masera, Torino, // massimo.masera@to.infn.it // Moved to STEER and adapted to ESD tracks: // F.Prino, Torino, prino@to.infn.it //----------------------------------------------------------------- //---- Root headers -------- #include #include //---- AliRoot headers ----- #include "AliStrLine.h" #include "AliVertexerTracks.h" #include "AliESD.h" #include "AliESDtrack.h" ClassImp(AliVertexerTracks) //---------------------------------------------------------------------------- AliVertexerTracks::AliVertexerTracks(): TObject(), fVert(), fCurrentVertex(0), fMinTracks(2), fMinITSClusters(5), fTrkArray(), fTrksToSkip(0), fNTrksToSkip(0), fDCAcut(0), fAlgo(1), fNSigma(3), fDebug(0) { // // Default constructor // SetVtxStart(); SetVtxStartSigma(); SetMinTracks(); SetMinITSClusters(); SetNSigmad0(); } //----------------------------------------------------------------------------- AliVertexerTracks::AliVertexerTracks(Double_t xStart, Double_t yStart): TObject(), fVert(), fCurrentVertex(0), fMinTracks(2), fMinITSClusters(5), fTrkArray(), fTrksToSkip(0), fNTrksToSkip(0), fDCAcut(0), fAlgo(1), fNSigma(3), fDebug(0) { // // Standard constructor // SetVtxStart(xStart,yStart); SetVtxStartSigma(); SetMinTracks(); SetMinITSClusters(); SetNSigmad0(); } //----------------------------------------------------------------------------- AliVertexerTracks::~AliVertexerTracks() { // Default Destructor // The objects pointed by the following pointer are not owned // by this class and are not deleted fCurrentVertex = 0; delete[] fTrksToSkip; } //--------------------------------------------------------------------------- void AliVertexerTracks::SetSkipTracks(Int_t n,Int_t *skipped) { // // Mark the tracks not ot be used in the vertex finding // fNTrksToSkip = n; fTrksToSkip = new Int_t[n]; for(Int_t i=0;iRelateToVertex(initVertex,field,100.); initVertex->GetXYZ(posVtx); initVertex->GetSigmaXYZ(sigmaVtx); covVtxXY = 0.; } else { // OptImpParCut==2 fCurrentVertex->GetSigmaXYZ(sigmaCurr); if((sigmaCurr[0]+sigmaCurr[1])<(fNominalSigma[0]+fNominalSigma[1])) { track->RelateToVertex(fCurrentVertex,field,100.); fCurrentVertex->GetXYZ(posVtx); fCurrentVertex->GetSigmaXYZ(sigmaVtx); fCurrentVertex->GetCovMatrix(covVtx); covVtxXY = covVtx[1]; } else { track->RelateToVertex(initVertex,field,100.); initVertex->GetXYZ(posVtx); initVertex->GetSigmaXYZ(sigmaVtx); covVtxXY = 0.; } } // select tracks with d0rphi < maxd0rphi track->GetImpactParameters(d0z0,covd0z0); track->GetXYZ(postrk); track->GetPxPyPz(momentum); pt = TMath::Sqrt(momentum[0]*momentum[0]+momentum[1]*momentum[1]); phitrk = TMath::ATan2(postrk[1]-posVtx[1],postrk[0]-posVtx[0]); sigmavtxphi = TMath::Sqrt(sigmaVtx[0]*sigmaVtx[0]* TMath::Cos(phitrk)*TMath::Cos(phitrk)+ sigmaVtx[1]*sigmaVtx[1]* TMath::Sin(phitrk)*TMath::Sin(phitrk)+ covVtxXY* TMath::Cos(phitrk)*TMath::Sin(phitrk)); sigma = TMath::Sqrt(Sigmad0rphi(pt)*Sigmad0rphi(pt)+ sigmavtxphi*sigmavtxphi); maxd0rphi = fNSigma*sigma; if(OptImpParCut==1) maxd0rphi *= 5.; if(fDebug) printf("trk %d; lab %d; |d0| = %f; cut = %f\n",i,track->GetLabel(),TMath::Abs(d0z0[0]),maxd0rphi); if(TMath::Abs(d0z0[0]) > maxd0rphi) { if(fDebug) printf(" rejected\n"); delete track; continue; } fTrkArray.AddLast(track); nTrks++; } delete initVertex; return nTrks; } //---------------------------------------------------------------------------- Double_t AliVertexerTracks::Sigmad0rphi(Double_t pt) const { // // Impact parameter resolution in rphi [cm] vs pt [GeV/c] // Double_t a_meas = 11.6; Double_t a_scatt = 65.8; Double_t b = 1.878; Double_t sigma = a_meas*a_meas+a_scatt*a_scatt/TMath::Power(pt,b); sigma = 0.0001*TMath::Sqrt(sigma); return sigma; } //---------------------------------------------------------------------------- AliVertex* AliVertexerTracks::VertexForSelectedTracks(TTree *trkTree) { // // Return vertex from tracks in trkTree // // get tracks and propagate them to initial vertex position Int_t nTrks = PrepareTracks(*trkTree,1); if(nTrks < fMinTracks) { printf("TooFewTracks\n"); Double_t vtx[3]={0,0,0}; fVert.SetXYZ(vtx); fVert.SetDispersion(999); fVert.SetNContributors(-5); fTrkArray.Delete(); return &fVert; } // Set initial vertex position from ESD if(fAlgo==1) StrLinVertexFinderMinDist(1); if(fAlgo==2) StrLinVertexFinderMinDist(0); if(fAlgo==3) HelixVertexFinder(); if(fAlgo==4) VertexFinder(1); if(fAlgo==5) VertexFinder(0); return &fVert; } //---------------------------------------------------------------------------- AliESDVertex* AliVertexerTracks::FindPrimaryVertex(const AliESD *esdEvent) { // // Primary vertex for current ESD event // (Two iterations: // 1st with 5*fNSigma*sigma(pt) cut w.r.t. to initial vertex; // 2nd with fNSigma*sigma(pt) cut w.r.t. to vertex found in 1st iteration) // All ESD tracks with inside the beam pipe are then propagated to found vertex // fCurrentVertex = 0; Int_t entr = (Int_t)esdEvent->GetNumberOfTracks(); TTree *trkTree = new TTree("TreeT","tracks"); AliESDtrack *esdTrack = 0; trkTree->Branch("tracks","AliESDtrack",&esdTrack); Bool_t skipThis; for(Int_t i=0; iGetTrack(i); esdTrack = new AliESDtrack(*et); if(!(esdTrack->GetStatus()&AliESDtrack::kITSin)) {delete esdTrack;continue;} if(!(esdTrack->GetStatus()&AliESDtrack::kITSrefit)) {delete esdTrack;continue;} Int_t nclus=esdTrack->GetNcls(0); // check number of clusters in ITS if(nclusFill(); delete esdTrack; } // delete esdTrack; // ITERATION 1 // propagate tracks to initVertex // preselect them (reject for |d0|>5*fNSigma*sigma w.r.t. initVertex) Int_t nTrks; nTrks = PrepareTracks(*trkTree,1); if(fDebug) printf(" tracks prepared - step 1: %d\n",nTrks); if(nTrks < fMinTracks) { printf("TooFewTracks\n"); TooFewTracks(esdEvent); if(fDebug) fCurrentVertex->PrintStatus(); fTrkArray.Delete(); delete trkTree; return fCurrentVertex; } // vertex finder switch (fAlgo) { case 1: StrLinVertexFinderMinDist(1); break; case 2: StrLinVertexFinderMinDist(0); break; case 3: HelixVertexFinder(); break; case 4: VertexFinder(1); break; case 5: VertexFinder(0); break; default: printf("Wrong algorithm\n"); break; } if(fDebug) printf(" vertex finding completed\n"); // vertex fitter VertexFitter(kTRUE); if(fDebug) printf(" vertex fit completed\n"); fTrkArray.Delete(); // ITERATION 2 // propagate tracks to best between initVertex and fCurrentVertex // preselect tracks (reject for |d0|>fNSigma*sigma w.r.t. best // between initVertex and fCurrentVertex) nTrks = PrepareTracks(*trkTree,2); delete trkTree; if(fDebug) printf(" tracks prepared - step 2: %d\n",nTrks); if(nTrks < fMinTracks) { printf("TooFewTracks\n"); TooFewTracks(esdEvent); if(fDebug) fCurrentVertex->PrintStatus(); fTrkArray.Delete(); return fCurrentVertex; } // vertex finder switch (fAlgo) { case 1: StrLinVertexFinderMinDist(1); break; case 2: StrLinVertexFinderMinDist(0); break; case 3: HelixVertexFinder(); break; case 4: VertexFinder(1); break; case 5: VertexFinder(0); break; default: printf("Wrong algorithm\n"); break; } if(fDebug) printf(" vertex finding completed\n"); // fitter VertexFitter(kTRUE); if(fDebug) printf(" vertex fit completed\n"); fTrkArray.Delete(); // take true pos from ESD Double_t tp[3]; esdEvent->GetVertex()->GetTruePos(tp); fCurrentVertex->SetTruePos(tp); if(fNominalSigma[0]>1.) { fCurrentVertex->SetTitle("VertexerTracksNoConstraint"); } else { fCurrentVertex->SetTitle("VertexerTracksWithConstraint"); } if(fDebug) fCurrentVertex->PrintStatus(); if(fTrksToSkip) delete [] fTrksToSkip; // propagate tracks to found vertex if(TMath::Sqrt(fCurrentVertex->GetXv()*fCurrentVertex->GetXv()+fCurrentVertex->GetYv()*fCurrentVertex->GetYv())<3.) { for(Int_t ii=0; iiGetTrack(ii); if(!et->GetStatus()&AliESDtrack::kITSin) continue; if(et->GetX()>3.) continue; et->RelateToVertex(fCurrentVertex,GetField(),100.); } } else { AliWarning("Found vertex outside beam pipe!"); } return fCurrentVertex; } //---------------------------------------------------------------------------- AliESDVertex* AliVertexerTracks::FindPrimaryVertexOld(const AliESD *esdEvent) { // // Primary vertex for current ESD event (one iteration with |d0rphi|<3cm) // fCurrentVertex = 0; Int_t entr = (Int_t)esdEvent->GetNumberOfTracks(); TTree *trkTree = new TTree("TreeT","tracks"); AliESDtrack *esdTrack = 0; trkTree->Branch("tracks","AliESDtrack",&esdTrack); for(Int_t i=0; iGetTrack(i); esdTrack = new AliESDtrack(*et); if(!(esdTrack->GetStatus()&AliESDtrack::kITSin)) {delete esdTrack;continue;} if(!(esdTrack->GetStatus()&AliESDtrack::kITSrefit)) {delete esdTrack;continue;} Int_t nclus=esdTrack->GetNcls(0); // check number of clusters in ITS if(nclusFill(); delete esdTrack; } // delete esdTrack; // preselect tracks and propagate them to initial vertex position Int_t nTrks = PrepareTracks(*trkTree,1); delete trkTree; if(fDebug) printf(" tracks prepared: %d\n",nTrks); if(nTrks < fMinTracks) { printf("TooFewTracks\n"); fCurrentVertex = new AliESDVertex(0.,0.,-1); fTrkArray.Delete(); return fCurrentVertex; } // VERTEX FINDER switch (fAlgo) { case 1: StrLinVertexFinderMinDist(1); break; case 2: StrLinVertexFinderMinDist(0); break; case 3: HelixVertexFinder(); break; case 4: VertexFinder(1); break; case 5: VertexFinder(0); break; default: printf("Wrong algorithm\n"); break; } // VERTEX FITTER VertexFitter(); if(fDebug) printf(" vertex fit completed\n"); Double_t tp[3]; esdEvent->GetVertex()->GetTruePos(tp); fCurrentVertex->SetTruePos(tp); fCurrentVertex->SetTitle("VertexerTracks"); fTrkArray.Delete(); return fCurrentVertex; } //---------------------------------------------------------------------------- AliVertex* AliVertexerTracks::VertexForSelectedTracks(TObjArray *trkArray) { // // Return vertex from array of tracks // // get tracks and propagate them to initial vertex position Int_t nTrks = trkArray->GetEntriesFast(); if(nTrks < fMinTracks) { printf("TooFewTracks\n"); Double_t vtx[3]={0,0,0}; fVert.SetXYZ(vtx); fVert.SetDispersion(999); fVert.SetNContributors(-5); fTrkArray.Delete(); return &fVert; } TTree *trkTree = new TTree("TreeT","tracks"); AliESDtrack *esdTrack = 0; trkTree->Branch("tracks","AliESDtrack",&esdTrack); for(Int_t i=0; iAt(i); trkTree->Fill(); } AliVertex *vtx = VertexForSelectedTracks(trkTree); delete trkTree; return vtx; } //--------------------------------------------------------------------------- void AliVertexerTracks::VertexFinder(Int_t optUseWeights) { // Get estimate of vertex position in (x,y) from tracks DCA Double_t initPos[3]; initPos[2] = 0.; for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i]; Int_t nacc = (Int_t)fTrkArray.GetEntriesFast(); Double_t aver[3]={0.,0.,0.}; Double_t aversq[3]={0.,0.,0.}; Double_t sigmasq[3]={0.,0.,0.}; Double_t sigma=0; Int_t ncombi = 0; AliESDtrack *track1; AliESDtrack *track2; Double_t pos[3],dir[3]; Double_t alpha,mindist; Double_t field=GetField(); for(Int_t i=0; iGetAlpha(); mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1]; track1->GetXYZAt(mindist,field,pos); track1->GetPxPyPzAt(mindist,field,dir); AliStrLine *line1 = new AliStrLine(pos,dir); // AliStrLine *line1 = new AliStrLine(); // track1->ApproximateHelixWithLine(mindist,field,line1); for(Int_t j=i+1; jGetAlpha(); mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1]; track2->GetXYZAt(mindist,field,pos); track2->GetPxPyPzAt(mindist,field,dir); AliStrLine *line2 = new AliStrLine(pos,dir); // AliStrLine *line2 = new AliStrLine(); // track2->ApproximateHelixWithLine(mindist,field,line2); Double_t distCA=line2->GetDCA(line1); if(fDCAcut<=0 || (fDCAcut>0&&distCACross(line1,crosspoint); if(retcode>=0){ ncombi++; for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj]; for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]); } } if(optUseWeights>0){ Int_t retcode = line1->CrossPoints(line2,pnt1,pnt2); if(retcode>=0){ Double_t alpha, cs, sn; alpha=track1->GetAlpha(); cs=TMath::Cos(alpha); sn=TMath::Sin(alpha); Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2(); alpha=track2->GetAlpha(); cs=TMath::Cos(alpha); sn=TMath::Sin(alpha); Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2(); Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2(); Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1; Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1; Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1; crosspoint[0]=wx1*pnt1[0] + wx2*pnt2[0]; crosspoint[1]=wy1*pnt1[1] + wy2*pnt2[1]; crosspoint[2]=wz1*pnt1[2] + wz2*pnt2[2]; ncombi++; for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj]; for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]); } } } delete line2; } delete line1; } if(ncombi>0){ for(Int_t jj=0;jj<3;jj++){ initPos[jj] = aver[jj]/ncombi; aversq[jj]/=ncombi; sigmasq[jj]=aversq[jj]-initPos[jj]*initPos[jj]; sigma+=sigmasq[jj]; } sigma=TMath::Sqrt(TMath::Abs(sigma)); } else { Warning("VertexFinder","Finder did not succed"); sigma=999; } fVert.SetXYZ(initPos); fVert.SetDispersion(sigma); fVert.SetNContributors(ncombi); } //--------------------------------------------------------------------------- void AliVertexerTracks::HelixVertexFinder() { // Get estimate of vertex position in (x,y) from tracks DCA Double_t initPos[3]; initPos[2] = 0.; for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i]; Double_t field=GetField(); Int_t nacc = (Int_t)fTrkArray.GetEntriesFast(); Double_t aver[3]={0.,0.,0.}; Double_t averquad[3]={0.,0.,0.}; Double_t sigmaquad[3]={0.,0.,0.}; Double_t sigma=0; Int_t ncombi = 0; AliESDtrack *track1; AliESDtrack *track2; Double_t distCA; Double_t x, par[5]; Double_t alpha, cs, sn; Double_t crosspoint[3]; for(Int_t i=0; iPropagateToDCA(track1,field); if(fDCAcut<=0 ||(fDCAcut>0&&distCAGetExternalParameters(x,par); alpha=track1->GetAlpha(); cs=TMath::Cos(alpha); sn=TMath::Sin(alpha); Double_t x1=x*cs - par[0]*sn; Double_t y1=x*sn + par[0]*cs; Double_t z1=par[1]; Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2(); track2->GetExternalParameters(x,par); alpha=track2->GetAlpha(); cs=TMath::Cos(alpha); sn=TMath::Sin(alpha); Double_t x2=x*cs - par[0]*sn; Double_t y2=x*sn + par[0]*cs; Double_t z2=par[1]; Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2(); Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2(); Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1; Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1; Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1; crosspoint[0]=wx1*x1 + wx2*x2; crosspoint[1]=wy1*y1 + wy2*y2; crosspoint[2]=wz1*z1 + wz2*z2; ncombi++; for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj]; for(Int_t jj=0;jj<3;jj++)averquad[jj]+=(crosspoint[jj]*crosspoint[jj]); } } } if(ncombi>0){ for(Int_t jj=0;jj<3;jj++){ initPos[jj] = aver[jj]/ncombi; averquad[jj]/=ncombi; sigmaquad[jj]=averquad[jj]-initPos[jj]*initPos[jj]; sigma+=sigmaquad[jj]; } sigma=TMath::Sqrt(TMath::Abs(sigma)); } else { Warning("HelixVertexFinder","Finder did not succed"); sigma=999; } fVert.SetXYZ(initPos); fVert.SetDispersion(sigma); fVert.SetNContributors(ncombi); } //--------------------------------------------------------------------------- void AliVertexerTracks::StrLinVertexFinderMinDist(Int_t optUseWeights){ // Calculate the point at minimum distance to prepared tracks Double_t initPos[3]; initPos[2] = 0.; Double_t sigma=0; for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i]; const Int_t knacc = (Int_t)fTrkArray.GetEntriesFast(); Double_t field=GetField(); AliESDtrack *track1; Double_t (*vectP0)[3]=new Double_t [knacc][3]; Double_t (*vectP1)[3]=new Double_t [knacc][3]; Double_t sum[3][3]; Double_t dsum[3]={0,0,0}; for(Int_t i=0;i<3;i++) for(Int_t j=0;j<3;j++)sum[i][j]=0; for(Int_t i=0; iGetAlpha(); Double_t mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1]; Double_t pos[3],dir[3]; track1->GetXYZAt(mindist,field,pos); track1->GetPxPyPzAt(mindist,field,dir); AliStrLine *line1 = new AliStrLine(pos,dir); // AliStrLine *line1 = new AliStrLine(); // track1->ApproximateHelixWithLine(mindist,field,line1); Double_t p0[3],cd[3]; line1->GetP0(p0); line1->GetCd(cd); Double_t p1[3]={p0[0]+cd[0],p0[1]+cd[1],p0[2]+cd[2]}; vectP0[i][0]=p0[0]; vectP0[i][1]=p0[1]; vectP0[i][2]=p0[2]; vectP1[i][0]=p1[0]; vectP1[i][1]=p1[1]; vectP1[i][2]=p1[2]; Double_t matr[3][3]; Double_t dknow[3]; if(optUseWeights==0)GetStrLinDerivMatrix(p0,p1,matr,dknow); if(optUseWeights==1){ Double_t sigmasq[3]; sigmasq[0]=track1->GetSigmaY2(); sigmasq[1]=track1->GetSigmaY2(); sigmasq[2]=track1->GetSigmaZ2(); GetStrLinDerivMatrix(p0,p1,sigmasq,matr,dknow); } for(Int_t iii=0;iii<3;iii++){ dsum[iii]+=dknow[iii]; for(Int_t lj=0;lj<3;lj++) sum[iii][lj]+=matr[iii][lj]; } delete line1; } Double_t vett[3][3]; Double_t det=GetDeterminant3X3(sum); if(det!=0){ for(Int_t zz=0;zz<3;zz++){ for(Int_t ww=0;ww<3;ww++){ for(Int_t kk=0;kk<3;kk++) vett[ww][kk]=sum[ww][kk]; } for(Int_t kk=0;kk<3;kk++) vett[kk][zz]=dsum[kk]; initPos[zz]=GetDeterminant3X3(vett)/det; } for(Int_t i=0; iGetAlpha(); xlStart = initPos[0]*TMath::Cos(alpha)+initPos[1]*TMath::Sin(alpha); t->AliExternalTrackParam::PropagateTo(xlStart,AliTracker::GetBz()); // to vtxSeed rotAngle = alpha; if(alpha<0.) rotAngle += 2.*TMath::Pi(); cosRot = TMath::Cos(rotAngle); sinRot = TMath::Sin(rotAngle); // vector of track global coordinates TMatrixD ri(3,1); ri(0,0) = t->GetX()*cosRot-t->GetY()*sinRot; ri(1,0) = t->GetX()*sinRot+t->GetY()*cosRot; ri(2,0) = t->GetZ(); // matrix to go from global (x,y,z) to local (y,z); TMatrixD qQi(2,3); qQi(0,0) = -sinRot; qQi(0,1) = cosRot; qQi(0,2) = 0.; qQi(1,0) = 0.; qQi(1,1) = 0.; qQi(1,2) = 1.; // covariance matrix of local (y,z) - inverted TMatrixD uUi(2,2); t->GetExternalCovariance(cc); uUi(0,0) = cc[0]; uUi(0,1) = cc[1]; uUi(1,0) = cc[1]; uUi(1,1) = cc[2]; // weights matrix: wWi = qQiT * uUiInv * qQi if(uUi.Determinant() <= 0.) continue; TMatrixD uUiInv(TMatrixD::kInverted,uUi); TMatrixD uUiInvQi(uUiInv,TMatrixD::kMult,qQi); TMatrixD wWi(qQi,TMatrixD::kTransposeMult,uUiInvQi); // track chi2 TMatrixD deltar = rv; deltar -= ri; TMatrixD wWideltar(wWi,TMatrixD::kMult,deltar); chi2i = deltar(0,0)*wWideltar(0,0)+ deltar(1,0)*wWideltar(1,0)+ deltar(2,0)*wWideltar(2,0); // add to total chi2 chi2 += chi2i; TMatrixD wWiri(wWi,TMatrixD::kMult,ri); sumWiri += wWiri; sumWi += wWi; nUsedTrks++; } // end loop on tracks if(nUsedTrks < fMinTracks) { failed=1; continue; } Double_t determinant = sumWi.Determinant(); //cerr<<" determinant: "<PrintStatus(); } return; } //--------------------------------------------------------------------------- void AliVertexerTracks::TooFewTracks(const AliESD* esdEvent) { // // When the number of tracks is < fMinTracks // // deal with vertices not found Double_t pos[3],err[3]; Int_t ncontr=0; pos[0] = fNominalPos[0]; err[0] = fNominalSigma[0]; pos[1] = fNominalPos[1]; err[1] = fNominalSigma[1]; pos[2] = esdEvent->GetVertex()->GetZv(); err[2] = esdEvent->GetVertex()->GetZRes(); if(err[0]>1. && esdEvent->GetVertex()->GetNContributors()<=0) ncontr = -1; // (x,y,z) = (0,0,0) if(err[0]>1. && esdEvent->GetVertex()->GetNContributors()>0) ncontr = -2; // (x,y,z) = (0,0,z_fromSPD) if(err[0]<1. && esdEvent->GetVertex()->GetNContributors()<=0) ncontr = -3; // (x,y,z) = (x_mean,y_mean,0) if(err[0]<1. && esdEvent->GetVertex()->GetNContributors()>0) ncontr = -4; // (x,y,z) = (x_mean,y_mean,z_fromSPD) fCurrentVertex = 0; fCurrentVertex = new AliESDVertex(pos,err); fCurrentVertex->SetNContributors(ncontr); Double_t tp[3]; esdEvent->GetVertex()->GetTruePos(tp); fCurrentVertex->SetTruePos(tp); fCurrentVertex->SetTitle("VertexerTracksWithConstraint"); if(ncontr==-1||ncontr==-2) fCurrentVertex->SetTitle("VertexerTracksNoConstraint"); return; }