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 tracks
18 // It accepts V2 and ESD tracks
20 // Origin: A.Dainese, Padova,
21 // andrea.dainese@pd.infn.it
23 // massimo.masera@to.infn.it
24 //-----------------------------------------------------------------
26 //---- standard headers ----
27 #include <Riostream.h>
28 //---- Root headers --------
33 //---- AliRoot headers -----
34 #include "AliITSStrLine.h"
35 #include "AliITStrackV2.h"
36 #include "AliESDVertex.h"
37 #include "AliITSVertexerTracks.h"
39 #include "AliESDtrack.h"
42 ClassImp(AliITSVertexerTracks)
45 //----------------------------------------------------------------------------
46 AliITSVertexerTracks::AliITSVertexerTracks():AliITSVertexer() {
48 // Default constructor
56 for(Int_t i=0; i<3; i++)fInitPos[i] = 0.;
58 //----------------------------------------------------------------------------
59 AliITSVertexerTracks::AliITSVertexerTracks(TFile *inFile,TFile *outFile,
62 Double_t xStart,Double_t yStart) {
64 // Standard constructor
72 SetVtxStart(xStart,yStart);
76 for(Int_t i=0; i<3; i++) fInitPos[i] = 0.;
79 //----------------------------------------------------------------------------
80 AliITSVertexerTracks::AliITSVertexerTracks(Double_t field, TString fn,
81 Double_t xStart,Double_t yStart)
84 // Alternative constructor
89 SetVtxStart(xStart,yStart);
93 for(Int_t i=0; i<3; i++) fInitPos[i] = 0.;
95 //______________________________________________________________________
96 AliITSVertexerTracks::AliITSVertexerTracks(const AliITSVertexerTracks &vtxr) : AliITSVertexer(vtxr) {
98 // Copies are not allowed. The method is protected to avoid misuse.
99 Error("AliITSVertexerTracks","Copy constructor not allowed\n");
102 //______________________________________________________________________
103 AliITSVertexerTracks& AliITSVertexerTracks::operator=(const AliITSVertexerTracks& /* vtxr */){
104 // Assignment operator
105 // Assignment is not allowed. The method is protected to avoid misuse.
106 Error("= operator","Assignment operator not allowed\n");
110 //-----------------------------------------------------------------------------
111 AliITSVertexerTracks::~AliITSVertexerTracks() {
112 // Default Destructor
113 // The objects poited by the following pointers are not owned
114 // by this class and are not deleted
120 //-----------------------------------------------------------------------------
121 Bool_t AliITSVertexerTracks::CheckField() const {
123 // Check if the conv. const. has been set
126 Double_t cc = t.GetConvConst();
127 Double_t field = 100./0.299792458/cc;
129 if(field<0.1 || field>0.6) {
130 printf("AliITSVertexerTracks::CheckField():\n ERROR: AliKalmanTrack::fConvConst not set\n Use AliKalmanTrack::SetConvConst() or AliITSVertexerTracks::SetField()\n");
133 printf("AliITSVertexerTracks::CheckField(): Using B = %3.1f T\n",field);
136 //---------------------------------------------------------------------------
137 void AliITSVertexerTracks::ComputeMaxChi2PerTrack(Int_t nTracks) {
139 // Max. contr. to the chi2 has been tuned as a function of multiplicity
141 if(nTracks < 7) { fMaxChi2PerTrack = 1.e6;
142 } else { fMaxChi2PerTrack = 100.; }
146 //---------------------------------------------------------------------------
147 void AliITSVertexerTracks::FindVertices() {
149 // Vertices for all events from fFirstEvent to fLastEvent
152 // Check if the conv. const. has been set
153 if(!CheckField()) return;
155 TDirectory *curdir = 0;
158 for(Int_t ev=fFirstEvent; ev<=fLastEvent; ev++) {
159 if(ev % 100 == 0 || fDebug) printf("--- Processing event %d of %d ---\n",ev,fLastEvent);
161 FindVertexForCurrentEvent(ev);
163 if(!fCurrentVertex) {
164 printf("AliITSVertexerTracks::FindVertixes(): no tracks tree for event %d\n",ev);
168 if(fDebug) fCurrentVertex->PrintStatus();
170 // write vertex to file
171 TString vtxName = "Vertex_";
173 fCurrentVertex->SetName(vtxName.Data());
174 fCurrentVertex->SetTitle("VertexerTracks");
175 //WriteCurrentVertex();
178 fCurrentVertex->Write();
181 } // loop over events
185 //---------------------------------------------------------------------------
186 void AliITSVertexerTracks::FindVerticesESD() {
188 // Vertices for all events from fFirstEvent to fLastEvent
191 // Check if the conv. const. has been set
192 if(!CheckField()) return;
194 TDirectory *curdir = 0;
198 TIter next(fInFile->GetListOfKeys());
199 // loop on events in file
200 while ((key=(TKey*)next())!=0) {
201 AliESD *esdEvent=(AliESD*)key->ReadObj();
203 printf("AliITSVertexerTracks::FindVerticesESD(): not an ESD!\n");
206 Int_t ev = (Int_t)esdEvent->GetEventNumber();
207 if(ev<fFirstEvent || ev>fLastEvent) { delete esdEvent; continue; }
208 if(ev % 100 == 0 || fDebug)
209 printf("--- Processing event %d of %d ---\n",ev,fLastEvent-fFirstEvent);
211 FindVertexForCurrentEvent(esdEvent);
213 if(!fCurrentVertex) {
214 printf("AliITSVertexerTracks::FindVertixesESD():\n no vertex for event %d\n",ev);
218 cout<<"VERTICE TROVATO\n";
219 if(fDebug) fCurrentVertex->PrintStatus();
221 // write the ESD to file
224 sprintf(ename,"%d",ev);
227 esdEvent->Write(ename,TObject::kOverwrite);
230 } // loop over events
234 //----------------------------------------------------------------------------
235 Int_t AliITSVertexerTracks::PrepareTracks(TTree &trkTree) {
237 // Propagate tracks to initial vertex position and store them in a TObjArray
239 Double_t maxd0rphi = 3.;
240 Double_t alpha,xlStart,d0rphi;
244 Int_t nEntries = (Int_t)trkTree.GetEntries();
246 if(!fTrkArray.IsEmpty()) fTrkArray.Clear();
247 fTrkArray.Expand(nEntries);
250 printf(" PrepareTracks()\n");
254 for(Int_t i=0; i<nEntries; i++) {
255 // check tracks to skip
257 for(Int_t j=0; j<fNTrksToSkip; j++) {
258 if(i==fTrksToSkip[j]) {
259 if(fDebug) printf("skipping track: %d\n",i);
263 if(skipThis) continue;
265 AliITStrackV2 *itstrack = new AliITStrackV2;
266 trkTree.SetBranchAddress("tracks",&itstrack);
270 // propagate track to vtxSeed
271 alpha = itstrack->GetAlpha();
272 xlStart = fNominalPos[0]*TMath::Cos(alpha)+fNominalPos[1]*TMath::Sin(alpha);
273 itstrack->PropagateTo(3.,0.0023,65.19); // to beam pipe (0.8 mm of Be)
274 itstrack->PropagateTo(xlStart,0.,0.); // to vtxSeed
276 // select tracks with d0rphi < maxd0rphi
277 d0rphi = TMath::Abs(itstrack->GetD(fNominalPos[0],fNominalPos[1]));
278 if(d0rphi > maxd0rphi) { delete itstrack; continue; }
280 fTrkArray.AddLast(itstrack);
285 if(fTrksToSkip) delete [] fTrksToSkip;
289 //----------------------------------------------------------------------------
290 void AliITSVertexerTracks::PrintStatus() const {
294 printf(" Initial position (%f,%f)\n",fNominalPos[0],fNominalPos[1]);
295 printf(" Vertex position after vertex finder (%f, %f, %f)\n",fInitPos[0],fInitPos[1],fInitPos[2]);
296 printf(" Number of tracks in array: %d\n",(Int_t)fTrkArray.GetEntriesFast());
297 printf(" Minimum # tracks required in fit: %d\n",fMinTracks);
301 //----------------------------------------------------------------------------
302 AliESDVertex* AliITSVertexerTracks::FindVertexForCurrentEvent(Int_t evnumb) {
304 // Vertex for current event
308 // get tree with tracks from input file
309 TString treeName = "TreeT_ITS_";
311 TTree *trkTree=(TTree*)fInFile->Get(treeName.Data());
312 if(!trkTree) return fCurrentVertex;
315 // get tracks and propagate them to initial vertex position
316 Int_t nTrks = PrepareTracks(*trkTree);
318 if(fDebug) printf(" tracks prepared: %d\n",nTrks);
319 if(nTrks < fMinTracks) { TooFewTracks(); return fCurrentVertex; }
325 ComputeMaxChi2PerTrack(nTrks);
327 if(fDebug) printf(" vertex fit completed\n");
329 return fCurrentVertex;
331 //----------------------------------------------------------------------------
332 AliESDVertex* AliITSVertexerTracks::FindVertexForCurrentEvent(AliESD *esdEvent)
335 // Vertex for current ESD event
338 Double_t vtx[3],cvtx[6];
340 // put tracks reco in ITS in a tree
341 Int_t entr = (Int_t)esdEvent->GetNumberOfTracks();
342 TTree *trkTree = new TTree("TreeT_ITS","its tracks");
343 AliITStrackV2 *itstrack = 0;
344 trkTree->Branch("tracks","AliITStrackV2",&itstrack,entr,0);
346 for(Int_t i=0; i<entr; i++) {
347 AliESDtrack *esdTrack = (AliESDtrack*)esdEvent->GetTrack(i);
348 if(!esdTrack->GetStatus()&AliESDtrack::kITSin)
349 { delete esdTrack; continue; }
350 itstrack = new AliITStrackV2(*esdTrack);
357 // preselect tracks and propagate them to initial vertex position
358 Int_t nTrks = PrepareTracks(*trkTree);
360 if(fDebug) printf(" tracks prepared: %d\n",nTrks);
361 if(nTrks < fMinTracks) { TooFewTracks(); return fCurrentVertex; }
363 // Set initial vertex position from ESD
364 esdEvent->GetVertex()->GetXYZ(vtx);
365 SetVtxStart(vtx[0],vtx[1]);
371 ComputeMaxChi2PerTrack(nTrks);
373 if(fDebug) printf(" vertex fit completed\n");
375 // store vertex information in ESD
376 fCurrentVertex->GetXYZ(vtx);
377 fCurrentVertex->GetCovMatrix(cvtx);
378 esdEvent->SetVertex(fCurrentVertex);
380 cout<<"Vertex: "<<vtx[0]<<", "<<vtx[1]<<", "<<vtx[2]<<endl;
381 return fCurrentVertex;
383 //---------------------------------------------------------------------------
384 void AliITSVertexerTracks::SetSkipTracks(Int_t n,Int_t *skipped) {
386 // Mark the tracks not ot be used in the vertex finding
389 fTrksToSkip = new Int_t[n];
390 for(Int_t i=0;i<n;i++) fTrksToSkip[i] = skipped[i];
393 //---------------------------------------------------------------------------
394 void AliITSVertexerTracks::TooFewTracks() {
396 // When the number of tracks is < fMinTracks the vertex is set to (0,0,0)
397 // and the number of tracks to -1
399 fCurrentVertex = new AliESDVertex(0.,0.,-1);
402 //---------------------------------------------------------------------------
403 void AliITSVertexerTracks::VertexFinder() {
405 // Get estimate of vertex position in (x,y) from tracks DCA
406 // Then this estimate is stored to the data member fInitPos
407 // (previous values are overwritten)
412 ******* TEMPORARY!!! FOR TEST ONLY!!! **********************************
414 fInitPos[0] = fNominalPos[0]+gRandom->Gaus(0.,0.0100); // 100 micron gaussian smearing
415 fInitPos[1] = fNominalPos[1]+gRandom->Gaus(0.,0.0100); // 100 micron gaussian smearing
419 for(Int_t i=0;i<2;i++)fInitPos[i]=fNominalPos[i];
421 Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
423 Double_t aver[3]={0.,0.,0.};
425 AliITStrackV2 *track1;
426 AliITStrackV2 *track2;
427 for(Int_t i=0; i<nacc; i++){
428 track1 = (AliITStrackV2*)fTrkArray.At(i);
431 track1->GetExternalParameters(xv,par);
432 cout<<"Track in position "<<i<<" xr= "<<xv<<endl;
433 for(Int_t ii=0;ii<5;ii++)cout<<par[ii]<<" ";
437 Double_t alpha = track1->GetAlpha();
438 Double_t azim = TMath::ASin(track1->GetSnp())+alpha;
439 Double_t theta = TMath::Pi()/2. - TMath::ATan(track1->GetTgl());
440 mom1[0] = TMath::Sin(theta)*TMath::Cos(azim);
441 mom1[1] = TMath::Sin(theta)*TMath::Sin(azim);
442 mom1[2] = TMath::Cos(theta);
445 Double_t mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
446 track1->GetGlobalXYZat(mindist,pos1[0],pos1[1],pos1[2]);
447 AliITSStrLine *line1 = new AliITSStrLine(pos1,mom1);
448 for(Int_t j=i+1; j<nacc; j++){
449 track2 = (AliITStrackV2*)fTrkArray.At(j);
451 alpha = track2->GetAlpha();
452 azim = TMath::ASin(track2->GetSnp())+alpha;
453 theta = TMath::Pi()/2. - TMath::ATan(track2->GetTgl());
454 mom2[0] = TMath::Sin(theta)*TMath::Cos(azim);
455 mom2[1] = TMath::Sin(theta)*TMath::Sin(azim);
456 mom2[2] = TMath::Cos(theta);
458 mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
459 track2->GetGlobalXYZat(mindist,pos2[0],pos2[1],pos2[2]);
460 AliITSStrLine *line2 = new AliITSStrLine(pos2,mom2);
461 Double_t crosspoint[3];
462 Int_t retcode = line2->Cross(line1,crosspoint);
464 if(fDebug>10)cout<<" i= "<<i<<", j= "<<j<<endl;
465 if(fDebug>10)cout<<"bad intersection\n";
466 line1->PrintStatus();
467 line2->PrintStatus();
471 for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
472 if(fDebug>10)cout<<" i= "<<i<<", j= "<<j<<endl;
473 if(fDebug>10)cout<<"\n Cross point: ";
474 if(fDebug>10)cout<<crosspoint[0]<<" "<<crosspoint[1]<<" "<<crosspoint[2]<<endl;
481 for(Int_t jj=0;jj<3;jj++)fInitPos[jj] = aver[jj]/ncombi;
484 Warning("VertexFinder","Finder did not succed");
488 //************************************************************************
492 //---------------------------------------------------------------------------
493 void AliITSVertexerTracks::VertexFitter() {
495 // The optimal estimate of the vertex position is given by a "weighted
496 // average of tracks positions"
497 // Original method: CMS Note 97/0051
500 printf(" VertexFitter(): start\n");
508 rv(0,0) = fInitPos[0];
509 rv(1,0) = fInitPos[1];
511 Double_t xlStart,alpha;
513 Double_t cosRot,sinRot;
517 Int_t arrEntries = (Int_t)fTrkArray.GetEntries();
518 AliITStrackV2 *t = 0;
521 Int_t *skipTrack = new Int_t[arrEntries];
522 for(i=0; i<arrEntries; i++) skipTrack[i]=0;
525 // 1st - first estimate of vtx using all tracks
526 // 2nd - apply cut on chi2 max per track
527 // 3rd - estimate of global chi2
528 for(step=0; step<3; step++) {
529 if(fDebug) printf(" step = %d\n",step);
533 TMatrixD sumWiri(3,1);
537 for(j=0; j<3; j++) sumWi(j,i) = 0.;
541 for(k=0; k<arrEntries; k++) {
542 if(skipTrack[k]) continue;
543 // get track from track array
544 t = (AliITStrackV2*)fTrkArray.At(k);
545 alpha = t->GetAlpha();
546 xlStart = fInitPos[0]*TMath::Cos(alpha)+fInitPos[1]*TMath::Sin(alpha);
547 t->PropagateTo(xlStart,0.,0.); // to vtxSeed
549 if(alpha<0.) rotAngle += 2.*TMath::Pi();
550 cosRot = TMath::Cos(rotAngle);
551 sinRot = TMath::Sin(rotAngle);
553 // vector of track global coordinates
555 ri(0,0) = t->GetX()*cosRot-t->GetY()*sinRot;
556 ri(1,0) = t->GetX()*sinRot+t->GetY()*cosRot;
559 // matrix to go from global (x,y,z) to local (y,z);
568 // covariance matrix of local (y,z) - inverted
570 t->GetExternalCovariance(cc);
576 // weights matrix: wWi = qQiT * uUiInv * qQi
577 if(uUi.Determinant() <= 0.) continue;
578 TMatrixD uUiInv(TMatrixD::kInverted,uUi);
579 TMatrixD uUiInvQi(uUiInv,TMatrixD::kMult,qQi);
580 TMatrixD wWi(qQi,TMatrixD::kTransposeMult,uUiInvQi);
583 TMatrixD deltar = rv; deltar -= ri;
584 TMatrixD wWideltar(wWi,TMatrixD::kMult,deltar);
585 chi2i = deltar(0,0)*wWideltar(0,0)+
586 deltar(1,0)*wWideltar(1,0)+
587 deltar(2,0)*wWideltar(2,0);
590 if(step==1 && chi2i > fMaxChi2PerTrack) {
598 TMatrixD wWiri(wWi,TMatrixD::kMult,ri);
604 } // end loop on tracks
606 if(nUsedTrks < fMinTracks) {
611 Double_t determinant = sumWi.Determinant();
612 //cerr<<" determinant: "<<determinant<<endl;
613 if(determinant < 100.) {
614 printf("det(V) = 0\n");
619 // inverted of weights matrix
620 TMatrixD invsumWi(TMatrixD::kInverted,sumWi);
623 // position of primary vertex
626 } // end loop on the 3 steps
636 Double_t position[3];
637 position[0] = rv(0,0);
638 position[1] = rv(1,0);
639 position[2] = rv(2,0);
640 Double_t covmatrix[6];
641 covmatrix[0] = vV(0,0);
642 covmatrix[1] = vV(0,1);
643 covmatrix[2] = vV(1,1);
644 covmatrix[3] = vV(0,2);
645 covmatrix[4] = vV(1,2);
646 covmatrix[5] = vV(2,2);
648 // store data in the vertex object
649 fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nUsedTrks);
652 printf(" VertexFitter(): finish\n");
653 printf(" rv = ( %f , %f , %f )\n\n",rv(0,0),rv(1,0),rv(2,0));
654 fCurrentVertex->PrintStatus();
659 //----------------------------------------------------------------------------
660 AliESDVertex *AliITSVertexerTracks::VertexOnTheFly(TTree &trkTree) {
662 // Return vertex from tracks in trkTree
664 if(fCurrentVertex) fCurrentVertex = 0;
666 // get tracks and propagate them to initial vertex position
667 Int_t nTrks = PrepareTracks(*(&trkTree));
668 if(fDebug) printf(" tracks prepared: %d\n",nTrks);
669 if(nTrks < fMinTracks) { TooFewTracks(); return fCurrentVertex; }
675 ComputeMaxChi2PerTrack(nTrks);
677 if(fDebug) printf(" vertex fit completed\n");
679 return fCurrentVertex;
681 //----------------------------------------------------------------------------