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