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