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