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