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