]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliVertexerTracks.cxx
AliITSVertexerTracks becomes AliVertexerTrack (M.Masera, F.Prino)
[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 //---- AliRoot headers -----
31 #include "AliStrLine.h"
32 #include "AliVertexerTracks.h"
33 #include "AliESDtrack.h"
34
35 ClassImp(AliVertexerTracks)
36
37
38 //----------------------------------------------------------------------------
39 AliVertexerTracks::AliVertexerTracks():
40  TObject(),fVert() 
41 {
42 //
43 // Default constructor
44 //
45   SetVtxStart();
46   SetMinTracks();
47   fDCAcut=0;
48   fAlgo=1;
49 }
50 //-----------------------------------------------------------------------------
51 AliVertexerTracks::AliVertexerTracks(Double_t xStart, Double_t yStart):
52  TObject(),fVert() 
53 {
54 //
55 // Standard constructor
56 //
57   SetVtxStart(xStart,yStart);
58   SetMinTracks();
59   fDCAcut=0;
60   fAlgo=1;
61 }
62 //-----------------------------------------------------------------------------
63 AliVertexerTracks::~AliVertexerTracks() {
64   // Default Destructor
65   // The objects poited by the following pointers are not owned
66   // by this class and are not deleted
67 }
68 //----------------------------------------------------------------------------
69 Int_t AliVertexerTracks::PrepareTracks(TTree &trkTree) {
70 //
71 // Propagate tracks to initial vertex position and store them in a TObjArray
72 //
73   Double_t alpha,xlStart;
74   Int_t    nTrks    = 0;
75
76   Int_t    nEntries = (Int_t)trkTree.GetEntries();
77   if(!fTrkArray.IsEmpty()) fTrkArray.Clear();
78   fTrkArray.Expand(nEntries);
79
80   for(Int_t i=0; i<nEntries; i++) {
81     AliESDtrack *track = new AliESDtrack; 
82     trkTree.SetBranchAddress("tracks",&track);
83     trkTree.GetEvent(i);
84
85     // propagate track to vtxSeed
86     alpha  = track->GetAlpha();
87     xlStart = fNominalPos[0]*TMath::Cos(alpha)+fNominalPos[1]*TMath::Sin(alpha);
88     track->PropagateTo(xlStart,GetField());   // to vtxSeed
89     fTrkArray.AddLast(track);
90     nTrks++; 
91   }
92   return nTrks;
93
94 //----------------------------------------------------------------------------
95 AliVertex* AliVertexerTracks::VertexForSelectedTracks(TTree *trkTree) {
96 //
97 // Return vertex from tracks in trkTree
98 //
99
100   // get tracks and propagate them to initial vertex position
101   Int_t nTrks = PrepareTracks(*trkTree);
102   if(nTrks < fMinTracks) {
103     printf("TooFewTracks\n");
104     Double_t vtx[3]={0,0,0};
105     fVert.SetXYZ(vtx);
106     fVert.SetDispersion(999);
107     fVert.SetNContributors(-5);
108     return &fVert;
109   }
110  
111   // Set initial vertex position from ESD
112   if(fAlgo==1)  StrLinVertexFinderMinDist(1);
113   if(fAlgo==2)  StrLinVertexFinderMinDist(0);
114   if(fAlgo==3)  HelixVertexFinder();
115   if(fAlgo==4)  VertexFinder(1);
116   if(fAlgo==5)  VertexFinder(0);
117   return &fVert;
118 }
119 //----------------------------------------------------------------------------
120 AliVertex* AliVertexerTracks::VertexForSelectedTracks(TObjArray *trkArray) {
121 //
122 // Return vertex from array of tracks
123 //
124
125   // get tracks and propagate them to initial vertex position
126   Int_t nTrks = trkArray->GetEntriesFast();
127   if(nTrks < fMinTracks) {
128     printf("TooFewTracks\n");
129     Double_t vtx[3]={0,0,0};
130     fVert.SetXYZ(vtx);
131     fVert.SetDispersion(999);
132     fVert.SetNContributors(-5);
133     return &fVert;
134   }
135   TTree *trkTree = new TTree("TreeT","tracks");
136   AliESDtrack *esdTrack = 0;
137   trkTree->Branch("tracks","AliESDtrack",&esdTrack);
138   for(Int_t i=0; i<nTrks; i++){
139     esdTrack = (AliESDtrack*)trkArray->At(i);
140     trkTree->Fill();
141   }
142     
143   return VertexForSelectedTracks(trkTree);
144 }
145
146 //---------------------------------------------------------------------------
147 void AliVertexerTracks::VertexFinder(Int_t optUseWeights) {
148
149   // Get estimate of vertex position in (x,y) from tracks DCA
150  
151   Double_t initPos[3];
152   initPos[2] = 0.;
153   for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
154   Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
155   Double_t aver[3]={0.,0.,0.};
156   Double_t aversq[3]={0.,0.,0.};
157   Double_t sigmasq[3]={0.,0.,0.};
158   Double_t sigma=0;
159   Int_t ncombi = 0;
160   AliESDtrack *track1;
161   AliESDtrack *track2;
162   Double_t alpha,mindist;
163   Double_t field=GetField();
164
165   for(Int_t i=0; i<nacc; i++){
166     track1 = (AliESDtrack*)fTrkArray.At(i);
167     alpha=track1->GetAlpha();
168     mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
169     AliStrLine *line1 = new AliStrLine();
170     track1->ApproximateHelixWithLine(mindist,field,line1);
171    
172     for(Int_t j=i+1; j<nacc; j++){
173       track2 = (AliESDtrack*)fTrkArray.At(j);
174       alpha=track2->GetAlpha();
175       mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
176       AliStrLine *line2 = new AliStrLine();
177       track2->ApproximateHelixWithLine(mindist,field,line2);
178       Double_t distCA=line2->GetDCA(line1);
179       if(fDCAcut<=0 || (fDCAcut>0&&distCA<fDCAcut)){
180         Double_t pnt1[3],pnt2[3],crosspoint[3];
181
182         if(optUseWeights<=0){
183           Int_t retcode = line2->Cross(line1,crosspoint);
184           if(retcode>=0){
185             ncombi++;
186             for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
187             for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]);
188           }
189         }
190         if(optUseWeights>0){
191           Int_t retcode = line1->CrossPoints(line2,pnt1,pnt2);
192           if(retcode>=0){
193             Double_t alpha, cs, sn;
194             alpha=track1->GetAlpha();
195             cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);   
196             Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2();
197             alpha=track2->GetAlpha();
198             cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
199             Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2();
200             Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2();
201             Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
202             Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
203             Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
204             crosspoint[0]=wx1*pnt1[0] + wx2*pnt2[0]; 
205             crosspoint[1]=wy1*pnt1[1] + wy2*pnt2[1]; 
206             crosspoint[2]=wz1*pnt1[2] + wz2*pnt2[2];
207           
208             ncombi++;
209             for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
210             for(Int_t jj=0;jj<3;jj++)aversq[jj]+=(crosspoint[jj]*crosspoint[jj]);
211           }
212         }
213       }
214       delete line2;
215     }
216     delete line1;
217   }
218   if(ncombi>0){
219     for(Int_t jj=0;jj<3;jj++){
220       initPos[jj] = aver[jj]/ncombi;
221       aversq[jj]/=ncombi;
222       sigmasq[jj]=aversq[jj]-initPos[jj]*initPos[jj];
223       sigma+=sigmasq[jj];
224     }
225     sigma=TMath::Sqrt(TMath::Abs(sigma));
226   }
227   else {
228     Warning("VertexFinder","Finder did not succed");
229     sigma=999;
230   }
231   fVert.SetXYZ(initPos);
232   fVert.SetDispersion(sigma);
233   fVert.SetNContributors(ncombi);
234 }
235 //---------------------------------------------------------------------------
236 void AliVertexerTracks::HelixVertexFinder() {
237
238   // Get estimate of vertex position in (x,y) from tracks DCA
239
240
241   Double_t initPos[3];
242   initPos[2] = 0.;
243   for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
244   Double_t field=GetField();
245
246   Int_t nacc = (Int_t)fTrkArray.GetEntriesFast();
247
248   Double_t aver[3]={0.,0.,0.};
249   Double_t averquad[3]={0.,0.,0.};
250   Double_t sigmaquad[3]={0.,0.,0.};
251   Double_t sigma=0;
252   Int_t ncombi = 0;
253   AliESDtrack *track1;
254   AliESDtrack *track2;
255   Double_t distCA;
256   Double_t x, par[5];
257   Double_t alpha, cs, sn;
258   Double_t crosspoint[3];
259   for(Int_t i=0; i<nacc; i++){
260     track1 = (AliESDtrack*)fTrkArray.At(i);
261     
262
263     for(Int_t j=i+1; j<nacc; j++){
264       track2 = (AliESDtrack*)fTrkArray.At(j);
265
266       distCA=track2->PropagateToDCA(track1,field);
267
268       if(fDCAcut<=0 ||(fDCAcut>0&&distCA<fDCAcut)){
269         track1->GetExternalParameters(x,par);
270         alpha=track1->GetAlpha();
271         cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
272         Double_t x1=x*cs - par[0]*sn;
273         Double_t y1=x*sn + par[0]*cs;
274         Double_t z1=par[1];
275         Double_t sx1=sn*sn*track1->GetSigmaY2(), sy1=cs*cs*track1->GetSigmaY2(); 
276         track2->GetExternalParameters(x,par);
277         alpha=track2->GetAlpha();
278         cs=TMath::Cos(alpha); sn=TMath::Sin(alpha);
279         Double_t x2=x*cs - par[0]*sn;
280         Double_t y2=x*sn + par[0]*cs;
281         Double_t z2=par[1];
282         Double_t sx2=sn*sn*track2->GetSigmaY2(), sy2=cs*cs*track2->GetSigmaY2();
283         Double_t sz1=track1->GetSigmaZ2(), sz2=track2->GetSigmaZ2();
284         Double_t wx1=sx2/(sx1+sx2), wx2=1.- wx1;
285         Double_t wy1=sy2/(sy1+sy2), wy2=1.- wy1;
286         Double_t wz1=sz2/(sz1+sz2), wz2=1.- wz1;
287         crosspoint[0]=wx1*x1 + wx2*x2; 
288         crosspoint[1]=wy1*y1 + wy2*y2; 
289         crosspoint[2]=wz1*z1 + wz2*z2;
290
291         ncombi++;
292         for(Int_t jj=0;jj<3;jj++)aver[jj]+=crosspoint[jj];
293         for(Int_t jj=0;jj<3;jj++)averquad[jj]+=(crosspoint[jj]*crosspoint[jj]);
294       }
295     }
296       
297   }
298   if(ncombi>0){
299     for(Int_t jj=0;jj<3;jj++){
300       initPos[jj] = aver[jj]/ncombi;
301       averquad[jj]/=ncombi;
302       sigmaquad[jj]=averquad[jj]-initPos[jj]*initPos[jj];
303       sigma+=sigmaquad[jj];
304     }
305     sigma=TMath::Sqrt(TMath::Abs(sigma));
306   }
307   else {
308     Warning("HelixVertexFinder","Finder did not succed");
309     sigma=999;
310   }
311   fVert.SetXYZ(initPos);
312   fVert.SetDispersion(sigma);
313   fVert.SetNContributors(ncombi);
314 }
315 //---------------------------------------------------------------------------
316 void AliVertexerTracks::StrLinVertexFinderMinDist(Int_t optUseWeights){
317
318   // Calculate the point at minimum distance to prepared tracks 
319   
320   Double_t initPos[3];
321   initPos[2] = 0.;
322   Double_t sigma=0;
323   for(Int_t i=0;i<2;i++)initPos[i]=fNominalPos[i];
324   const Int_t knacc = (Int_t)fTrkArray.GetEntriesFast();
325   Double_t field=GetField();
326
327   AliESDtrack *track1;
328   Double_t (*vectP0)[3]=new Double_t [knacc][3];
329   Double_t (*vectP1)[3]=new Double_t [knacc][3];
330   
331   Double_t sum[3][3];
332   Double_t dsum[3]={0,0,0};
333   for(Int_t i=0;i<3;i++)
334     for(Int_t j=0;j<3;j++)sum[i][j]=0;
335   for(Int_t i=0; i<knacc; i++){
336     track1 = (AliESDtrack*)fTrkArray.At(i);
337     Double_t alpha=track1->GetAlpha();
338     Double_t mindist = TMath::Cos(alpha)*fNominalPos[0]+TMath::Sin(alpha)*fNominalPos[1];
339     AliStrLine *line1 = new AliStrLine();
340     track1->ApproximateHelixWithLine(mindist,field,line1);
341
342     Double_t p0[3],cd[3];
343     line1->GetP0(p0);
344     line1->GetCd(cd);
345     Double_t p1[3]={p0[0]+cd[0],p0[1]+cd[1],p0[2]+cd[2]};
346     vectP0[i][0]=p0[0];
347     vectP0[i][1]=p0[1];
348     vectP0[i][2]=p0[2];
349     vectP1[i][0]=p1[0];
350     vectP1[i][1]=p1[1];
351     vectP1[i][2]=p1[2];
352     
353     Double_t matr[3][3];
354     Double_t dknow[3];
355     if(optUseWeights==0)GetStrLinDerivMatrix(p0,p1,matr,dknow);
356     if(optUseWeights==1){
357       Double_t sigmasq[3];
358       sigmasq[0]=track1->GetSigmaY2();
359       sigmasq[1]=track1->GetSigmaY2();
360       sigmasq[2]=track1->GetSigmaZ2();
361       GetStrLinDerivMatrix(p0,p1,sigmasq,matr,dknow);
362     }
363
364     for(Int_t iii=0;iii<3;iii++){
365       dsum[iii]+=dknow[iii]; 
366       for(Int_t lj=0;lj<3;lj++) sum[iii][lj]+=matr[iii][lj];
367     }
368     delete line1;
369   }
370   
371   Double_t vett[3][3];
372   Double_t det=GetDeterminant3X3(sum);
373   
374    if(det!=0){
375      for(Int_t zz=0;zz<3;zz++){
376        for(Int_t ww=0;ww<3;ww++){
377          for(Int_t kk=0;kk<3;kk++) vett[ww][kk]=sum[ww][kk];
378        }
379        for(Int_t kk=0;kk<3;kk++) vett[kk][zz]=dsum[kk];
380        initPos[zz]=GetDeterminant3X3(vett)/det;
381      }
382
383
384      for(Int_t i=0; i<knacc; i++){
385        Double_t p0[3]={0,0,0},p1[3]={0,0,0};
386        for(Int_t ii=0;ii<3;ii++){
387          p0[ii]=vectP0[i][ii];
388          p1[ii]=vectP1[i][ii];
389        }
390        sigma+=GetStrLinMinDist(p0,p1,initPos);
391      }
392
393      sigma=TMath::Sqrt(sigma);
394    }else{
395     Warning("StrLinVertexFinderMinDist","Finder did not succed");
396     sigma=999;
397   }
398   delete vectP0;
399   delete vectP1;
400   fVert.SetXYZ(initPos);
401   fVert.SetDispersion(sigma);
402   fVert.SetNContributors(knacc);
403 }
404 //_______________________________________________________________________
405 Double_t AliVertexerTracks::GetDeterminant3X3(Double_t matr[][3]){
406   //
407   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];
408  return det;
409 }
410 //____________________________________________________________________________
411 void AliVertexerTracks::GetStrLinDerivMatrix(Double_t *p0,Double_t *p1,Double_t m[][3],Double_t *d){
412
413   //
414   Double_t x12=p0[0]-p1[0];
415   Double_t y12=p0[1]-p1[1];
416   Double_t z12=p0[2]-p1[2];
417   Double_t kk=x12*x12+y12*y12+z12*z12;
418   m[0][0]=2-2/kk*x12*x12;
419   m[0][1]=-2/kk*x12*y12;
420   m[0][2]=-2/kk*x12*z12;
421   m[1][0]=-2/kk*x12*y12;
422   m[1][1]=2-2/kk*y12*y12;
423   m[1][2]=-2/kk*y12*z12;
424   m[2][0]=-2/kk*x12*z12;
425   m[2][1]=-2*y12*z12;
426   m[2][2]=2-2/kk*z12*z12;
427   d[0]=2*p0[0]-2/kk*p0[0]*x12*x12-2/kk*p0[2]*x12*z12-2/kk*p0[1]*x12*y12;
428   d[1]=2*p0[1]-2/kk*p0[1]*y12*y12-2/kk*p0[0]*x12*y12-2/kk*p0[2]*z12*y12;
429   d[2]=2*p0[2]-2/kk*p0[2]*z12*z12-2/kk*p0[0]*x12*z12-2/kk*p0[1]*z12*y12;
430
431 }
432 //____________________________________________________________________________
433 void AliVertexerTracks::GetStrLinDerivMatrix(Double_t *p0,Double_t *p1,Double_t *sigmasq,Double_t m[][3],Double_t *d){
434   //
435   Double_t x12=p1[0]-p0[0];
436   Double_t y12=p1[1]-p0[1];
437   Double_t z12=p1[2]-p0[2];
438
439   Double_t den= x12*x12*sigmasq[1]*sigmasq[2]+y12*y12*sigmasq[0]*sigmasq[2]+z12*z12*sigmasq[0]*sigmasq[1];
440
441   Double_t kk= 2*(x12*x12/sigmasq[0]+y12*y12/sigmasq[1]+z12*z12/sigmasq[2]);
442
443   Double_t cc[3];
444   cc[0]=-x12/sigmasq[0];
445   cc[1]=-y12/sigmasq[1];
446   cc[2]=-z12/sigmasq[2];
447
448   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;
449
450   Double_t ss= -p0[0]*cc[0]-p0[1]*cc[1]-p0[2]*cc[2];
451
452   Double_t aa[3];
453   aa[0]=x12*sigmasq[1]*sigmasq[2]/den;
454   aa[1]=y12*sigmasq[0]*sigmasq[2]/den;
455   aa[2]=z12*sigmasq[0]*sigmasq[1]/den;
456
457   m[0][0]=aa[0]*(aa[0]*kk+2*cc[0])+2*cc[0]*aa[0]+2/sigmasq[0];
458   m[0][1]=aa[1]*(aa[0]*kk+2*cc[0])+2*cc[1]*aa[0];
459   m[0][2]=aa[2]*(aa[0]*kk+2*cc[0])+2*cc[2]*aa[0];
460
461   m[1][0]=aa[0]*(aa[1]*kk+2*cc[1])+2*cc[0]*aa[1];
462   m[1][1]=aa[1]*(aa[1]*kk+2*cc[1])+2*cc[1]*aa[1]+2/sigmasq[1];
463   m[1][2]=aa[2]*(aa[1]*kk+2*cc[1])+2*cc[2]*aa[1];
464
465   m[2][0]=aa[0]*(aa[2]*kk+2*cc[2])+2*cc[0]*aa[2];
466   m[2][1]=aa[1]*(aa[2]*kk+2*cc[2])+2*cc[1]*aa[2];
467   m[2][2]=aa[2]*(aa[2]*kk+2*cc[2])+2*cc[2]*aa[2]+2/sigmasq[2];
468
469   d[0]=-ww*(aa[0]*kk+2*cc[0])-2*ss*aa[0]+2*p0[0]/sigmasq[0];
470   d[1]=-ww*(aa[1]*kk+2*cc[1])-2*ss*aa[1]+2*p0[1]/sigmasq[1];
471   d[2]=-ww*(aa[2]*kk+2*cc[2])-2*ss*aa[2]+2*p0[2]/sigmasq[2];
472
473 }
474 //_____________________________________________________________________________
475 Double_t AliVertexerTracks::GetStrLinMinDist(Double_t *p0,Double_t *p1,Double_t *x0){
476   //
477   Double_t x12=p0[0]-p1[0];
478   Double_t y12=p0[1]-p1[1];
479   Double_t z12=p0[2]-p1[2];
480   Double_t x10=p0[0]-x0[0];
481   Double_t y10=p0[1]-x0[1];
482   Double_t z10=p0[2]-x0[2];
483   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);
484 }