]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliVertexerTracks.cxx
Bit added to indicate if this track was used to fit the primary vertex.
[u/mrichter/AliRoot.git] / STEER / AliVertexerTracks.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2003, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 //-----------------------------------------------------------------
17 //    Implementation of the vertexer from ESD tracks
18 //
19 // Origin: AliITSVertexerTracks
20 //         A.Dainese, Padova, 
21 //         andrea.dainese@pd.infn.it
22 //         M.Masera,  Torino, 
23 //         massimo.masera@to.infn.it
24 // Moved to STEER and adapted to ESD tracks: 
25 //          F.Prino,  Torino, prino@to.infn.it
26 //-----------------------------------------------------------------
27
28 //---- Root headers --------
29 #include <TTree.h>
30 #include <TMatrixD.h>
31 //---- AliRoot headers -----
32 #include "AliStrLine.h"
33 #include "AliVertexerTracks.h"
34 #include "AliESD.h"
35 #include "AliESDtrack.h"
36
37 ClassImp(AliVertexerTracks)
38
39
40 //----------------------------------------------------------------------------
41 AliVertexerTracks::AliVertexerTracks():
42 TObject(),
43 fVert(),
44 fCurrentVertex(0),
45 fConstraint(kFALSE),
46 fOnlyFitter(kFALSE),
47 fMinTracks(1),
48 fMinITSClusters(5),
49 fTrkArray(),
50 fTrksToSkip(0),
51 fNTrksToSkip(0),
52 fDCAcut(0),
53 fAlgo(1),
54 fNSigma(3),
55 fMaxd0z0(0.5),
56 fITSin(kTRUE),
57 fITSrefit(kTRUE),
58 fDebug(0)
59 {
60 //
61 // Default constructor
62 //
63   SetVtxStart();
64   SetVtxStartSigma();
65   SetMinTracks();
66   SetMinITSClusters();
67   SetNSigmad0(); 
68   SetMaxd0z0(); 
69 }
70 //-----------------------------------------------------------------------------
71 AliVertexerTracks::AliVertexerTracks(Double_t xStart, Double_t yStart):
72 TObject(),
73 fVert(),
74 fCurrentVertex(0),
75 fConstraint(kFALSE),
76 fOnlyFitter(kFALSE),
77 fMinTracks(1),
78 fMinITSClusters(5),
79 fTrkArray(),
80 fTrksToSkip(0),
81 fNTrksToSkip(0),
82 fDCAcut(0),
83 fAlgo(1),
84 fNSigma(3),
85 fMaxd0z0(0.5),
86 fITSin(kTRUE),
87 fITSrefit(kTRUE),
88 fDebug(0)
89 {
90 //
91 // Standard constructor
92 //
93   SetVtxStart(xStart,yStart);
94   SetVtxStartSigma();
95   SetMinTracks();
96   SetMinITSClusters();
97   SetNSigmad0(); 
98   SetMaxd0z0(); 
99 }
100 //-----------------------------------------------------------------------------
101 AliVertexerTracks::~AliVertexerTracks() 
102 {
103   // Default Destructor
104   // The objects pointed by the following pointer are not owned
105   // by this class and are not deleted
106   fCurrentVertex = 0;
107   delete[] fTrksToSkip;
108 }
109 //----------------------------------------------------------------------------
110 AliESDVertex* AliVertexerTracks::FindPrimaryVertex(const AliESD *esdEvent)
111 {
112 //
113 // Primary vertex for current ESD event
114 // (Two iterations: 
115 //  1st with 5*fNSigma*sigma cut w.r.t. to initial vertex
116 //      + cut on sqrt(d0d0+z0z0) if fConstraint=kFALSE  
117 //  2nd with fNSigma*sigma cut w.r.t. to vertex found in 1st iteration) 
118 // All ESD tracks with inside the beam pipe are then propagated to found vertex
119 //
120   fCurrentVertex = 0;
121
122   // read tracks from ESD
123   Int_t nTrksTot = (Int_t)esdEvent->GetNumberOfTracks();
124   if (nTrksTot<=0) {
125     if(fDebug) printf("TooFewTracks\n");
126     TooFewTracks(esdEvent);
127     return fCurrentVertex;
128   } 
129
130   TTree *trkTree = new TTree("TreeT","tracks");
131   AliESDtrack *esdTrack = 0;
132   trkTree->Branch("tracks","AliESDtrack",&esdTrack);
133
134   Bool_t   skipThis;
135   for(Int_t i=0; i<nTrksTot; i++) {
136     AliESDtrack *et = esdEvent->GetTrack(i);
137     esdTrack = new AliESDtrack(*et);
138     // check tracks to skip
139     skipThis = kFALSE;
140     for(Int_t j=0; j<fNTrksToSkip; j++) { 
141       if(et->GetID()==fTrksToSkip[j]) {
142         if(fDebug) printf("skipping track: %d\n",i);
143         skipThis = kTRUE;
144       }
145     }
146     if(skipThis) {delete esdTrack;continue;}
147     if(fITSin) {
148       if(!(esdTrack->GetStatus()&AliESDtrack::kITSin)) {delete esdTrack;continue;}
149       if(fITSrefit && !(esdTrack->GetStatus()&AliESDtrack::kITSrefit)) {delete esdTrack;continue;}
150       Int_t nclus=esdTrack->GetNcls(0); // check number of clusters in ITS
151       if(nclus<fMinITSClusters) {delete esdTrack;continue;}
152     }
153     trkTree->Fill();
154     delete esdTrack;
155   }
156   
157   // If fConstraint=kFALSE
158   // run VertexFinder(1) to get rough estimate of initVertex (x,y)
159   if(!fConstraint) {
160     // fill fTrkArray, for VertexFinder()
161     if(!fTrkArray.IsEmpty()) fTrkArray.Delete();
162     PrepareTracks(*trkTree,0);
163     Double_t cutsave = fDCAcut;  fDCAcut = 0.1; // 1 mm
164     VertexFinder(1); // using weights, cutting dca < fDCAcut
165     fDCAcut = cutsave;
166     fTrkArray.Delete();
167     if(fVert.GetNContributors()>0) {
168       fVert.GetXYZ(fNominalPos);
169       fNominalPos[0] = fVert.GetXv();
170       fNominalPos[1] = fVert.GetYv();
171       fNominalPos[2] = fVert.GetZv();
172       if(fDebug) printf("No mean vertex: VertexFinder gives (%f, %f, %f)\n",fNominalPos[0],fNominalPos[1],fNominalPos[2]);
173     } else {
174       fNominalPos[0] = 0.;
175       fNominalPos[1] = 0.;
176       fNominalPos[2] = 0.;
177       if(fDebug) printf("No mean vertex and VertexFinder failed\n");
178     }
179   }
180   
181   // TWO ITERATIONS:
182   //
183   // ITERATION 1
184   // propagate tracks to fNominalPos vertex
185   // preselect them:
186   // if(constraint) reject for |d0|>5*fNSigma*sigma w.r.t. fNominal... vertex
187   // else  reject for |d0|\oplus|z0| > 5 mm w.r.t. fNominal... vertex
188   // ITERATION 2
189   // propagate tracks to best between initVertex and fCurrentVertex
190   // preselect tracks (reject for |d0|>fNSigma*sigma w.r.t. best 
191   //                   between initVertex and fCurrentVertex) 
192   for(Int_t iter=0; iter<2; iter++) {
193     if(fOnlyFitter && iter==0) continue; 
194     Int_t nTrksPrep = PrepareTracks(*trkTree,iter+1);
195     if(fDebug) printf(" tracks prepared - iteration %d: %d\n",iter+1,nTrksPrep);
196     if(nTrksPrep < fMinTracks) {
197       if(fDebug) printf("TooFewTracks\n");
198       TooFewTracks(esdEvent);
199       if(fDebug) fCurrentVertex->PrintStatus();
200       fTrkArray.Delete();
201       delete trkTree;
202       return fCurrentVertex; 
203     }
204
205     // vertex finder
206     if(!fOnlyFitter) {
207       if(nTrksPrep==1){
208         if(fDebug) printf("Just one track\n");
209         OneTrackVertFinder();
210       }else{
211         switch (fAlgo) {
212         case 1: StrLinVertexFinderMinDist(1); break;
213         case 2: StrLinVertexFinderMinDist(0); break;
214         case 3: HelixVertexFinder();          break;
215         case 4: VertexFinder(1);              break;
216         case 5: VertexFinder(0);              break;
217         default: printf("Wrong algorithm\n"); break;  
218         }
219       }
220       if(fDebug) printf(" Vertex finding completed\n");
221     }
222
223     // vertex fitter
224     VertexFitter(fConstraint);
225     if(fDebug) printf(" Vertex fit completed\n");
226     if(iter==0) fTrkArray.Delete();
227   } // end loop on the two iterations
228
229
230   // take true pos from SPD vertex in ESD and write it in tracks' vertex
231   Double_t tp[3];
232   esdEvent->GetVertex()->GetTruePos(tp);
233   fCurrentVertex->SetTruePos(tp);
234
235   if(fConstraint) {
236     if(fOnlyFitter) {
237       fCurrentVertex->SetTitle("VertexerTracksWithConstraintOnlyFitter");
238     } else {
239       fCurrentVertex->SetTitle("VertexerTracksWithConstraint");
240     }
241   } else {
242     fCurrentVertex->SetTitle("VertexerTracksNoConstraint");
243   }
244
245     
246   // set indices of used tracks
247   UShort_t *indices = 0;
248   AliESDtrack *ett = 0;
249   if(fCurrentVertex->GetNContributors()>0) {
250     indices = new UShort_t[fCurrentVertex->GetNContributors()];
251     for(Int_t jj=0;jj<(Int_t)fTrkArray.GetEntriesFast();jj++) {
252       ett = (AliESDtrack*)fTrkArray.At(jj);
253       indices[jj] = (UShort_t)ett->GetID();
254     }
255     fCurrentVertex->SetIndices(fCurrentVertex->GetNContributors(),indices);
256   }
257   delete [] indices;
258
259   delete trkTree;
260   
261   fTrkArray.Delete();
262
263   if(fTrksToSkip) delete [] fTrksToSkip;
264
265
266   if(fDebug) fCurrentVertex->PrintStatus();
267   if(fDebug) fCurrentVertex->PrintIndices();
268
269   
270   return fCurrentVertex;
271 }
272 //------------------------------------------------------------------------
273 Double_t AliVertexerTracks::GetDeterminant3X3(Double_t matr[][3])
274 {
275   //
276   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];
277  return det;
278 }
279 //-------------------------------------------------------------------------
280 void AliVertexerTracks::GetStrLinDerivMatrix(Double_t *p0,Double_t *p1,Double_t (*m)[3],Double_t *d)
281 {
282   //
283   Double_t x12=p0[0]-p1[0];
284   Double_t y12=p0[1]-p1[1];
285   Double_t z12=p0[2]-p1[2];
286   Double_t kk=x12*x12+y12*y12+z12*z12;
287   m[0][0]=2-2/kk*x12*x12;
288   m[0][1]=-2/kk*x12*y12;
289   m[0][2]=-2/kk*x12*z12;
290   m[1][0]=-2/kk*x12*y12;
291   m[1][1]=2-2/kk*y12*y12;
292   m[1][2]=-2/kk*y12*z12;
293   m[2][0]=-2/kk*x12*z12;
294   m[2][1]=-2*y12*z12;
295   m[2][2]=2-2/kk*z12*z12;
296   d[0]=2*p0[0]-2/kk*p0[0]*x12*x12-2/kk*p0[2]*x12*z12-2/kk*p0[1]*x12*y12;
297   d[1]=2*p0[1]-2/kk*p0[1]*y12*y12-2/kk*p0[0]*x12*y12-2/kk*p0[2]*z12*y12;
298   d[2]=2*p0[2]-2/kk*p0[2]*z12*z12-2/kk*p0[0]*x12*z12-2/kk*p0[1]*z12*y12;
299
300 }
301 //--------------------------------------------------------------------------  
302 void AliVertexerTracks::GetStrLinDerivMatrix(Double_t *p0,Double_t *p1,Double_t *sigmasq,Double_t (*m)[3],Double_t *d)
303 {
304   //
305   Double_t x12=p1[0]-p0[0];
306   Double_t y12=p1[1]-p0[1];
307   Double_t z12=p1[2]-p0[2];
308
309   Double_t den= x12*x12*sigmasq[1]*sigmasq[2]+y12*y12*sigmasq[0]*sigmasq[2]+z12*z12*sigmasq[0]*sigmasq[1];
310
311   Double_t kk= 2*(x12*x12/sigmasq[0]+y12*y12/sigmasq[1]+z12*z12/sigmasq[2]);
312
313   Double_t cc[3];
314   cc[0]=-x12/sigmasq[0];
315   cc[1]=-y12/sigmasq[1];
316   cc[2]=-z12/sigmasq[2];
317
318   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;
319
320   Double_t ss= -p0[0]*cc[0]-p0[1]*cc[1]-p0[2]*cc[2];
321
322   Double_t aa[3];
323   aa[0]=x12*sigmasq[1]*sigmasq[2]/den;
324   aa[1]=y12*sigmasq[0]*sigmasq[2]/den;
325   aa[2]=z12*sigmasq[0]*sigmasq[1]/den;
326
327   m[0][0]=aa[0]*(aa[0]*kk+2*cc[0])+2*cc[0]*aa[0]+2/sigmasq[0];
328   m[0][1]=aa[1]*(aa[0]*kk+2*cc[0])+2*cc[1]*aa[0];
329   m[0][2]=aa[2]*(aa[0]*kk+2*cc[0])+2*cc[2]*aa[0];
330
331   m[1][0]=aa[0]*(aa[1]*kk+2*cc[1])+2*cc[0]*aa[1];
332   m[1][1]=aa[1]*(aa[1]*kk+2*cc[1])+2*cc[1]*aa[1]+2/sigmasq[1];
333   m[1][2]=aa[2]*(aa[1]*kk+2*cc[1])+2*cc[2]*aa[1];
334
335   m[2][0]=aa[0]*(aa[2]*kk+2*cc[2])+2*cc[0]*aa[2];
336   m[2][1]=aa[1]*(aa[2]*kk+2*cc[2])+2*cc[1]*aa[2];
337   m[2][2]=aa[2]*(aa[2]*kk+2*cc[2])+2*cc[2]*aa[2]+2/sigmasq[2];
338
339   d[0]=-ww*(aa[0]*kk+2*cc[0])-2*ss*aa[0]+2*p0[0]/sigmasq[0];
340   d[1]=-ww*(aa[1]*kk+2*cc[1])-2*ss*aa[1]+2*p0[1]/sigmasq[1];
341   d[2]=-ww*(aa[2]*kk+2*cc[2])-2*ss*aa[2]+2*p0[2]/sigmasq[2];
342
343   }
344 //--------------------------------------------------------------------------   
345 Double_t AliVertexerTracks::GetStrLinMinDist(Double_t *p0,Double_t *p1,Double_t *x0)
346 {
347   //
348   Double_t x12=p0[0]-p1[0];
349   Double_t y12=p0[1]-p1[1];
350   Double_t z12=p0[2]-p1[2];
351   Double_t x10=p0[0]-x0[0];
352   Double_t y10=p0[1]-x0[1];
353   Double_t z10=p0[2]-x0[2];
354   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);
355 }
356 //---------------------------------------------------------------------------
357 void AliVertexerTracks::OneTrackVertFinder() 
358 {
359   // find vertex for events with 1 track, using DCA to nominal beam axis
360   if(fDebug) printf("Number of prepared tracks =%d - Call OneTrackVertFinder",fTrkArray.GetEntries());
361   AliESDtrack *track1;
362   track1 = (AliESDtrack*)fTrkArray.At(0);
363   Double_t field=GetField();
364   Double_t alpha=track1->GetAlpha();
365   Double_t mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
366   Double_t pos[3],dir[3]; 
367   track1->GetXYZAt(mindist,field,pos);
368   track1->GetPxPyPzAt(mindist,field,dir);
369   AliStrLine *line1 = new AliStrLine(pos,dir);
370   Double_t p1[3]={fNominalPos[0],fNominalPos[1],0.}; 
371   Double_t p2[3]={fNominalPos[0],fNominalPos[1],10.}; 
372   AliStrLine *zeta=new AliStrLine(p1,p2,kTRUE);
373   Double_t crosspoint[3]={0.,0.,0.};
374   Double_t sigma=999.;
375   Int_t nContrib=-1;
376   Int_t retcode = zeta->Cross(line1,crosspoint);
377   if(retcode>=0){
378     sigma=line1->GetDistFromPoint(crosspoint);
379     nContrib=1;
380   }
381   delete zeta;
382   delete line1;
383   fVert.SetXYZ(crosspoint);
384   fVert.SetDispersion(sigma);
385   fVert.SetNContributors(nContrib);  
386 }
387 //---------------------------------------------------------------------------
388 void AliVertexerTracks::HelixVertexFinder() 
389 {
390   // Get estimate of vertex position in (x,y) from tracks DCA
391
392
393   Double_t initPos[3];
394   initPos[2] = 0.;
395   for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
396   Double_t field=GetField();
397
398   Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
399
400   Double_t aver[3]={0.,0.,0.};
401   Double_t averquad[3]={0.,0.,0.};
402   Double_t sigmaquad[3]={0.,0.,0.};
403   Double_t sigma=0;
404   Int_t ncombi = 0;
405   AliESDtrack *track1;
406   AliESDtrack *track2;
407   Double_t distCA;
408   Double_t x, par[5];
409   Double_t alpha, cs, sn;
410   Double_t crosspoint[3];
411   for(Int_t i=0; i<nacc; i++){
412     track1 = (AliESDtrack*)fTrkArray.At(i);
413     
414
415     for(Int_t j=i+1; j<nacc; j++){
416       track2 = (AliESDtrack*)fTrkArray.At(j);
417
418       distCA=track2->PropagateToDCA(track1,field);
419       if(fDCAcut<=0 ||(fDCAcut>0&&distCA<fDCAcut)){
420         track1->GetExternalParameters(x,par);
421         alpha=track1->GetAlpha();
422         cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
423         Double_t x1=x*cs - par[0]*sn;
424         Double_t y1=x*sn + par[0]*cs;
425         Double_t z1=par[1];
426         Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2(); 
427         track2->GetExternalParameters(x,par);
428         alpha=track2->GetAlpha();
429         cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
430         Double_t x2=x*cs - par[0]*sn;
431         Double_t y2=x*sn + par[0]*cs;
432         Double_t z2=par[1];
433         Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2();
434         Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2();
435         Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
436         Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
437         Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
438         crosspoint[0]=wx1*x1 + wx2*x2; 
439         crosspoint[1]=wy1*y1 + wy2*y2; 
440         crosspoint[2]=wz1*z1 + wz2*z2;
441
442         ncombi++;
443         for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
444         for(Int_t jj=0;jj<3;jj++)averquad[jj]+=(crosspoint[jj]*crosspoint[jj]);
445       }
446     }
447       
448   }
449   if(ncombi>0){
450     for(Int_t jj=0;jj<3;jj++){
451       initPos[jj] = aver[jj]/ncombi;
452       averquad[jj]/=ncombi;
453       sigmaquad[jj]=averquad[jj]-initPos[jj]*initPos[jj];
454       sigma+=sigmaquad[jj];
455     }
456     sigma=TMath::Sqrt(TMath::Abs(sigma));
457   }
458   else {
459     Warning("HelixVertexFinder","Finder did not succed");
460     sigma=999;
461   }
462   fVert.SetXYZ(initPos);
463   fVert.SetDispersion(sigma);
464   fVert.SetNContributors(ncombi);
465 }
466 //----------------------------------------------------------------------------
467 Int_t AliVertexerTracks::PrepareTracks(TTree &trkTree,Int_t optImpParCut) 
468 {
469 //
470 // Propagate tracks to initial vertex position and store them in a TObjArray
471 //
472   Double_t maxd0rphi; 
473   Double_t maxd0z0 = fMaxd0z0; // default is 5 mm  
474   Int_t    nTrks    = 0;
475   Double_t sigmaCurr[3];
476   Double_t normdistx,normdisty;
477   Float_t  d0z0[2],covd0z0[3]; 
478   Double_t sigma;
479   Double_t field=GetField();
480
481   AliESDVertex *initVertex = new AliESDVertex(fNominalPos,fNominalCov,1,1);
482
483   Int_t    nEntries = (Int_t)trkTree.GetEntries();
484   if(!fTrkArray.IsEmpty()) fTrkArray.Delete();
485
486   if(fDebug) {
487     printf(" PrepareTracks()\n");
488   }
489
490   for(Int_t i=0; i<nEntries; i++) {
491     AliESDtrack *track = new AliESDtrack; 
492     trkTree.SetBranchAddress("tracks",&track);
493     trkTree.GetEvent(i);
494
495     // propagate track to vertex
496     if(optImpParCut<=1 || fOnlyFitter) { // optImpParCut==1 or 0
497       track->RelateToVertex(initVertex,field,100.);
498     } else {              // optImpParCut==2
499       fCurrentVertex->GetSigmaXYZ(sigmaCurr);
500       normdistx = TMath::Abs(fCurrentVertex->GetXv()-fNominalPos[0])/TMath::Sqrt(sigmaCurr[0]*sigmaCurr[0]+fNominalCov[0]);
501       normdisty = TMath::Abs(fCurrentVertex->GetYv()-fNominalPos[1])/TMath::Sqrt(sigmaCurr[1]*sigmaCurr[1]+fNominalCov[2]);
502       if(normdistx < 3. && normdisty < 3. &&
503          (sigmaCurr[0]+sigmaCurr[1])<(TMath::Sqrt(fNominalCov[0])+TMath::Sqrt(fNominalCov[2]))) {
504         track->RelateToVertex(fCurrentVertex,field,100.);
505       } else {
506         track->RelateToVertex(initVertex,field,100.);
507       }
508     }
509
510     track->GetImpactParameters(d0z0,covd0z0);
511     sigma = TMath::Sqrt(covd0z0[0]);
512     maxd0rphi = fNSigma*sigma;
513     if(optImpParCut==1) maxd0rphi *= 5.;
514
515
516
517     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);
518
519     // during iterations 1 and 2, if fConstraint=kFALSE,
520     // select tracks with d0oplusz0 < maxd0z0
521     if(optImpParCut>=1 && !fConstraint && nEntries>=3 && 
522        fVert.GetNContributors()>0) {
523       if(TMath::Sqrt(d0z0[0]*d0z0[0]+d0z0[1]*d0z0[1]) > maxd0z0) { 
524         if(fDebug) printf("     rejected\n");
525         delete track; continue; 
526       }
527     }
528
529     // select tracks with d0rphi < maxd0rphi
530     if(optImpParCut>0 && TMath::Abs(d0z0[0]) > maxd0rphi) { 
531       if(fDebug) printf("     rejected\n");
532       delete track; continue; 
533     }
534
535     fTrkArray.AddLast(track);
536     nTrks++; 
537   }
538
539   delete initVertex;
540
541   return nTrks;
542
543 //---------------------------------------------------------------------------
544 AliESDVertex* AliVertexerTracks::RemoveTracksFromVertex(AliESDVertex *inVtx,
545                                                         TTree *trksTree,
546                                                         Float_t *diamondxy) 
547 {
548 //
549 // Removes tracks in trksTree from fit of inVtx
550 //
551
552   if(!strstr(inVtx->GetTitle(),"VertexerTracksWithConstraintOnlyFitter"))
553     printf("WARNING: result of tracks' removal will be only approximately correct\n");
554
555   TMatrixD rv(3,1);
556   rv(0,0) = inVtx->GetXv();
557   rv(1,0) = inVtx->GetYv();
558   rv(2,0) = inVtx->GetZv();
559   TMatrixD vV(3,3);
560   Double_t cov[6];
561   inVtx->GetCovMatrix(cov);
562   vV(0,0) = cov[0];
563   vV(0,1) = cov[1]; vV(1,0) = cov[1];
564   vV(1,1) = cov[2];
565   vV(0,2) = cov[3]; vV(2,0) = cov[3];
566   vV(1,2) = cov[4]; vV(2,1) = cov[4]; 
567   vV(2,2) = cov[5];
568
569   TMatrixD sumWi(TMatrixD::kInverted,vV);
570   TMatrixD sumWiri(sumWi,TMatrixD::kMult,rv);
571
572   Int_t nUsedTrks = inVtx->GetNContributors();
573   Double_t chi2 = inVtx->GetChi2();
574
575   AliESDtrack *track = 0;
576   trksTree->SetBranchAddress("tracks",&track);
577   Int_t ntrks = trksTree->GetEntries();
578   for(Int_t i=0;i<ntrks;i++) {
579     trksTree->GetEvent(i);
580     if(!inVtx->UsesTrack(track->GetID())) {
581       printf("track %d was not used in vertex fit\n",track->GetID());
582       continue;
583     }
584     Double_t alpha = track->GetAlpha();
585     Double_t xl = diamondxy[0]*TMath::Cos(alpha)+diamondxy[1]*TMath::Sin(alpha);
586     track->AliExternalTrackParam::PropagateTo(xl,AliTracker::GetBz()); 
587     // vector of track global coordinates
588     TMatrixD ri(3,1);
589     // covariance matrix of ri
590     TMatrixD wWi(3,3);
591     
592     // get space point from track
593     if(!TrackToPoint(track,ri,wWi)) continue;
594
595     TMatrixD wWiri(wWi,TMatrixD::kMult,ri); 
596
597     sumWi -= wWi;
598     sumWiri -= wWiri;
599
600     // track chi2
601     TMatrixD deltar = rv; deltar -= ri;
602     TMatrixD wWideltar(wWi,TMatrixD::kMult,deltar);
603     Double_t chi2i = deltar(0,0)*wWideltar(0,0)+
604                      deltar(1,0)*wWideltar(1,0)+
605                      deltar(2,0)*wWideltar(2,0);
606     // remove from total chi2
607     chi2 -= chi2i;
608
609     nUsedTrks--;
610     if(nUsedTrks<2) {
611       printf("Trying to remove too many tracks!\n");
612       return 0x0;
613     }
614   }
615
616   TMatrixD rvnew(3,1);
617   TMatrixD vVnew(3,3);
618
619   // new inverted of weights matrix
620   TMatrixD invsumWi(TMatrixD::kInverted,sumWi);
621   vVnew = invsumWi;
622   // new position of primary vertex
623   rvnew.Mult(vVnew,sumWiri);
624
625   Double_t position[3];
626   position[0] = rvnew(0,0);
627   position[1] = rvnew(1,0);
628   position[2] = rvnew(2,0);
629   cov[0] = vVnew(0,0);
630   cov[1] = vVnew(0,1);
631   cov[2] = vVnew(1,1);
632   cov[3] = vVnew(0,2);
633   cov[4] = vVnew(1,2);
634   cov[5] = vVnew(2,2);
635   
636   // store data in the vertex object
637   AliESDVertex *outVtx = new AliESDVertex(position,cov,chi2,nUsedTrks);
638   outVtx->SetTitle(inVtx->GetTitle());
639   Double_t tp[3];
640   inVtx->GetTruePos(tp);
641   outVtx->SetTruePos(tp);
642   UShort_t *inindices = inVtx->GetIndices();
643   UShort_t *outindices = new UShort_t[outVtx->GetNContributors()];
644   Int_t j=0;
645   Bool_t copyindex;
646   for(Int_t k=0; k<inVtx->GetNIndices(); k++) {
647     copyindex=kTRUE;
648     for(Int_t l=0; l<ntrks; l++) {
649       trksTree->GetEvent(l);
650       if(inindices[k]==track->GetID()) copyindex=kFALSE;
651     }
652     if(copyindex) {
653       outindices[j] = inindices[k]; j++;
654     }
655   }
656   outVtx->SetIndices(outVtx->GetNContributors(),outindices);
657   delete [] outindices;
658
659   if(fDebug) {
660     printf("Vertex before removing tracks:\n");
661     inVtx->PrintStatus();
662     inVtx->PrintIndices();
663     printf("Vertex after removing tracks:\n");
664     outVtx->PrintStatus();
665     outVtx->PrintIndices();
666   }
667
668   return outVtx;
669 }
670 //---------------------------------------------------------------------------
671 void AliVertexerTracks::SetSkipTracks(Int_t n,Int_t *skipped) 
672 {
673 //
674 // Mark the tracks not to be used in the vertex reconstruction.
675 // Tracks are identified by AliESDtrack::GetID()
676 //
677   fNTrksToSkip = n;  fTrksToSkip = new Int_t[n]; 
678   for(Int_t i=0;i<n;i++) fTrksToSkip[i] = skipped[i]; 
679   return; 
680 }
681 //---------------------------------------------------------------------------
682 void  AliVertexerTracks::SetVtxStart(AliESDVertex *vtx) 
683
684 //
685 // Set initial vertex knowledge
686 //
687   vtx->GetXYZ(fNominalPos);
688   vtx->GetCovMatrix(fNominalCov);
689   SetConstraintOn();
690   return; 
691 }
692 //---------------------------------------------------------------------------
693 void AliVertexerTracks::StrLinVertexFinderMinDist(Int_t optUseWeights)
694 {
695   // Calculate the point at minimum distance to prepared tracks 
696   
697   Double_t initPos[3];
698   initPos[2] = 0.;
699   Double_t sigma=0;
700   for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
701   const Int_t knacc = (Int_t)fTrkArray.GetEntriesFast();
702   Double_t field=GetField();
703
704   AliESDtrack *track1;
705   Double_t (*vectP0)[3]=new Double_t [knacc][3];
706   Double_t (*vectP1)[3]=new Double_t [knacc][3];
707   
708   Double_t sum[3][3];
709   Double_t dsum[3]={0,0,0};
710   for(Int_t i=0;i<3;i++)
711     for(Int_t j=0;j<3;j++)sum[i][j]=0;
712   for(Int_t i=0; i<knacc; i++){
713     track1 = (AliESDtrack*)fTrkArray.At(i);
714     Double_t alpha=track1->GetAlpha();
715     Double_t mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
716     Double_t pos[3],dir[3]; 
717     track1->GetXYZAt(mindist,field,pos);
718     track1->GetPxPyPzAt(mindist,field,dir);
719     AliStrLine *line1 = new AliStrLine(pos,dir); 
720     //    AliStrLine *line1 = new AliStrLine();
721     //    track1->ApproximateHelixWithLine(mindist,field,line1);
722
723     Double_t p0[3],cd[3];
724     line1->GetP0(p0);
725     line1->GetCd(cd);
726     Double_t p1[3]={p0[0]+cd[0],p0[1]+cd[1],p0[2]+cd[2]};
727     vectP0[i][0]=p0[0];
728     vectP0[i][1]=p0[1];
729     vectP0[i][2]=p0[2];
730     vectP1[i][0]=p1[0];
731     vectP1[i][1]=p1[1];
732     vectP1[i][2]=p1[2];
733     
734     Double_t matr[3][3];
735     Double_t dknow[3];
736     if(optUseWeights==0)GetStrLinDerivMatrix(p0,p1,matr,dknow);
737     if(optUseWeights==1){
738       Double_t sigmasq[3];
739       sigmasq[0]=track1->GetSigmaY2();
740       sigmasq[1]=track1->GetSigmaY2();
741       sigmasq[2]=track1->GetSigmaZ2();
742       GetStrLinDerivMatrix(p0,p1,sigmasq,matr,dknow);
743     }
744
745     for(Int_t iii=0;iii<3;iii++){
746       dsum[iii]+=dknow[iii]; 
747       for(Int_t lj=0;lj<3;lj++) sum[iii][lj]+=matr[iii][lj];
748     }
749     delete line1;
750   }
751   
752   Double_t vett[3][3];
753   Double_t det=GetDeterminant3X3(sum);
754   
755   if(det!=0){
756     for(Int_t zz=0;zz<3;zz++){
757       for(Int_t ww=0;ww<3;ww++){
758         for(Int_t kk=0;kk<3;kk++) vett[ww][kk]=sum[ww][kk];
759       }
760       for(Int_t kk=0;kk<3;kk++) vett[kk][zz]=dsum[kk];
761       initPos[zz]=GetDeterminant3X3(vett)/det;
762     }
763
764
765     for(Int_t i=0; i<knacc; i++){
766       Double_t p0[3]={0,0,0},p1[3]={0,0,0};
767       for(Int_t ii=0;ii<3;ii++){
768         p0[ii]=vectP0[i][ii];
769         p1[ii]=vectP1[i][ii];
770       }
771       sigma+=GetStrLinMinDist(p0,p1,initPos);
772     }
773
774     sigma=TMath::Sqrt(sigma);
775   }else{
776     Warning("StrLinVertexFinderMinDist","Finder did not succed");
777     sigma=999;
778   }
779   delete vectP0;
780   delete vectP1;
781   fVert.SetXYZ(initPos);
782   fVert.SetDispersion(sigma);
783   fVert.SetNContributors(knacc);
784 }
785 //---------------------------------------------------------------------------
786 Bool_t AliVertexerTracks::TrackToPoint(AliESDtrack *t,
787                                        TMatrixD &ri,TMatrixD &wWi) const 
788 {
789 //
790 // Extract from the AliESDtrack the global coordinates ri and covariance matrix
791 // wWi of the space point that it represents (to be used in VertexFitter())
792 //
793
794   
795   Double_t rotAngle = t->GetAlpha();
796   if(rotAngle<0.) rotAngle += 2.*TMath::Pi();
797   Double_t cosRot = TMath::Cos(rotAngle);
798   Double_t sinRot = TMath::Sin(rotAngle);
799
800   ri(0,0) = t->GetX()*cosRot-t->GetY()*sinRot;
801   ri(1,0) = t->GetX()*sinRot+t->GetY()*cosRot;
802   ri(2,0) = t->GetZ();
803
804   // matrix to go from global (x,y,z) to local (y,z);
805   TMatrixD qQi(2,3);
806   qQi(0,0) = -sinRot;
807   qQi(0,1) = cosRot;
808   qQi(0,2) = 0.;
809   qQi(1,0) = 0.;
810   qQi(1,1) = 0.;
811   qQi(1,2) = 1.;
812
813   // covariance matrix of local (y,z) - inverted
814   TMatrixD uUi(2,2);
815   Double_t cc[15];
816   t->GetExternalCovariance(cc);
817   uUi(0,0) = cc[0];
818   uUi(0,1) = cc[1];
819   uUi(1,0) = cc[1];
820   uUi(1,1) = cc[2];
821   if(uUi.Determinant() <= 0.) return kFALSE;
822   TMatrixD uUiInv(TMatrixD::kInverted,uUi);
823   
824   // weights matrix: wWi = qQiT * uUiInv * qQi
825   TMatrixD uUiInvQi(uUiInv,TMatrixD::kMult,qQi);
826   TMatrixD m(qQi,TMatrixD::kTransposeMult,uUiInvQi);
827   wWi = m;
828
829   return kTRUE;
830
831 //---------------------------------------------------------------------------
832 void AliVertexerTracks::TooFewTracks(const AliESD* esdEvent) 
833 {
834 //
835 // When the number of tracks is < fMinTracks
836 //
837
838   // deal with vertices not found
839   Double_t pos[3],err[3];
840   Int_t    ncontr=0;
841   pos[0] = fNominalPos[0];
842   err[0] = TMath::Sqrt(fNominalCov[0]);
843   pos[1] = fNominalPos[1];
844   err[1] = TMath::Sqrt(fNominalCov[2]);
845   pos[2] = esdEvent->GetVertex()->GetZv();
846   err[2] = esdEvent->GetVertex()->GetZRes();
847   if(err[0]>1. && esdEvent->GetVertex()->GetNContributors()<=0) 
848     ncontr = -1; // (x,y,z) = (0,0,0) 
849   if(err[0]>1. && esdEvent->GetVertex()->GetNContributors()>0) 
850     ncontr = -2; // (x,y,z) = (0,0,z_fromSPD) 
851   if(err[0]<1. && esdEvent->GetVertex()->GetNContributors()<=0) 
852     ncontr = -3; // (x,y,z) = (x_mean,y_mean,0)
853   if(err[0]<1. && esdEvent->GetVertex()->GetNContributors()>0) 
854     ncontr = -4; // (x,y,z) = (x_mean,y_mean,z_fromSPD)
855   fCurrentVertex = 0;
856   fCurrentVertex = new AliESDVertex(pos,err);
857   fCurrentVertex->SetNContributors(ncontr);
858
859   Double_t tp[3];
860   esdEvent->GetVertex()->GetTruePos(tp);
861   fCurrentVertex->SetTruePos(tp);
862   if(fConstraint) {
863     fCurrentVertex->SetTitle("VertexerTracksWithConstraint");
864   } else {
865     fCurrentVertex->SetTitle("VertexerTracksNoConstraint");
866   }
867
868   return;
869 }
870 //---------------------------------------------------------------------------
871 void AliVertexerTracks::VertexFinder(Int_t optUseWeights) 
872 {
873
874   // Get estimate of vertex position in (x,y) from tracks DCA
875  
876   Double_t initPos[3];
877   initPos[2] = 0.;
878   for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
879   Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
880   Double_t aver[3]={0.,0.,0.};
881   Double_t aversq[3]={0.,0.,0.};
882   Double_t sigmasq[3]={0.,0.,0.};
883   Double_t sigma=0;
884   Int_t ncombi = 0;
885   AliESDtrack *track1;
886   AliESDtrack *track2;
887   Double_t pos[3],dir[3]; 
888   Double_t alpha,mindist;
889   Double_t field=GetField();
890
891   for(Int_t i=0; i<nacc; i++){
892     track1 = (AliESDtrack*)fTrkArray.At(i);
893     alpha=track1->GetAlpha();
894     mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
895     track1->GetXYZAt(mindist,field,pos);
896     track1->GetPxPyPzAt(mindist,field,dir);
897     AliStrLine *line1 = new AliStrLine(pos,dir); 
898
899    //    AliStrLine *line1 = new AliStrLine();
900    //    track1->ApproximateHelixWithLine(mindist,field,line1);
901    
902     for(Int_t j=i+1; j<nacc; j++){
903       track2 = (AliESDtrack*)fTrkArray.At(j);
904       alpha=track2->GetAlpha();
905       mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
906       track2->GetXYZAt(mindist,field,pos);
907       track2->GetPxPyPzAt(mindist,field,dir);
908       AliStrLine *line2 = new AliStrLine(pos,dir); 
909     //      AliStrLine *line2 = new AliStrLine();
910     //  track2->ApproximateHelixWithLine(mindist,field,line2);
911       Double_t distCA=line2->GetDCA(line1);
912       //printf("%d   %d   %f\n",i,j,distCA);
913        if(fDCAcut<=0 || (fDCAcut>0&&distCA<fDCAcut)){
914         Double_t pnt1[3],pnt2[3],crosspoint[3];
915
916         if(optUseWeights<=0){
917           Int_t retcode = line2->Cross(line1,crosspoint);
918           if(retcode>=0){
919             ncombi++;
920             for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
921             for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]);
922           }
923         }
924         if(optUseWeights>0){
925           Int_t retcode = line1->CrossPoints(line2,pnt1,pnt2);
926           if(retcode>=0){
927             Double_t cs, sn;
928             alpha=track1->GetAlpha();
929             cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);   
930             Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2();
931             alpha=track2->GetAlpha();
932             cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
933             Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2();
934             Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2();
935             Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
936             Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
937             Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
938             crosspoint[0]=wx1*pnt1[0] + wx2*pnt2[0]; 
939             crosspoint[1]=wy1*pnt1[1] + wy2*pnt2[1]; 
940             crosspoint[2]=wz1*pnt1[2] + wz2*pnt2[2];
941           
942             ncombi++;
943             for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
944             for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]);
945           }
946         }
947       }
948       delete line2;
949     }
950     delete line1;
951   }
952   if(ncombi>0){
953     for(Int_t jj=0;jj<3;jj++){
954       initPos[jj] = aver[jj]/ncombi;
955       //printf("%f\n",initPos[jj]);
956       aversq[jj]/=ncombi;
957       sigmasq[jj]=aversq[jj]-initPos[jj]*initPos[jj];
958       sigma+=sigmasq[jj];
959     }
960     sigma=TMath::Sqrt(TMath::Abs(sigma));
961   }
962   else {
963     Warning("VertexFinder","Finder did not succed");
964     sigma=999;
965   }
966   fVert.SetXYZ(initPos);
967   fVert.SetDispersion(sigma);
968   fVert.SetNContributors(ncombi);
969 }
970 //---------------------------------------------------------------------------
971 void AliVertexerTracks::VertexFitter(Bool_t useConstraint) 
972 {
973 //
974 // The optimal estimate of the vertex position is given by a "weighted 
975 // average of tracks positions"
976 // Original method: V. Karimaki, CMS Note 97/0051
977 //
978   Double_t initPos[3];
979   fVert.GetXYZ(initPos);
980   Int_t arrEntries = (Int_t)fTrkArray.GetEntries();
981   if(arrEntries==1) useConstraint=kTRUE;
982   if(fDebug) { 
983     printf(" VertexFitter(): start\n");
984     printf(" Number of tracks in array: %d\n",(Int_t)fTrkArray.GetEntriesFast());
985     printf(" Minimum # tracks required in fit: %d\n",fMinTracks);
986     printf("Vertex position after finder: %f,%f,%f\n",initPos[0],initPos[1],initPos[2]);
987     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])); 
988   }
989
990   Int_t i,j,k,step=0;
991   TMatrixD rv(3,1);
992   TMatrixD vV(3,3);
993   rv(0,0) = initPos[0];
994   rv(1,0) = initPos[1];
995   rv(2,0) = 0.;
996   Double_t xlStart,alpha;
997   Int_t nUsedTrks;
998   Double_t chi2,chi2i,chi2b;
999   AliESDtrack *t = 0;
1000   Int_t failed = 0;
1001
1002   // initial vertex covariance matrix
1003   TMatrixD vVb(3,3);
1004   vVb(0,0) = fNominalCov[0];
1005   vVb(0,1) = fNominalCov[1];
1006   vVb(0,2) = 0.;
1007   vVb(1,0) = fNominalCov[1];
1008   vVb(1,1) = fNominalCov[2];
1009   vVb(1,2) = 0.;
1010   vVb(2,0) = 0.;
1011   vVb(2,1) = 0.;
1012   vVb(2,2) = fNominalCov[5];
1013   TMatrixD vVbInv(TMatrixD::kInverted,vVb);
1014   TMatrixD rb(3,1);
1015   rb(0,0) = fNominalPos[0];
1016   rb(1,0) = fNominalPos[1];
1017   rb(2,0) = fNominalPos[2];
1018   TMatrixD vVbInvrb(vVbInv,TMatrixD::kMult,rb);
1019
1020
1021   // 2 steps:
1022   // 1st - estimate of vtx using all tracks
1023   // 2nd - estimate of global chi2
1024   for(step=0; step<2; step++) {
1025     if(fDebug) printf(" step = %d\n",step);
1026     chi2 = 0.;
1027     nUsedTrks = 0;
1028
1029     TMatrixD sumWiri(3,1);
1030     TMatrixD sumWi(3,3);
1031     for(i=0; i<3; i++) {
1032       sumWiri(i,0) = 0.;
1033       for(j=0; j<3; j++) sumWi(j,i) = 0.;
1034     }
1035
1036     // mean vertex constraint
1037     if(useConstraint) {
1038       for(i=0;i<3;i++) {
1039         sumWiri(i,0) += vVbInvrb(i,0);
1040         for(k=0;k<3;k++) sumWi(i,k) += vVbInv(i,k);
1041       }
1042       // chi2
1043       TMatrixD deltar = rv; deltar -= rb;
1044       TMatrixD vVbInvdeltar(vVbInv,TMatrixD::kMult,deltar);
1045       chi2b = deltar(0,0)*vVbInvdeltar(0,0)+
1046               deltar(1,0)*vVbInvdeltar(1,0)+
1047               deltar(2,0)*vVbInvdeltar(2,0);
1048       chi2 += chi2b;
1049     }
1050
1051
1052     // loop on tracks  
1053     for(k=0; k<arrEntries; k++) {
1054       // get track from track array
1055       t = (AliESDtrack*)fTrkArray.At(k);
1056       alpha = t->GetAlpha();
1057       xlStart = initPos[0]*TMath::Cos(alpha)+initPos[1]*TMath::Sin(alpha);
1058       // to vtxSeed (from finder)
1059       t->AliExternalTrackParam::PropagateTo(xlStart,AliTracker::GetBz());   
1060  
1061
1062       // vector of track global coordinates
1063       TMatrixD ri(3,1);
1064       // covariance matrix of ri
1065       TMatrixD wWi(3,3);
1066
1067       // get space point from track
1068       if(!TrackToPoint(t,ri,wWi)) continue;
1069
1070       TMatrixD wWiri(wWi,TMatrixD::kMult,ri); 
1071
1072       // track chi2
1073       TMatrixD deltar = rv; deltar -= ri;
1074       TMatrixD wWideltar(wWi,TMatrixD::kMult,deltar);
1075       chi2i = deltar(0,0)*wWideltar(0,0)+
1076               deltar(1,0)*wWideltar(1,0)+
1077               deltar(2,0)*wWideltar(2,0);
1078
1079       // add to total chi2
1080       chi2 += chi2i;
1081
1082       sumWiri += wWiri;
1083       sumWi   += wWi;
1084
1085       nUsedTrks++;
1086     } // end loop on tracks
1087
1088     if(nUsedTrks < fMinTracks) {
1089       failed=1;
1090       continue;
1091     }
1092
1093     Double_t determinant = sumWi.Determinant();
1094     //cerr<<" determinant: "<<determinant<<endl;
1095     if(determinant < 100.)  { 
1096       printf("det(V) = 0\n");       
1097       failed=1;
1098       continue;
1099     }
1100
1101     // inverted of weights matrix
1102     TMatrixD invsumWi(TMatrixD::kInverted,sumWi);
1103     vV = invsumWi;
1104      
1105     // position of primary vertex
1106     rv.Mult(vV,sumWiri);
1107
1108   } // end loop on the 2 steps
1109
1110   // delete t;
1111
1112   if(failed) { 
1113     if(fDebug) printf("TooFewTracks\n");
1114     fCurrentVertex = new AliESDVertex(0.,0.,-1);
1115     return; 
1116   }
1117
1118   Double_t position[3];
1119   position[0] = rv(0,0);
1120   position[1] = rv(1,0);
1121   position[2] = rv(2,0);
1122   Double_t covmatrix[6];
1123   covmatrix[0] = vV(0,0);
1124   covmatrix[1] = vV(0,1);
1125   covmatrix[2] = vV(1,1);
1126   covmatrix[3] = vV(0,2);
1127   covmatrix[4] = vV(1,2);
1128   covmatrix[5] = vV(2,2);
1129   
1130   // store data in the vertex object
1131   fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nUsedTrks);
1132
1133   if(fDebug) {
1134     printf(" VertexFitter(): finish\n");
1135     printf(" rv = ( %f , %f , %f )\n\n",rv(0,0),rv(1,0),rv(2,0));
1136     fCurrentVertex->PrintStatus();
1137   }
1138
1139   return;
1140 }
1141 //----------------------------------------------------------------------------
1142 AliESDVertex* AliVertexerTracks::VertexForSelectedTracks(TTree *trkTree,Bool_t optUseFitter,Bool_t optPropagate) 
1143 {
1144 //
1145 // Return vertex from tracks in trkTree
1146 //
1147
1148   SetConstraintOff();
1149
1150   // get tracks and propagate them to initial vertex position
1151   Int_t nTrksPrep = PrepareTracks(*trkTree,0);
1152   if(nTrksPrep <  TMath::Max(2,fMinTracks) ) {
1153     if(fDebug) printf("TooFewTracks\n");
1154     fCurrentVertex = new AliESDVertex(0.,0.,-1);
1155     return fCurrentVertex;
1156   }
1157  
1158   switch (fAlgo) {
1159     case 1: StrLinVertexFinderMinDist(1); break;
1160     case 2: StrLinVertexFinderMinDist(0); break;
1161     case 3: HelixVertexFinder();          break;
1162     case 4: VertexFinder(1);              break;
1163     case 5: VertexFinder(0);              break;
1164     default: printf("Wrong algorithm\n"); break;  
1165   }
1166   
1167   if(fDebug) printf(" Vertex finding completed\n");
1168
1169   // vertex fitter
1170   if(optUseFitter){
1171     VertexFitter(fConstraint);
1172     if(fDebug) printf(" Vertex fit completed\n");
1173   }else{
1174     Double_t position[3]={fVert.GetXv(),fVert.GetYv(),fVert.GetZv()};
1175     Double_t covmatrix[6]={0.,0.,0.,0.,0.,0.};
1176     Double_t chi2=99999.;
1177     Int_t    nUsedTrks=fVert.GetNContributors();
1178     fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nUsedTrks);    
1179   }
1180   fCurrentVertex->SetDispersion(fVert.GetDispersion());
1181
1182
1183   // set indices of used tracks and propagate track to found vertex
1184   UShort_t *indices = 0;
1185   AliESDtrack *eta = 0;
1186   if(fCurrentVertex->GetNContributors()>0) {
1187     indices = new UShort_t[fCurrentVertex->GetNContributors()];
1188     for(Int_t jj=0;jj<(Int_t)fTrkArray.GetEntriesFast();jj++) {
1189       eta = (AliESDtrack*)fTrkArray.At(jj);
1190       indices[jj] = (UShort_t)eta->GetID();
1191       if(optPropagate&&optUseFitter){
1192         if(TMath::Sqrt(fCurrentVertex->GetXv()*fCurrentVertex->GetXv()+fCurrentVertex->GetYv()*fCurrentVertex->GetYv())<3.) {
1193           eta->RelateToVertex(fCurrentVertex,GetField(),100.);
1194           if(fDebug) printf("Track %d propagated to found vertex\n",jj);
1195         }else{
1196           AliWarning("Found vertex outside beam pipe!");
1197         }
1198       }
1199     }
1200     fCurrentVertex->SetIndices(fCurrentVertex->GetNContributors(),indices);
1201   }
1202   delete [] indices;
1203   fTrkArray.Delete();
1204   
1205   return fCurrentVertex;
1206 }
1207 //----------------------------------------------------------------------------
1208 AliESDVertex* AliVertexerTracks::VertexForSelectedTracks(TObjArray *trkArray,Bool_t optUseFitter, Bool_t optPropagate) 
1209 {
1210 //
1211 // Return vertex from array of tracks
1212 //
1213
1214   // get tracks and propagate them to initial vertex position
1215   Int_t nTrks = trkArray->GetEntriesFast();
1216   if(nTrks < TMath::Max(2,fMinTracks) ) {
1217     if(fDebug) printf("TooFewTracks\n");
1218     fCurrentVertex = new AliESDVertex(0.,0.,-1);
1219     return fCurrentVertex;
1220   }
1221   TTree *trkTree = new TTree("TreeT","tracks");
1222   AliESDtrack *esdTrack = 0;
1223   trkTree->Branch("tracks","AliESDtrack",&esdTrack);
1224   for(Int_t i=0; i<nTrks; i++){
1225     esdTrack = (AliESDtrack*)trkArray->At(i);
1226     trkTree->Fill();
1227   }
1228     
1229   AliESDVertex *vtx =  VertexForSelectedTracks(trkTree,optUseFitter,optPropagate);
1230   delete trkTree;
1231   return vtx;
1232 }
1233 //--------------------------------------------------------------------------