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