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