]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliVertexerTracks.cxx
Calling the V0 and cascade finders from AliReconstruction, and the possibility to...
[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 fMinTracks(2),
46 fMinITSClusters(5),
47 fTrkArray(),
48 fTrksToSkip(0),
49 fNTrksToSkip(0),
50 fDCAcut(0),
51 fAlgo(1),
52 fNSigma(3),
53 fITSin(kTRUE),
54 fITSrefit(kTRUE),
55 fDebug(0)
56 {
57 //
58 // Default constructor
59 //
60   SetVtxStart();
61   SetVtxStartSigma();
62   SetMinTracks();
63   SetMinITSClusters();
64   SetNSigmad0(); 
65 }
66 //-----------------------------------------------------------------------------
67 AliVertexerTracks::AliVertexerTracks(Double_t xStart, Double_t yStart):
68 TObject(),
69 fVert(),
70 fCurrentVertex(0),
71 fMinTracks(2),
72 fMinITSClusters(5),
73 fTrkArray(),
74 fTrksToSkip(0),
75 fNTrksToSkip(0),
76 fDCAcut(0),
77 fAlgo(1),
78 fNSigma(3),
79 fITSin(kTRUE),
80 fITSrefit(kTRUE),
81 fDebug(0)
82 {
83 //
84 // Standard constructor
85 //
86   SetVtxStart(xStart,yStart);
87   SetVtxStartSigma();
88   SetMinTracks();
89   SetMinITSClusters();
90   SetNSigmad0(); 
91 }
92 //-----------------------------------------------------------------------------
93 AliVertexerTracks::~AliVertexerTracks() {
94   // Default Destructor
95   // The objects pointed by the following pointer are not owned
96   // by this class and are not deleted
97   fCurrentVertex = 0;
98   delete[] fTrksToSkip;
99 }
100 //----------------------------------------------------------------------------
101 AliESDVertex* AliVertexerTracks::FindPrimaryVertex(const AliESD *esdEvent)
102 {
103 //
104 // Primary vertex for current ESD event
105 // (Two iterations: 
106 //  1st with 5*fNSigma*sigma(pt) cut w.r.t. to initial vertex; 
107 //  2nd with fNSigma*sigma(pt) cut w.r.t. to vertex found in 1st iteration) 
108 // All ESD tracks with inside the beam pipe are then propagated to found vertex
109 //
110   fCurrentVertex = 0;
111
112   Int_t nTrksTot = (Int_t)esdEvent->GetNumberOfTracks();
113   TTree *trkTree = new TTree("TreeT","tracks");
114   AliESDtrack *esdTrack = 0;
115   trkTree->Branch("tracks","AliESDtrack",&esdTrack);
116
117   Bool_t   skipThis;
118   for(Int_t i=0; i<nTrksTot; i++) {
119     // check tracks to skip
120     skipThis = kFALSE;
121     for(Int_t j=0; j<fNTrksToSkip; j++) { 
122       if(i==fTrksToSkip[j]) {
123         if(fDebug) printf("skipping track: %d\n",i);
124         skipThis = kTRUE;
125       }
126     }
127     AliESDtrack *et = esdEvent->GetTrack(i);
128     esdTrack = new AliESDtrack(*et);
129     if(skipThis) {delete esdTrack;continue;}
130     if(fITSin) {
131       if(!(esdTrack->GetStatus()&AliESDtrack::kITSin)) {delete esdTrack;continue;}
132       if(fITSrefit && !(esdTrack->GetStatus()&AliESDtrack::kITSrefit)) {delete esdTrack;continue;}
133       Int_t nclus=esdTrack->GetNcls(0); // check number of clusters in ITS
134       if(nclus<fMinITSClusters) {delete esdTrack;continue;}
135     }
136     trkTree->Fill();
137     delete esdTrack;
138   }
139  
140
141   // ITERATION 1
142   // propagate tracks to initVertex
143   // preselect them  (reject for |d0|>5*fNSigma*sigma w.r.t. initVertex)
144   Int_t nTrksPrep;
145   nTrksPrep = PrepareTracks(*trkTree,1);
146   if(fDebug) printf(" tracks prepared - step 1: %d\n",nTrksPrep);
147   if(nTrksPrep < fMinTracks) {
148     if(fDebug) printf("TooFewTracks\n");
149     TooFewTracks(esdEvent);
150     if(fDebug) fCurrentVertex->PrintStatus();
151     fTrkArray.Delete();
152     delete trkTree;
153     return fCurrentVertex; 
154   }
155
156   // vertex finder
157   switch (fAlgo) {
158     case 1: StrLinVertexFinderMinDist(1); break;
159     case 2: StrLinVertexFinderMinDist(0); break;
160     case 3: HelixVertexFinder();          break;
161     case 4: VertexFinder(1);              break;
162     case 5: VertexFinder(0);              break;
163     default: printf("Wrong algorithm\n"); break;  
164   }
165   if(fDebug) printf(" vertex finding completed\n");
166
167   // vertex fitter
168   VertexFitter(kTRUE);
169   if(fDebug) printf(" vertex fit completed\n");
170   fTrkArray.Delete();
171   // ITERATION 2
172   // propagate tracks to best between initVertex and fCurrentVertex
173   // preselect tracks (reject for |d0|>fNSigma*sigma w.r.t. best 
174   //                   between initVertex and fCurrentVertex) 
175   nTrksPrep = PrepareTracks(*trkTree,2);
176   delete trkTree;
177   if(fDebug) printf(" tracks prepared - step 2: %d\n",nTrksPrep);
178   if(nTrksPrep < fMinTracks) {
179     if(fDebug) printf("TooFewTracks\n");
180     TooFewTracks(esdEvent);
181     if(fDebug) fCurrentVertex->PrintStatus();
182     fTrkArray.Delete();
183     return fCurrentVertex; 
184   }
185
186   // vertex finder
187   switch (fAlgo) {
188     case 1: StrLinVertexFinderMinDist(1); break;
189     case 2: StrLinVertexFinderMinDist(0); break;
190     case 3: HelixVertexFinder();          break;
191     case 4: VertexFinder(1);              break;
192     case 5: VertexFinder(0);              break;
193     default: printf("Wrong algorithm\n"); break;  
194   }
195   if(fDebug) printf(" vertex finding completed\n");
196
197   // fitter
198   VertexFitter(kTRUE);
199   if(fDebug) printf(" vertex fit completed\n");
200
201
202   // take true pos from ESD
203   Double_t tp[3];
204   esdEvent->GetVertex()->GetTruePos(tp);
205   fCurrentVertex->SetTruePos(tp);
206   if(fNominalCov[0]>1.) {
207     fCurrentVertex->SetTitle("VertexerTracksNoConstraint");
208   } else {
209     fCurrentVertex->SetTitle("VertexerTracksWithConstraint");
210   }
211
212   
213   // propagate tracks to found vertex
214   if(TMath::Sqrt(fCurrentVertex->GetXv()*fCurrentVertex->GetXv()+fCurrentVertex->GetYv()*fCurrentVertex->GetYv())<3.) {
215     for(Int_t ii=0; ii<nTrksTot; ii++) {
216       AliESDtrack *et = esdEvent->GetTrack(ii);
217       if(!(et->GetStatus()&AliESDtrack::kITSin)) continue;
218       if(et->GetX()>3.) continue;
219       et->RelateToVertex(fCurrentVertex,GetField(),100.);
220     }
221   } else {
222     AliWarning("Found vertex outside beam pipe!");
223   }
224
225   // set indices of used tracks
226   UShort_t *indices = 0;
227   AliESDtrack *ett = 0;
228   if(fCurrentVertex->GetNContributors()>0) {
229     indices = new UShort_t[fCurrentVertex->GetNContributors()];
230     for(Int_t jj=0;jj<(Int_t)fTrkArray.GetEntriesFast();jj++) {
231       ett = (AliESDtrack*)fTrkArray.At(jj);
232       indices[jj] = (UShort_t)ett->GetID();
233     }
234     fCurrentVertex->SetIndices(fCurrentVertex->GetNContributors(),indices);
235   }
236   delete [] indices;
237   fTrkArray.Delete();
238
239   if(fTrksToSkip) delete [] fTrksToSkip;
240
241
242   if(fDebug) fCurrentVertex->PrintStatus();
243   if(fDebug) fCurrentVertex->PrintIndices();
244  
245   return fCurrentVertex;
246 }
247 //------------------------------------------------------------------------
248 Double_t AliVertexerTracks::GetDeterminant3X3(Double_t matr[][3]){
249   //
250   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];
251  return det;
252 }
253 //-------------------------------------------------------------------------
254 void AliVertexerTracks::GetStrLinDerivMatrix(Double_t *p0,Double_t *p1,Double_t (*m)[3],Double_t *d){
255
256   //
257   Double_t x12=p0[0]-p1[0];
258   Double_t y12=p0[1]-p1[1];
259   Double_t z12=p0[2]-p1[2];
260   Double_t kk=x12*x12+y12*y12+z12*z12;
261   m[0][0]=2-2/kk*x12*x12;
262   m[0][1]=-2/kk*x12*y12;
263   m[0][2]=-2/kk*x12*z12;
264   m[1][0]=-2/kk*x12*y12;
265   m[1][1]=2-2/kk*y12*y12;
266   m[1][2]=-2/kk*y12*z12;
267   m[2][0]=-2/kk*x12*z12;
268   m[2][1]=-2*y12*z12;
269   m[2][2]=2-2/kk*z12*z12;
270   d[0]=2*p0[0]-2/kk*p0[0]*x12*x12-2/kk*p0[2]*x12*z12-2/kk*p0[1]*x12*y12;
271   d[1]=2*p0[1]-2/kk*p0[1]*y12*y12-2/kk*p0[0]*x12*y12-2/kk*p0[2]*z12*y12;
272   d[2]=2*p0[2]-2/kk*p0[2]*z12*z12-2/kk*p0[0]*x12*z12-2/kk*p0[1]*z12*y12;
273
274 }
275 //--------------------------------------------------------------------------  
276 void AliVertexerTracks::GetStrLinDerivMatrix(Double_t *p0,Double_t *p1,Double_t *sigmasq,Double_t (*m)[3],Double_t *d){
277   //
278   Double_t x12=p1[0]-p0[0];
279   Double_t y12=p1[1]-p0[1];
280   Double_t z12=p1[2]-p0[2];
281
282   Double_t den= x12*x12*sigmasq[1]*sigmasq[2]+y12*y12*sigmasq[0]*sigmasq[2]+z12*z12*sigmasq[0]*sigmasq[1];
283
284   Double_t kk= 2*(x12*x12/sigmasq[0]+y12*y12/sigmasq[1]+z12*z12/sigmasq[2]);
285
286   Double_t cc[3];
287   cc[0]=-x12/sigmasq[0];
288   cc[1]=-y12/sigmasq[1];
289   cc[2]=-z12/sigmasq[2];
290
291   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;
292
293   Double_t ss= -p0[0]*cc[0]-p0[1]*cc[1]-p0[2]*cc[2];
294
295   Double_t aa[3];
296   aa[0]=x12*sigmasq[1]*sigmasq[2]/den;
297   aa[1]=y12*sigmasq[0]*sigmasq[2]/den;
298   aa[2]=z12*sigmasq[0]*sigmasq[1]/den;
299
300   m[0][0]=aa[0]*(aa[0]*kk+2*cc[0])+2*cc[0]*aa[0]+2/sigmasq[0];
301   m[0][1]=aa[1]*(aa[0]*kk+2*cc[0])+2*cc[1]*aa[0];
302   m[0][2]=aa[2]*(aa[0]*kk+2*cc[0])+2*cc[2]*aa[0];
303
304   m[1][0]=aa[0]*(aa[1]*kk+2*cc[1])+2*cc[0]*aa[1];
305   m[1][1]=aa[1]*(aa[1]*kk+2*cc[1])+2*cc[1]*aa[1]+2/sigmasq[1];
306   m[1][2]=aa[2]*(aa[1]*kk+2*cc[1])+2*cc[2]*aa[1];
307
308   m[2][0]=aa[0]*(aa[2]*kk+2*cc[2])+2*cc[0]*aa[2];
309   m[2][1]=aa[1]*(aa[2]*kk+2*cc[2])+2*cc[1]*aa[2];
310   m[2][2]=aa[2]*(aa[2]*kk+2*cc[2])+2*cc[2]*aa[2]+2/sigmasq[2];
311
312   d[0]=-ww*(aa[0]*kk+2*cc[0])-2*ss*aa[0]+2*p0[0]/sigmasq[0];
313   d[1]=-ww*(aa[1]*kk+2*cc[1])-2*ss*aa[1]+2*p0[1]/sigmasq[1];
314   d[2]=-ww*(aa[2]*kk+2*cc[2])-2*ss*aa[2]+2*p0[2]/sigmasq[2];
315
316   }
317 //--------------------------------------------------------------------------   
318 Double_t AliVertexerTracks::GetStrLinMinDist(Double_t *p0,Double_t *p1,Double_t *x0){
319   //
320   Double_t x12=p0[0]-p1[0];
321   Double_t y12=p0[1]-p1[1];
322   Double_t z12=p0[2]-p1[2];
323   Double_t x10=p0[0]-x0[0];
324   Double_t y10=p0[1]-x0[1];
325   Double_t z10=p0[2]-x0[2];
326   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);
327 }
328 //---------------------------------------------------------------------------
329 void AliVertexerTracks::HelixVertexFinder() {
330
331   // Get estimate of vertex position in (x,y) from tracks DCA
332
333
334   Double_t initPos[3];
335   initPos[2] = 0.;
336   for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
337   Double_t field=GetField();
338
339   Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
340
341   Double_t aver[3]={0.,0.,0.};
342   Double_t averquad[3]={0.,0.,0.};
343   Double_t sigmaquad[3]={0.,0.,0.};
344   Double_t sigma=0;
345   Int_t ncombi = 0;
346   AliESDtrack *track1;
347   AliESDtrack *track2;
348   Double_t distCA;
349   Double_t x, par[5];
350   Double_t alpha, cs, sn;
351   Double_t crosspoint[3];
352   for(Int_t i=0; i<nacc; i++){
353     track1 = (AliESDtrack*)fTrkArray.At(i);
354     
355
356     for(Int_t j=i+1; j<nacc; j++){
357       track2 = (AliESDtrack*)fTrkArray.At(j);
358
359       distCA=track2->PropagateToDCA(track1,field);
360
361       if(fDCAcut<=0 ||(fDCAcut>0&&distCA<fDCAcut)){
362         track1->GetExternalParameters(x,par);
363         alpha=track1->GetAlpha();
364         cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
365         Double_t x1=x*cs - par[0]*sn;
366         Double_t y1=x*sn + par[0]*cs;
367         Double_t z1=par[1];
368         Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2(); 
369         track2->GetExternalParameters(x,par);
370         alpha=track2->GetAlpha();
371         cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
372         Double_t x2=x*cs - par[0]*sn;
373         Double_t y2=x*sn + par[0]*cs;
374         Double_t z2=par[1];
375         Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2();
376         Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2();
377         Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
378         Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
379         Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
380         crosspoint[0]=wx1*x1 + wx2*x2; 
381         crosspoint[1]=wy1*y1 + wy2*y2; 
382         crosspoint[2]=wz1*z1 + wz2*z2;
383
384         ncombi++;
385         for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
386         for(Int_t jj=0;jj<3;jj++)averquad[jj]+=(crosspoint[jj]*crosspoint[jj]);
387       }
388     }
389       
390   }
391   if(ncombi>0){
392     for(Int_t jj=0;jj<3;jj++){
393       initPos[jj] = aver[jj]/ncombi;
394       averquad[jj]/=ncombi;
395       sigmaquad[jj]=averquad[jj]-initPos[jj]*initPos[jj];
396       sigma+=sigmaquad[jj];
397     }
398     sigma=TMath::Sqrt(TMath::Abs(sigma));
399   }
400   else {
401     Warning("HelixVertexFinder","Finder did not succed");
402     sigma=999;
403   }
404   fVert.SetXYZ(initPos);
405   fVert.SetDispersion(sigma);
406   fVert.SetNContributors(ncombi);
407 }
408 //----------------------------------------------------------------------------
409 Int_t AliVertexerTracks::PrepareTracks(TTree &trkTree,Int_t OptImpParCut) {
410 //
411 // Propagate tracks to initial vertex position and store them in a TObjArray
412 //
413   Double_t maxd0rphi = 3.;  
414   Int_t    nTrks    = 0;
415   Double_t sigmaCurr[3];
416   Float_t  d0z0[2],covd0z0[3]; 
417   Double_t sigma;
418   Double_t field=GetField();
419
420   AliESDVertex *initVertex = new AliESDVertex(fNominalPos,fNominalCov,1,1);
421
422   Int_t    nEntries = (Int_t)trkTree.GetEntries();
423   if(!fTrkArray.IsEmpty()) fTrkArray.Delete();
424
425   if(fDebug) {
426     printf(" PrepareTracks()\n");
427   }
428
429   for(Int_t i=0; i<nEntries; i++) {
430     AliESDtrack *track = new AliESDtrack; 
431     trkTree.SetBranchAddress("tracks",&track);
432     trkTree.GetEvent(i);
433
434     // propagate track to vertex
435     if(OptImpParCut==1) { // OptImpParCut==1
436       track->RelateToVertex(initVertex,field,100.);
437     } else {              // OptImpParCut==2
438       fCurrentVertex->GetSigmaXYZ(sigmaCurr);
439       if((sigmaCurr[0]+sigmaCurr[1])<(TMath::Sqrt(fNominalCov[0])+TMath::Sqrt(fNominalCov[2]))) {
440         track->RelateToVertex(fCurrentVertex,field,100.);
441       } else {
442         track->RelateToVertex(initVertex,field,100.);
443       }
444     }
445
446     // select tracks with d0rphi < maxd0rphi
447     track->GetImpactParameters(d0z0,covd0z0);
448     sigma = TMath::Sqrt(covd0z0[0]);
449     maxd0rphi = fNSigma*sigma;
450     if(OptImpParCut==1) maxd0rphi *= 5.;
451
452     if(fDebug) printf("trk %d; lab %d; |d0| = %f;  cut = %f\n",i,track->GetLabel(),TMath::Abs(d0z0[0]),maxd0rphi);
453     if(TMath::Abs(d0z0[0]) > maxd0rphi) { 
454       if(fDebug) printf("     rejected\n");
455       delete track; continue; 
456     }
457
458     fTrkArray.AddLast(track);
459     nTrks++; 
460   }
461
462   delete initVertex;
463
464   return nTrks;
465
466 //---------------------------------------------------------------------------
467 void AliVertexerTracks::SetSkipTracks(Int_t n,Int_t *skipped) {
468 //
469 // Mark the tracks not ot be used in the vertex finding
470 //
471   fNTrksToSkip = n;
472   fTrksToSkip = new Int_t[n]; 
473   for(Int_t i=0;i<n;i++) fTrksToSkip[i] = skipped[i]; 
474   return; 
475 }
476 //---------------------------------------------------------------------------
477 void  AliVertexerTracks::SetVtxStart(AliESDVertex *vtx) { 
478 //
479 // Set initial vertex knowledge
480 //
481   vtx->GetXYZ(fNominalPos);
482   vtx->GetCovMatrix(fNominalCov);
483   return; 
484 }
485 //---------------------------------------------------------------------------
486 void AliVertexerTracks::StrLinVertexFinderMinDist(Int_t optUseWeights){
487
488   // Calculate the point at minimum distance to prepared tracks 
489   
490   Double_t initPos[3];
491   initPos[2] = 0.;
492   Double_t sigma=0;
493   for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
494   const Int_t knacc = (Int_t)fTrkArray.GetEntriesFast();
495   Double_t field=GetField();
496
497   AliESDtrack *track1;
498   Double_t (*vectP0)[3]=new Double_t [knacc][3];
499   Double_t (*vectP1)[3]=new Double_t [knacc][3];
500   
501   Double_t sum[3][3];
502   Double_t dsum[3]={0,0,0};
503   for(Int_t i=0;i<3;i++)
504     for(Int_t j=0;j<3;j++)sum[i][j]=0;
505   for(Int_t i=0; i<knacc; i++){
506     track1 = (AliESDtrack*)fTrkArray.At(i);
507     Double_t alpha=track1->GetAlpha();
508     Double_t mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
509     Double_t pos[3],dir[3]; 
510     track1->GetXYZAt(mindist,field,pos);
511     track1->GetPxPyPzAt(mindist,field,dir);
512     AliStrLine *line1 = new AliStrLine(pos,dir); 
513     //    AliStrLine *line1 = new AliStrLine();
514     //    track1->ApproximateHelixWithLine(mindist,field,line1);
515
516     Double_t p0[3],cd[3];
517     line1->GetP0(p0);
518     line1->GetCd(cd);
519     Double_t p1[3]={p0[0]+cd[0],p0[1]+cd[1],p0[2]+cd[2]};
520     vectP0[i][0]=p0[0];
521     vectP0[i][1]=p0[1];
522     vectP0[i][2]=p0[2];
523     vectP1[i][0]=p1[0];
524     vectP1[i][1]=p1[1];
525     vectP1[i][2]=p1[2];
526     
527     Double_t matr[3][3];
528     Double_t dknow[3];
529     if(optUseWeights==0)GetStrLinDerivMatrix(p0,p1,matr,dknow);
530     if(optUseWeights==1){
531       Double_t sigmasq[3];
532       sigmasq[0]=track1->GetSigmaY2();
533       sigmasq[1]=track1->GetSigmaY2();
534       sigmasq[2]=track1->GetSigmaZ2();
535       GetStrLinDerivMatrix(p0,p1,sigmasq,matr,dknow);
536     }
537
538     for(Int_t iii=0;iii<3;iii++){
539       dsum[iii]+=dknow[iii]; 
540       for(Int_t lj=0;lj<3;lj++) sum[iii][lj]+=matr[iii][lj];
541     }
542     delete line1;
543   }
544   
545   Double_t vett[3][3];
546   Double_t det=GetDeterminant3X3(sum);
547   
548    if(det!=0){
549      for(Int_t zz=0;zz<3;zz++){
550        for(Int_t ww=0;ww<3;ww++){
551          for(Int_t kk=0;kk<3;kk++) vett[ww][kk]=sum[ww][kk];
552        }
553        for(Int_t kk=0;kk<3;kk++) vett[kk][zz]=dsum[kk];
554        initPos[zz]=GetDeterminant3X3(vett)/det;
555      }
556
557
558      for(Int_t i=0; i<knacc; i++){
559        Double_t p0[3]={0,0,0},p1[3]={0,0,0};
560        for(Int_t ii=0;ii<3;ii++){
561          p0[ii]=vectP0[i][ii];
562          p1[ii]=vectP1[i][ii];
563        }
564        sigma+=GetStrLinMinDist(p0,p1,initPos);
565      }
566
567      sigma=TMath::Sqrt(sigma);
568    }else{
569     Warning("StrLinVertexFinderMinDist","Finder did not succed");
570     sigma=999;
571   }
572   delete vectP0;
573   delete vectP1;
574   fVert.SetXYZ(initPos);
575   fVert.SetDispersion(sigma);
576   fVert.SetNContributors(knacc);
577 }
578 //---------------------------------------------------------------------------
579 void AliVertexerTracks::TooFewTracks(const AliESD* esdEvent) {
580 //
581 // When the number of tracks is < fMinTracks
582 //
583
584   // deal with vertices not found
585   Double_t pos[3],err[3];
586   Int_t    ncontr=0;
587   pos[0] = fNominalPos[0];
588   err[0] = TMath::Sqrt(fNominalCov[0]);
589   pos[1] = fNominalPos[1];
590   err[1] = TMath::Sqrt(fNominalCov[2]);
591   pos[2] = esdEvent->GetVertex()->GetZv();
592   err[2] = esdEvent->GetVertex()->GetZRes();
593   if(err[0]>1. && esdEvent->GetVertex()->GetNContributors()<=0) 
594     ncontr = -1; // (x,y,z) = (0,0,0) 
595   if(err[0]>1. && esdEvent->GetVertex()->GetNContributors()>0) 
596     ncontr = -2; // (x,y,z) = (0,0,z_fromSPD) 
597   if(err[0]<1. && esdEvent->GetVertex()->GetNContributors()<=0) 
598     ncontr = -3; // (x,y,z) = (x_mean,y_mean,0)
599   if(err[0]<1. && esdEvent->GetVertex()->GetNContributors()>0) 
600     ncontr = -4; // (x,y,z) = (x_mean,y_mean,z_fromSPD)
601   fCurrentVertex = 0;
602   fCurrentVertex = new AliESDVertex(pos,err);
603   fCurrentVertex->SetNContributors(ncontr);
604
605   Double_t tp[3];
606   esdEvent->GetVertex()->GetTruePos(tp);
607   fCurrentVertex->SetTruePos(tp);
608   fCurrentVertex->SetTitle("VertexerTracksWithConstraint");
609   if(ncontr==-1||ncontr==-2) 
610     fCurrentVertex->SetTitle("VertexerTracksNoConstraint");
611
612   return;
613 }
614 //---------------------------------------------------------------------------
615 void AliVertexerTracks::VertexFinder(Int_t optUseWeights) {
616
617   // Get estimate of vertex position in (x,y) from tracks DCA
618  
619   Double_t initPos[3];
620   initPos[2] = 0.;
621   for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
622   Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
623   Double_t aver[3]={0.,0.,0.};
624   Double_t aversq[3]={0.,0.,0.};
625   Double_t sigmasq[3]={0.,0.,0.};
626   Double_t sigma=0;
627   Int_t ncombi = 0;
628   AliESDtrack *track1;
629   AliESDtrack *track2;
630   Double_t pos[3],dir[3]; 
631   Double_t alpha,mindist;
632   Double_t field=GetField();
633
634   for(Int_t i=0; i<nacc; i++){
635     track1 = (AliESDtrack*)fTrkArray.At(i);
636     alpha=track1->GetAlpha();
637     mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
638     track1->GetXYZAt(mindist,field,pos);
639     track1->GetPxPyPzAt(mindist,field,dir);
640     AliStrLine *line1 = new AliStrLine(pos,dir); 
641
642    //    AliStrLine *line1 = new AliStrLine();
643    //    track1->ApproximateHelixWithLine(mindist,field,line1);
644    
645     for(Int_t j=i+1; j<nacc; j++){
646       track2 = (AliESDtrack*)fTrkArray.At(j);
647       alpha=track2->GetAlpha();
648       mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
649       track2->GetXYZAt(mindist,field,pos);
650       track2->GetPxPyPzAt(mindist,field,dir);
651       AliStrLine *line2 = new AliStrLine(pos,dir); 
652     //      AliStrLine *line2 = new AliStrLine();
653     //  track2->ApproximateHelixWithLine(mindist,field,line2);
654       Double_t distCA=line2->GetDCA(line1);
655       if(fDCAcut<=0 || (fDCAcut>0&&distCA<fDCAcut)){
656         Double_t pnt1[3],pnt2[3],crosspoint[3];
657
658         if(optUseWeights<=0){
659           Int_t retcode = line2->Cross(line1,crosspoint);
660           if(retcode>=0){
661             ncombi++;
662             for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
663             for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]);
664           }
665         }
666         if(optUseWeights>0){
667           Int_t retcode = line1->CrossPoints(line2,pnt1,pnt2);
668           if(retcode>=0){
669             Double_t alpha, cs, sn;
670             alpha=track1->GetAlpha();
671             cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);   
672             Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2();
673             alpha=track2->GetAlpha();
674             cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
675             Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2();
676             Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2();
677             Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
678             Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
679             Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
680             crosspoint[0]=wx1*pnt1[0] + wx2*pnt2[0]; 
681             crosspoint[1]=wy1*pnt1[1] + wy2*pnt2[1]; 
682             crosspoint[2]=wz1*pnt1[2] + wz2*pnt2[2];
683           
684             ncombi++;
685             for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
686             for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]);
687           }
688         }
689       }
690       delete line2;
691     }
692     delete line1;
693   }
694   if(ncombi>0){
695     for(Int_t jj=0;jj<3;jj++){
696       initPos[jj] = aver[jj]/ncombi;
697       aversq[jj]/=ncombi;
698       sigmasq[jj]=aversq[jj]-initPos[jj]*initPos[jj];
699       sigma+=sigmasq[jj];
700     }
701     sigma=TMath::Sqrt(TMath::Abs(sigma));
702   }
703   else {
704     Warning("VertexFinder","Finder did not succed");
705     sigma=999;
706   }
707   fVert.SetXYZ(initPos);
708   fVert.SetDispersion(sigma);
709   fVert.SetNContributors(ncombi);
710 }
711 //---------------------------------------------------------------------------
712 void AliVertexerTracks::VertexFitter(Bool_t useNominalVtx) {
713 //
714 // The optimal estimate of the vertex position is given by a "weighted 
715 // average of tracks positions"
716 // Original method: CMS Note 97/0051
717 //
718   Double_t initPos[3];
719   fVert.GetXYZ(initPos);
720   if(fDebug) { 
721     printf(" VertexFitter(): start\n");
722     printf(" Number of tracks in array: %d\n",(Int_t)fTrkArray.GetEntriesFast());
723     printf(" Minimum # tracks required in fit: %d\n",fMinTracks);
724     printf("Vertex position after finder: %f,%f,%f\n",initPos[0],initPos[1],initPos[2]);
725     if(useNominalVtx) 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])); 
726   }
727
728   Int_t i,j,k,step=0;
729   TMatrixD rv(3,1);
730   TMatrixD vV(3,3);
731   rv(0,0) = initPos[0];
732   rv(1,0) = initPos[1];
733   rv(2,0) = 0.;
734   Double_t xlStart,alpha;
735   Double_t rotAngle;
736   Double_t cosRot,sinRot;
737   Double_t cc[15];
738   Int_t nUsedTrks;
739   Double_t chi2,chi2i;
740   Int_t arrEntries = (Int_t)fTrkArray.GetEntries();
741   AliESDtrack *t = 0;
742   Int_t failed = 0;
743
744   // initial vertex covariance matrix
745   TMatrixD vVb(3,3);
746   vVb(0,0) = fNominalCov[0];
747   vVb(0,1) = fNominalCov[1];
748   vVb(0,2) = 0.;
749   vVb(1,0) = fNominalCov[1];
750   vVb(1,1) = fNominalCov[2];
751   vVb(1,2) = 0.;
752   vVb(2,0) = 0.;
753   vVb(2,1) = 0.;
754   vVb(2,2) = fNominalCov[5];
755   TMatrixD vVbInv(TMatrixD::kInverted,vVb);
756   TMatrixD rb(3,1);
757   rb(0,0) = fNominalPos[0];
758   rb(1,0) = fNominalPos[1];
759   rb(2,0) = fNominalPos[2];
760   TMatrixD vVbInvrb(vVbInv,TMatrixD::kMult,rb);
761
762
763   // 2 steps:
764   // 1st - estimate of vtx using all tracks
765   // 2nd - estimate of global chi2
766   for(step=0; step<2; step++) {
767     if(fDebug) printf(" step = %d\n",step);
768     chi2 = 0.;
769     nUsedTrks = 0;
770
771     TMatrixD sumWiri(3,1);
772     TMatrixD sumWi(3,3);
773     for(i=0; i<3; i++) {
774       sumWiri(i,0) = 0.;
775       for(j=0; j<3; j++) sumWi(j,i) = 0.;
776     }
777
778
779     if(useNominalVtx) {
780       for(i=0;i<3;i++) {
781         sumWiri(i,0) += vVbInvrb(i,0);
782         for(k=0;k<3;k++) sumWi(i,k) += vVbInv(i,k);
783       }
784     }
785
786
787     // loop on tracks  
788     for(k=0; k<arrEntries; k++) {
789       // get track from track array
790       t = (AliESDtrack*)fTrkArray.At(k);
791       alpha = t->GetAlpha();
792       xlStart = initPos[0]*TMath::Cos(alpha)+initPos[1]*TMath::Sin(alpha);
793       t->AliExternalTrackParam::PropagateTo(xlStart,AliTracker::GetBz());   // to vtxSeed
794       rotAngle = alpha;
795       if(alpha<0.) rotAngle += 2.*TMath::Pi();
796       cosRot = TMath::Cos(rotAngle);
797       sinRot = TMath::Sin(rotAngle);
798       
799       // vector of track global coordinates
800       TMatrixD ri(3,1);
801       ri(0,0) = t->GetX()*cosRot-t->GetY()*sinRot;
802       ri(1,0) = t->GetX()*sinRot+t->GetY()*cosRot;
803       ri(2,0) = t->GetZ();
804
805       // matrix to go from global (x,y,z) to local (y,z);
806       TMatrixD qQi(2,3);
807       qQi(0,0) = -sinRot;
808       qQi(0,1) = cosRot;
809       qQi(0,2) = 0.;
810       qQi(1,0) = 0.;
811       qQi(1,1) = 0.;
812       qQi(1,2) = 1.;
813
814       // covariance matrix of local (y,z) - inverted
815       TMatrixD uUi(2,2);
816       t->GetExternalCovariance(cc);
817       uUi(0,0) = cc[0];
818       uUi(0,1) = cc[1];
819       uUi(1,0) = cc[1];
820       uUi(1,1) = cc[2];
821
822       // weights matrix: wWi = qQiT * uUiInv * qQi
823       if(uUi.Determinant() <= 0.) continue;
824       TMatrixD uUiInv(TMatrixD::kInverted,uUi);
825       TMatrixD uUiInvQi(uUiInv,TMatrixD::kMult,qQi);
826       TMatrixD wWi(qQi,TMatrixD::kTransposeMult,uUiInvQi);
827
828       // track chi2
829       TMatrixD deltar = rv; deltar -= ri;
830       TMatrixD wWideltar(wWi,TMatrixD::kMult,deltar);
831       chi2i = deltar(0,0)*wWideltar(0,0)+
832               deltar(1,0)*wWideltar(1,0)+
833               deltar(2,0)*wWideltar(2,0);
834
835
836       // add to total chi2
837       chi2 += chi2i;
838
839       TMatrixD wWiri(wWi,TMatrixD::kMult,ri); 
840
841       sumWiri += wWiri;
842       sumWi   += wWi;
843
844       nUsedTrks++;
845     } // end loop on tracks
846
847     if(nUsedTrks < fMinTracks) {
848       failed=1;
849       continue;
850     }
851
852     Double_t determinant = sumWi.Determinant();
853     //cerr<<" determinant: "<<determinant<<endl;
854     if(determinant < 100.)  { 
855       printf("det(V) = 0\n");       
856       failed=1;
857       continue;
858     }
859
860     // inverted of weights matrix
861     TMatrixD invsumWi(TMatrixD::kInverted,sumWi);
862     vV = invsumWi;
863      
864     // position of primary vertex
865     rv.Mult(vV,sumWiri);
866
867   } // end loop on the 2 steps
868
869   // delete t;
870
871   if(failed) { 
872     if(fDebug) printf("TooFewTracks\n");
873     fCurrentVertex = new AliESDVertex(0.,0.,-1);
874     return; 
875   }
876
877   Double_t position[3];
878   position[0] = rv(0,0);
879   position[1] = rv(1,0);
880   position[2] = rv(2,0);
881   Double_t covmatrix[6];
882   covmatrix[0] = vV(0,0);
883   covmatrix[1] = vV(0,1);
884   covmatrix[2] = vV(1,1);
885   covmatrix[3] = vV(0,2);
886   covmatrix[4] = vV(1,2);
887   covmatrix[5] = vV(2,2);
888   
889   // store data in the vertex object
890   fCurrentVertex = new AliESDVertex(position,covmatrix,chi2,nUsedTrks);
891
892   if(fDebug) {
893     printf(" VertexFitter(): finish\n");
894     printf(" rv = ( %f , %f , %f )\n\n",rv(0,0),rv(1,0),rv(2,0));
895     fCurrentVertex->PrintStatus();
896   }
897
898   return;
899 }
900 //----------------------------------------------------------------------------
901 AliVertex* AliVertexerTracks::VertexForSelectedTracks(TTree *trkTree) {
902 //
903 // Return vertex from tracks in trkTree
904 //
905
906   // get tracks and propagate them to initial vertex position
907   Int_t nTrksPrep = PrepareTracks(*trkTree,1);
908   if(nTrksPrep < fMinTracks) {
909     if(fDebug) printf("TooFewTracks\n");
910     Double_t vtx[3]={0,0,0};
911     fVert.SetXYZ(vtx);
912     fVert.SetDispersion(999);
913     fVert.SetNContributors(-5);
914     fTrkArray.Delete();
915     return &fVert;
916   }
917  
918   // Set initial vertex position from ESD
919   if(fAlgo==1)  StrLinVertexFinderMinDist(1);
920   if(fAlgo==2)  StrLinVertexFinderMinDist(0);
921   if(fAlgo==3)  HelixVertexFinder();
922   if(fAlgo==4)  VertexFinder(1);
923   if(fAlgo==5)  VertexFinder(0);
924
925   
926   // set indices of used tracks
927   UShort_t *indices = 0;
928   AliESDtrack *eta = 0;
929   if(fVert.GetNContributors()>0) {
930     indices = new UShort_t[fVert.GetNContributors()];
931     for(Int_t jj=0;jj<(Int_t)fTrkArray.GetEntriesFast();jj++) {
932       eta = (AliESDtrack*)fTrkArray.At(jj);
933       indices[jj] = (UShort_t)eta->GetID();
934     }
935     fVert.SetIndices(fVert.GetNContributors(),indices);
936   }
937   delete [] indices;
938   fTrkArray.Delete();
939   
940   return &fVert;
941 }
942 //----------------------------------------------------------------------------
943 AliVertex* AliVertexerTracks::VertexForSelectedTracks(TObjArray *trkArray) {
944 //
945 // Return vertex from array of tracks
946 //
947
948   // get tracks and propagate them to initial vertex position
949   Int_t nTrks = trkArray->GetEntriesFast();
950   if(nTrks < fMinTracks) {
951     if(fDebug) printf("TooFewTracks\n");
952     Double_t vtx[3]={0,0,0};
953     fVert.SetXYZ(vtx);
954     fVert.SetDispersion(999);
955     fVert.SetNContributors(-5);
956     fTrkArray.Delete();
957     return &fVert;
958   }
959   TTree *trkTree = new TTree("TreeT","tracks");
960   AliESDtrack *esdTrack = 0;
961   trkTree->Branch("tracks","AliESDtrack",&esdTrack);
962   for(Int_t i=0; i<nTrks; i++){
963     esdTrack = (AliESDtrack*)trkArray->At(i);
964     trkTree->Fill();
965   }
966     
967   AliVertex *vtx =  VertexForSelectedTracks(trkTree);
968   delete trkTree;
969   return vtx;
970 }
971 //--------------------------------------------------------------------------