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