]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliVertexerTracks.cxx
New version of AliVertexerTracks (F.Prino)
[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(),fVert(),fCurrentVertex(0)
43 {
44 //
45 // Default constructor
46 //
47   SetVtxStart();
48   SetMinTracks();
49   fDCAcut=0;
50   fTrksToSkip = 0;
51   fNTrksToSkip = 0;
52   fDebug=0;
53   fAlgo=1;
54 }
55 //-----------------------------------------------------------------------------
56 AliVertexerTracks::AliVertexerTracks(Double_t xStart, Double_t yStart):
57  TObject(),fVert(),fCurrentVertex(0)
58 {
59 //
60 // Standard constructor
61 //
62   SetVtxStart(xStart,yStart);
63   SetMinTracks();
64   fTrksToSkip = 0;
65   fNTrksToSkip = 0;
66   fDCAcut=0;
67   fDebug=0;
68   fAlgo=1;
69 }
70 //-----------------------------------------------------------------------------
71 AliVertexerTracks::~AliVertexerTracks() {
72   // Default Destructor
73   // The objects poited by the following pointers are not owned
74   // by this class and are not deleted
75   fCurrentVertex = 0;
76 }
77 //---------------------------------------------------------------------------
78 void AliVertexerTracks::SetSkipTracks(Int_t n,Int_t *skipped) {
79 //
80 // Mark the tracks not ot be used in the vertex finding
81 //
82   fNTrksToSkip = n;
83   fTrksToSkip = new Int_t[n]; 
84   for(Int_t i=0;i<n;i++) fTrksToSkip[i] = skipped[i]; 
85   return; 
86 }
87 //---------------------------------------------------------------------------
88 void AliVertexerTracks::ComputeMaxChi2PerTrack(Int_t nTracks) {
89 //
90 // Max. contr. to the chi2 has been tuned as a function of multiplicity
91 //
92   if(nTracks < 7) { fMaxChi2PerTrack = 1.e6;
93   } else { fMaxChi2PerTrack = 100.; }
94
95   return;
96 }
97 //----------------------------------------------------------------------------
98 Int_t AliVertexerTracks::PrepareTracks(TTree &trkTree, Int_t OptImpParCut) {
99 //
100 // Propagate tracks to initial vertex position and store them in a TObjArray
101 //
102   Double_t maxd0rphi = 3.;  
103   Double_t alpha,xlStart,d0rphi;
104   Int_t    nTrks    = 0;
105   Bool_t   skipThis;
106   Double_t field=GetField();
107
108   Int_t    nEntries = (Int_t)trkTree.GetEntries();
109   if(!fTrkArray.IsEmpty()) fTrkArray.Clear();
110   fTrkArray.Expand(nEntries);
111
112   if(fDebug) {
113     printf(" PrepareTracks()\n");
114   }
115
116   for(Int_t i=0; i<nEntries; i++) {
117     // check tracks to skip
118     skipThis = kFALSE;
119     for(Int_t j=0; j<fNTrksToSkip; j++) { 
120       if(i==fTrksToSkip[j]) {
121         if(fDebug) printf("skipping track: %d\n",i);
122         skipThis = kTRUE;
123       }
124     }
125     if(skipThis) continue;
126
127     AliESDtrack *track = new AliESDtrack; 
128     trkTree.SetBranchAddress("tracks",&track);
129     trkTree.GetEvent(i);
130
131     // propagate track to vtxSeed
132     alpha  = track->GetAlpha();
133     xlStart = fNominalPos[0]*TMath::Cos(alpha)+fNominalPos[1]*TMath::Sin(alpha);
134     track->PropagateTo(xlStart,field);   // to vtxSeed
135     d0rphi = TMath::Abs(track->GetD(fNominalPos[0],fNominalPos[1],field));
136     if(OptImpParCut==1 && d0rphi > maxd0rphi) { delete track; continue; }
137     if(OptImpParCut==2) { 
138       // to be replaced by new Andrea's selection
139       if(d0rphi > maxd0rphi){ 
140         delete track; 
141         continue; 
142       }
143       // end to be replaced by new Andrea's selection
144     }
145     fTrkArray.AddLast(track);
146     nTrks++; 
147     if(fDebug) printf(" :-) nTrks=%d, d0rphi=%f\n",nTrks,d0rphi);
148   }
149   if(fTrksToSkip) delete [] fTrksToSkip;
150   return nTrks;
151
152 //----------------------------------------------------------------------------
153 AliVertex* AliVertexerTracks::VertexForSelectedTracks(TTree *trkTree) {
154 //
155 // Return vertex from tracks in trkTree
156 //
157
158   // get tracks and propagate them to initial vertex position
159   Int_t nTrks = PrepareTracks(*trkTree,1);
160   if(nTrks < fMinTracks) {
161     printf("TooFewTracks\n");
162     Double_t vtx[3]={0,0,0};
163     fVert.SetXYZ(vtx);
164     fVert.SetDispersion(999);
165     fVert.SetNContributors(-5);
166     return &fVert;
167   }
168  
169   // Set initial vertex position from ESD
170   if(fAlgo==1)  StrLinVertexFinderMinDist(1);
171   if(fAlgo==2)  StrLinVertexFinderMinDist(0);
172   if(fAlgo==3)  HelixVertexFinder();
173   if(fAlgo==4)  VertexFinder(1);
174   if(fAlgo==5)  VertexFinder(0);
175   return &fVert;
176 }
177
178 //----------------------------------------------------------------------------
179 AliESDVertex* AliVertexerTracks::FindPrimaryVertex(const AliESD *esdEvent)
180 {
181 //
182 // Primary vertex for current ESD event
183 //
184   fCurrentVertex = 0;
185
186   Int_t entr = (Int_t)esdEvent->GetNumberOfTracks();
187   TTree *trkTree = new TTree("TreeT","tracks");
188   AliESDtrack *esdTrack = 0;
189   trkTree->Branch("tracks","AliESDtrack",&esdTrack);
190
191   for(Int_t i=0; i<entr; i++) {
192     AliESDtrack *et = esdEvent->GetTrack(i);
193     esdTrack = new AliESDtrack(*et);
194     if(!esdTrack->GetStatus()&AliESDtrack::kITSin) continue;
195     if(!esdTrack->GetStatus()&AliESDtrack::kITSrefit) continue;
196     Int_t nclus=esdTrack->GetNcls(0); // check number of clusters in ITS
197     if(nclus<6) continue;
198
199     trkTree->Fill();
200   }
201   delete esdTrack;
202
203   // preselect tracks and propagate them to initial vertex position
204   Int_t nTrks = PrepareTracks(*trkTree,2);
205   delete trkTree;
206   if(fDebug) printf(" tracks prepared: %d\n",nTrks);
207   if(nTrks < fMinTracks) {
208     printf("TooFewTracks\n");
209     fCurrentVertex = new AliESDVertex(0.,0.,-1);
210     return fCurrentVertex; 
211   }
212
213   // VERTEX FINDER
214   switch (fAlgo) {
215     case 1: StrLinVertexFinderMinDist(1); break;
216     case 2: StrLinVertexFinderMinDist(0); break;
217     case 3: HelixVertexFinder();          break;
218     case 4: VertexFinder(1);              break;
219     case 5: VertexFinder(0);              break;
220     default: printf("Wrong algorithm\n"); break;  
221   }
222
223   // VERTEX FITTER
224   ComputeMaxChi2PerTrack(nTrks);
225   VertexFitter();
226   if(fDebug) printf(" vertex fit completed\n");
227
228
229   Double_t tp[3];
230   esdEvent->GetVertex()->GetTruePos(tp);
231   fCurrentVertex->SetTruePos(tp);
232   fCurrentVertex->SetTitle("VertexerTracks");
233
234   fTrkArray.Clear();
235   return fCurrentVertex;
236 }
237 //----------------------------------------------------------------------------
238 AliVertex* AliVertexerTracks::VertexForSelectedTracks(TObjArray *trkArray) {
239 //
240 // Return vertex from array of tracks
241 //
242
243   // get tracks and propagate them to initial vertex position
244   Int_t nTrks = trkArray->GetEntriesFast();
245   if(nTrks < fMinTracks) {
246     printf("TooFewTracks\n");
247     Double_t vtx[3]={0,0,0};
248     fVert.SetXYZ(vtx);
249     fVert.SetDispersion(999);
250     fVert.SetNContributors(-5);
251     return &fVert;
252   }
253   TTree *trkTree = new TTree("TreeT","tracks");
254   AliESDtrack *esdTrack = 0;
255   trkTree->Branch("tracks","AliESDtrack",&esdTrack);
256   for(Int_t i=0; i<nTrks; i++){
257     esdTrack = (AliESDtrack*)trkArray->At(i);
258     trkTree->Fill();
259   }
260     
261   AliVertex *vtx =  VertexForSelectedTracks(trkTree);
262   delete trkTree;
263   return vtx;
264 }
265
266 //---------------------------------------------------------------------------
267 void AliVertexerTracks::VertexFinder(Int_t optUseWeights) {
268
269   // Get estimate of vertex position in (x,y) from tracks DCA
270  
271   Double_t initPos[3];
272   initPos[2] = 0.;
273   for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
274   Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
275   Double_t aver[3]={0.,0.,0.};
276   Double_t aversq[3]={0.,0.,0.};
277   Double_t sigmasq[3]={0.,0.,0.};
278   Double_t sigma=0;
279   Int_t ncombi = 0;
280   AliESDtrack *track1;
281   AliESDtrack *track2;
282   Double_t pos[3],dir[3]; 
283   Double_t alpha,mindist;
284   Double_t field=GetField();
285
286   for(Int_t i=0; i<nacc; i++){
287     track1 = (AliESDtrack*)fTrkArray.At(i);
288     alpha=track1->GetAlpha();
289     mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
290     track1->GetXYZAt(mindist,field,pos);
291     track1->GetPxPyPzAt(mindist,field,dir);
292     AliStrLine *line1 = new AliStrLine(pos,dir); 
293
294    //    AliStrLine *line1 = new AliStrLine();
295    //    track1->ApproximateHelixWithLine(mindist,field,line1);
296    
297     for(Int_t j=i+1; j<nacc; j++){
298       track2 = (AliESDtrack*)fTrkArray.At(j);
299       alpha=track2->GetAlpha();
300       mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
301       track2->GetXYZAt(mindist,field,pos);
302       track2->GetPxPyPzAt(mindist,field,dir);
303       AliStrLine *line2 = new AliStrLine(pos,dir); 
304     //      AliStrLine *line2 = new AliStrLine();
305     //  track2->ApproximateHelixWithLine(mindist,field,line2);
306       Double_t distCA=line2->GetDCA(line1);
307       if(fDCAcut<=0 || (fDCAcut>0&&distCA<fDCAcut)){
308         Double_t pnt1[3],pnt2[3],crosspoint[3];
309
310         if(optUseWeights<=0){
311           Int_t retcode = line2->Cross(line1,crosspoint);
312           if(retcode>=0){
313             ncombi++;
314             for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
315             for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]);
316           }
317         }
318         if(optUseWeights>0){
319           Int_t retcode = line1->CrossPoints(line2,pnt1,pnt2);
320           if(retcode>=0){
321             Double_t alpha, cs, sn;
322             alpha=track1->GetAlpha();
323             cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);   
324             Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2();
325             alpha=track2->GetAlpha();
326             cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
327             Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2();
328             Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2();
329             Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
330             Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
331             Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
332             crosspoint[0]=wx1*pnt1[0] + wx2*pnt2[0]; 
333             crosspoint[1]=wy1*pnt1[1] + wy2*pnt2[1]; 
334             crosspoint[2]=wz1*pnt1[2] + wz2*pnt2[2];
335           
336             ncombi++;
337             for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
338             for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]);
339           }
340         }
341       }
342       delete line2;
343     }
344     delete line1;
345   }
346   if(ncombi>0){
347     for(Int_t jj=0;jj<3;jj++){
348       initPos[jj] = aver[jj]/ncombi;
349       aversq[jj]/=ncombi;
350       sigmasq[jj]=aversq[jj]-initPos[jj]*initPos[jj];
351       sigma+=sigmasq[jj];
352     }
353     sigma=TMath::Sqrt(TMath::Abs(sigma));
354   }
355   else {
356     Warning("VertexFinder","Finder did not succed");
357     sigma=999;
358   }
359   fVert.SetXYZ(initPos);
360   fVert.SetDispersion(sigma);
361   fVert.SetNContributors(ncombi);
362 }
363 //---------------------------------------------------------------------------
364 void AliVertexerTracks::HelixVertexFinder() {
365
366   // Get estimate of vertex position in (x,y) from tracks DCA
367
368
369   Double_t initPos[3];
370   initPos[2] = 0.;
371   for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
372   Double_t field=GetField();
373
374   Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
375
376   Double_t aver[3]={0.,0.,0.};
377   Double_t averquad[3]={0.,0.,0.};
378   Double_t sigmaquad[3]={0.,0.,0.};
379   Double_t sigma=0;
380   Int_t ncombi = 0;
381   AliESDtrack *track1;
382   AliESDtrack *track2;
383   Double_t distCA;
384   Double_t x, par[5];
385   Double_t alpha, cs, sn;
386   Double_t crosspoint[3];
387   for(Int_t i=0; i<nacc; i++){
388     track1 = (AliESDtrack*)fTrkArray.At(i);
389     
390
391     for(Int_t j=i+1; j<nacc; j++){
392       track2 = (AliESDtrack*)fTrkArray.At(j);
393
394       distCA=track2->PropagateToDCA(track1,field);
395
396       if(fDCAcut<=0 ||(fDCAcut>0&&distCA<fDCAcut)){
397         track1->GetExternalParameters(x,par);
398         alpha=track1->GetAlpha();
399         cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
400         Double_t x1=x*cs - par[0]*sn;
401         Double_t y1=x*sn + par[0]*cs;
402         Double_t z1=par[1];
403         Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2(); 
404         track2->GetExternalParameters(x,par);
405         alpha=track2->GetAlpha();
406         cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
407         Double_t x2=x*cs - par[0]*sn;
408         Double_t y2=x*sn + par[0]*cs;
409         Double_t z2=par[1];
410         Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2();
411         Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2();
412         Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
413         Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
414         Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
415         crosspoint[0]=wx1*x1 + wx2*x2; 
416         crosspoint[1]=wy1*y1 + wy2*y2; 
417         crosspoint[2]=wz1*z1 + wz2*z2;
418
419         ncombi++;
420         for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
421         for(Int_t jj=0;jj<3;jj++)averquad[jj]+=(crosspoint[jj]*crosspoint[jj]);
422       }
423     }
424       
425   }
426   if(ncombi>0){
427     for(Int_t jj=0;jj<3;jj++){
428       initPos[jj] = aver[jj]/ncombi;
429       averquad[jj]/=ncombi;
430       sigmaquad[jj]=averquad[jj]-initPos[jj]*initPos[jj];
431       sigma+=sigmaquad[jj];
432     }
433     sigma=TMath::Sqrt(TMath::Abs(sigma));
434   }
435   else {
436     Warning("HelixVertexFinder","Finder did not succed");
437     sigma=999;
438   }
439   fVert.SetXYZ(initPos);
440   fVert.SetDispersion(sigma);
441   fVert.SetNContributors(ncombi);
442 }
443 //---------------------------------------------------------------------------
444 void AliVertexerTracks::StrLinVertexFinderMinDist(Int_t optUseWeights){
445
446   // Calculate the point at minimum distance to prepared tracks 
447   
448   Double_t initPos[3];
449   initPos[2] = 0.;
450   Double_t sigma=0;
451   for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
452   const Int_t knacc = (Int_t)fTrkArray.GetEntriesFast();
453   Double_t field=GetField();
454
455   AliESDtrack *track1;
456   Double_t (*vectP0)[3]=new Double_t [knacc][3];
457   Double_t (*vectP1)[3]=new Double_t [knacc][3];
458   
459   Double_t sum[3][3];
460   Double_t dsum[3]={0,0,0};
461   for(Int_t i=0;i<3;i++)
462     for(Int_t j=0;j<3;j++)sum[i][j]=0;
463   for(Int_t i=0; i<knacc; i++){
464     track1 = (AliESDtrack*)fTrkArray.At(i);
465     Double_t alpha=track1->GetAlpha();
466     Double_t mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
467     Double_t pos[3],dir[3]; 
468     track1->GetXYZAt(mindist,field,pos);
469     track1->GetPxPyPzAt(mindist,field,dir);
470     AliStrLine *line1 = new AliStrLine(pos,dir); 
471     //    AliStrLine *line1 = new AliStrLine();
472     //    track1->ApproximateHelixWithLine(mindist,field,line1);
473
474     Double_t p0[3],cd[3];
475     line1->GetP0(p0);
476     line1->GetCd(cd);
477     Double_t p1[3]={p0[0]+cd[0],p0[1]+cd[1],p0[2]+cd[2]};
478     vectP0[i][0]=p0[0];
479     vectP0[i][1]=p0[1];
480     vectP0[i][2]=p0[2];
481     vectP1[i][0]=p1[0];
482     vectP1[i][1]=p1[1];
483     vectP1[i][2]=p1[2];
484     
485     Double_t matr[3][3];
486     Double_t dknow[3];
487     if(optUseWeights==0)GetStrLinDerivMatrix(p0,p1,matr,dknow);
488     if(optUseWeights==1){
489       Double_t sigmasq[3];
490       sigmasq[0]=track1->GetSigmaY2();
491       sigmasq[1]=track1->GetSigmaY2();
492       sigmasq[2]=track1->GetSigmaZ2();
493       GetStrLinDerivMatrix(p0,p1,sigmasq,matr,dknow);
494     }
495
496     for(Int_t iii=0;iii<3;iii++){
497       dsum[iii]+=dknow[iii]; 
498       for(Int_t lj=0;lj<3;lj++) sum[iii][lj]+=matr[iii][lj];
499     }
500     delete line1;
501   }
502   
503   Double_t vett[3][3];
504   Double_t det=GetDeterminant3X3(sum);
505   
506    if(det!=0){
507      for(Int_t zz=0;zz<3;zz++){
508        for(Int_t ww=0;ww<3;ww++){
509          for(Int_t kk=0;kk<3;kk++) vett[ww][kk]=sum[ww][kk];
510        }
511        for(Int_t kk=0;kk<3;kk++) vett[kk][zz]=dsum[kk];
512        initPos[zz]=GetDeterminant3X3(vett)/det;
513      }
514
515
516      for(Int_t i=0; i<knacc; i++){
517        Double_t p0[3]={0,0,0},p1[3]={0,0,0};
518        for(Int_t ii=0;ii<3;ii++){
519          p0[ii]=vectP0[i][ii];
520          p1[ii]=vectP1[i][ii];
521        }
522        sigma+=GetStrLinMinDist(p0,p1,initPos);
523      }
524
525      sigma=TMath::Sqrt(sigma);
526    }else{
527     Warning("StrLinVertexFinderMinDist","Finder did not succed");
528     sigma=999;
529   }
530   delete vectP0;
531   delete vectP1;
532   fVert.SetXYZ(initPos);
533   fVert.SetDispersion(sigma);
534   fVert.SetNContributors(knacc);
535 }
536 //_______________________________________________________________________
537 Double_t AliVertexerTracks::GetDeterminant3X3(Double_t matr[][3]){
538   //
539   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];
540  return det;
541 }
542 //____________________________________________________________________________
543 void AliVertexerTracks::GetStrLinDerivMatrix(Double_t *p0,Double_t *p1,Double_t m[][3],Double_t *d){
544
545   //
546   Double_t x12=p0[0]-p1[0];
547   Double_t y12=p0[1]-p1[1];
548   Double_t z12=p0[2]-p1[2];
549   Double_t kk=x12*x12+y12*y12+z12*z12;
550   m[0][0]=2-2/kk*x12*x12;
551   m[0][1]=-2/kk*x12*y12;
552   m[0][2]=-2/kk*x12*z12;
553   m[1][0]=-2/kk*x12*y12;
554   m[1][1]=2-2/kk*y12*y12;
555   m[1][2]=-2/kk*y12*z12;
556   m[2][0]=-2/kk*x12*z12;
557   m[2][1]=-2*y12*z12;
558   m[2][2]=2-2/kk*z12*z12;
559   d[0]=2*p0[0]-2/kk*p0[0]*x12*x12-2/kk*p0[2]*x12*z12-2/kk*p0[1]*x12*y12;
560   d[1]=2*p0[1]-2/kk*p0[1]*y12*y12-2/kk*p0[0]*x12*y12-2/kk*p0[2]*z12*y12;
561   d[2]=2*p0[2]-2/kk*p0[2]*z12*z12-2/kk*p0[0]*x12*z12-2/kk*p0[1]*z12*y12;
562
563 }
564 //____________________________________________________________________________
565 void AliVertexerTracks::GetStrLinDerivMatrix(Double_t *p0,Double_t *p1,Double_t *sigmasq,Double_t m[][3],Double_t *d){
566   //
567   Double_t x12=p1[0]-p0[0];
568   Double_t y12=p1[1]-p0[1];
569   Double_t z12=p1[2]-p0[2];
570
571   Double_t den= x12*x12*sigmasq[1]*sigmasq[2]+y12*y12*sigmasq[0]*sigmasq[2]+z12*z12*sigmasq[0]*sigmasq[1];
572
573   Double_t kk= 2*(x12*x12/sigmasq[0]+y12*y12/sigmasq[1]+z12*z12/sigmasq[2]);
574
575   Double_t cc[3];
576   cc[0]=-x12/sigmasq[0];
577   cc[1]=-y12/sigmasq[1];
578   cc[2]=-z12/sigmasq[2];
579
580   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;
581
582   Double_t ss= -p0[0]*cc[0]-p0[1]*cc[1]-p0[2]*cc[2];
583
584   Double_t aa[3];
585   aa[0]=x12*sigmasq[1]*sigmasq[2]/den;
586   aa[1]=y12*sigmasq[0]*sigmasq[2]/den;
587   aa[2]=z12*sigmasq[0]*sigmasq[1]/den;
588
589   m[0][0]=aa[0]*(aa[0]*kk+2*cc[0])+2*cc[0]*aa[0]+2/sigmasq[0];
590   m[0][1]=aa[1]*(aa[0]*kk+2*cc[0])+2*cc[1]*aa[0];
591   m[0][2]=aa[2]*(aa[0]*kk+2*cc[0])+2*cc[2]*aa[0];
592
593   m[1][0]=aa[0]*(aa[1]*kk+2*cc[1])+2*cc[0]*aa[1];
594   m[1][1]=aa[1]*(aa[1]*kk+2*cc[1])+2*cc[1]*aa[1]+2/sigmasq[1];
595   m[1][2]=aa[2]*(aa[1]*kk+2*cc[1])+2*cc[2]*aa[1];
596
597   m[2][0]=aa[0]*(aa[2]*kk+2*cc[2])+2*cc[0]*aa[2];
598   m[2][1]=aa[1]*(aa[2]*kk+2*cc[2])+2*cc[1]*aa[2];
599   m[2][2]=aa[2]*(aa[2]*kk+2*cc[2])+2*cc[2]*aa[2]+2/sigmasq[2];
600
601   d[0]=-ww*(aa[0]*kk+2*cc[0])-2*ss*aa[0]+2*p0[0]/sigmasq[0];
602   d[1]=-ww*(aa[1]*kk+2*cc[1])-2*ss*aa[1]+2*p0[1]/sigmasq[1];
603   d[2]=-ww*(aa[2]*kk+2*cc[2])-2*ss*aa[2]+2*p0[2]/sigmasq[2];
604
605 }
606 //_____________________________________________________________________________
607 Double_t AliVertexerTracks::GetStrLinMinDist(Double_t *p0,Double_t *p1,Double_t *x0){
608   //
609   Double_t x12=p0[0]-p1[0];
610   Double_t y12=p0[1]-p1[1];
611   Double_t z12=p0[2]-p1[2];
612   Double_t x10=p0[0]-x0[0];
613   Double_t y10=p0[1]-x0[1];
614   Double_t z10=p0[2]-x0[2];
615   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);
616 }
617
618 //---------------------------------------------------------------------------
619 void AliVertexerTracks::VertexFitter() {
620 //
621 // The optimal estimate of the vertex position is given by a "weighted 
622 // average of tracks positions"
623 // Original method: CMS Note 97/0051
624 //
625   Double_t initPos[3];
626   fVert.GetXYZ(initPos);
627   if(fDebug) { 
628     printf(" VertexFitter(): start\n");
629     printf(" Number of tracks in array: %d\n",(Int_t)fTrkArray.GetEntriesFast());
630     printf(" Minimum # tracks required in fit: %d\n",fMinTracks);
631     printf("Vertex position after finder: %f,%f,%f\n",initPos[0],initPos[1],initPos[2]);
632   }
633
634   Int_t i,j,k,step=0;
635   TMatrixD rv(3,1);
636   TMatrixD vV(3,3);
637   rv(0,0) = initPos[0];
638   rv(1,0) = initPos[1];
639   rv(2,0) = 0.;
640   Double_t xlStart,alpha;
641   Double_t rotAngle;
642   Double_t cosRot,sinRot;
643   Double_t cc[15];
644   Int_t nUsedTrks;
645   Double_t chi2,chi2i;
646   Int_t arrEntries = (Int_t)fTrkArray.GetEntries();
647   AliESDtrack *t = 0;
648   Int_t failed = 0;
649
650   Int_t *skipTrack = new Int_t[arrEntries];
651   for(i=0; i<arrEntries; i++) skipTrack[i]=0;
652
653   // 3 steps:
654   // 1st - first estimate of vtx using all tracks
655   // 2nd - apply cut on chi2 max per track
656   // 3rd - estimate of global chi2
657   for(step=0; step<3; step++) {
658     if(fDebug) printf(" step = %d\n",step);
659     chi2 = 0.;
660     nUsedTrks = 0;
661
662     TMatrixD sumWiri(3,1);
663     TMatrixD sumWi(3,3);
664     for(i=0; i<3; i++) {
665       sumWiri(i,0) = 0.;
666       for(j=0; j<3; j++) sumWi(j,i) = 0.;
667     }
668
669     // loop on tracks  
670     for(k=0; k<arrEntries; k++) {
671       if(skipTrack[k]) continue;
672       // get track from track array
673       t = (AliESDtrack*)fTrkArray.At(k);
674       alpha = t->GetAlpha();
675       xlStart = initPos[0]*TMath::Cos(alpha)+initPos[1]*TMath::Sin(alpha);
676       t->PropagateTo(xlStart,AliTracker::GetBz());   // to vtxSeed
677       rotAngle = alpha;
678       if(alpha<0.) rotAngle += 2.*TMath::Pi();
679       cosRot = TMath::Cos(rotAngle);
680       sinRot = TMath::Sin(rotAngle);
681       
682       // vector of track global coordinates
683       TMatrixD ri(3,1);
684       ri(0,0) = t->GetX()*cosRot-t->GetY()*sinRot;
685       ri(1,0) = t->GetX()*sinRot+t->GetY()*cosRot;
686       ri(2,0) = t->GetZ();
687
688       // matrix to go from global (x,y,z) to local (y,z);
689       TMatrixD qQi(2,3);
690       qQi(0,0) = -sinRot;
691       qQi(0,1) = cosRot;
692       qQi(0,2) = 0.;
693       qQi(1,0) = 0.;
694       qQi(1,1) = 0.;
695       qQi(1,2) = 1.;
696
697       // covariance matrix of local (y,z) - inverted
698       TMatrixD uUi(2,2);
699       t->GetExternalCovariance(cc);
700       uUi(0,0) = cc[0];
701       uUi(0,1) = cc[1];
702       uUi(1,0) = cc[1];
703       uUi(1,1) = cc[2];
704
705       // weights matrix: wWi = qQiT * uUiInv * qQi
706       if(uUi.Determinant() <= 0.) continue;
707       TMatrixD uUiInv(TMatrixD::kInverted,uUi);
708       TMatrixD uUiInvQi(uUiInv,TMatrixD::kMult,qQi);
709       TMatrixD wWi(qQi,TMatrixD::kTransposeMult,uUiInvQi);
710
711       // track chi2
712       TMatrixD deltar = rv; deltar -= ri;
713       TMatrixD wWideltar(wWi,TMatrixD::kMult,deltar);
714       chi2i = deltar(0,0)*wWideltar(0,0)+
715               deltar(1,0)*wWideltar(1,0)+
716               deltar(2,0)*wWideltar(2,0);
717
718
719       if(step==1 && chi2i > fMaxChi2PerTrack) {
720         skipTrack[k] = 1;
721         continue;
722       }
723
724       // add to total chi2
725       chi2 += chi2i;
726
727       TMatrixD wWiri(wWi,TMatrixD::kMult,ri); 
728
729       sumWiri += wWiri;
730       sumWi   += wWi;
731
732       nUsedTrks++;
733     } // end loop on tracks
734
735     if(nUsedTrks < fMinTracks) {
736       failed=1;
737       continue;
738     }
739
740     Double_t determinant = sumWi.Determinant();
741     //cerr<<" determinant: "<<determinant<<endl;
742     if(determinant < 100.)  { 
743       printf("det(V) = 0\n");       
744       failed=1;
745       continue;
746     }
747
748     // inverted of weights matrix
749     TMatrixD invsumWi(TMatrixD::kInverted,sumWi);
750     vV = invsumWi;
751      
752     // position of primary vertex
753     rv.Mult(vV,sumWiri);
754
755   } // end loop on the 3 steps
756
757   delete [] skipTrack;
758   delete t;
759
760   if(failed) { 
761     printf("TooFewTracks\n");
762     fCurrentVertex = new AliESDVertex(0.,0.,-1);
763     return; 
764   }
765
766   Double_t position[3];
767   position[0] = rv(0,0);
768   position[1] = rv(1,0);
769   position[2] = rv(2,0);
770   Double_t covmatrix[6];
771   covmatrix[0] = vV(0,0);
772   covmatrix[1] = vV(0,1);
773   covmatrix[2] = vV(1,1);
774   covmatrix[3] = vV(0,2);
775   covmatrix[4] = vV(1,2);
776   covmatrix[5] = vV(2,2);
777   
778   // store data in the vertex object
779   fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nUsedTrks);
780
781   if(fDebug) {
782     printf(" VertexFitter(): finish\n");
783     printf(" rv = ( %f , %f , %f )\n\n",rv(0,0),rv(1,0),rv(2,0));
784     fCurrentVertex->PrintStatus();
785   }
786
787   return;
788 }