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