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