]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliExternalTrackParam.cxx
Moving the AliStrLine and AliTracker classes from libESD.so back to libSTEER.so....
[u/mrichter/AliRoot.git] / STEER / AliExternalTrackParam.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 // Implementation of the external track parameterisation class.              //
21 //                                                                           //
22 // This parameterisation is used to exchange tracks between the detectors.   //
23 // A set of functions returning the position and the momentum of tracks      //
24 // in the global coordinate system as well as the track impact parameters    //
25 // are implemented.
26 // Origin: I.Belikov, CERN, Jouri.Belikov@cern.ch                            //
27 ///////////////////////////////////////////////////////////////////////////////
28 #include "AliExternalTrackParam.h"
29 #include "AliKalmanTrack.h"
30 #include "AliESDVertex.h"
31
32
33 ClassImp(AliExternalTrackParam)
34
35 //_____________________________________________________________________________
36 AliExternalTrackParam::AliExternalTrackParam() :
37   fX(0),
38   fAlpha(0)
39 {
40   //
41   // default constructor
42   //
43   for (Int_t i = 0; i < 5; i++) fP[i] = 0;
44   for (Int_t i = 0; i < 15; i++) fC[i] = 0;
45 }
46
47 //_____________________________________________________________________________
48 AliExternalTrackParam::AliExternalTrackParam(Double_t x, Double_t alpha, 
49                                              const Double_t param[5], 
50                                              const Double_t covar[15]) :
51   fX(x),
52   fAlpha(alpha)
53 {
54   //
55   // create external track parameters from given arguments
56   //
57   for (Int_t i = 0; i < 5; i++)  fP[i] = param[i];
58   for (Int_t i = 0; i < 15; i++) fC[i] = covar[i];
59 }
60
61 //_____________________________________________________________________________
62 AliExternalTrackParam::AliExternalTrackParam(const AliKalmanTrack& track) :
63   fAlpha(track.GetAlpha())
64 {
65   //
66   //
67   track.GetExternalParameters(fX,fP);
68   track.GetExternalCovariance(fC);
69 }
70
71 //_____________________________________________________________________________
72 void AliExternalTrackParam::Set(const AliKalmanTrack& track) {
73   //
74   //
75   fAlpha=track.GetAlpha();
76   track.GetExternalParameters(fX,fP);
77   track.GetExternalCovariance(fC);
78 }
79
80 //_____________________________________________________________________________
81 void AliExternalTrackParam::Reset() {
82   fX=fAlpha=0.;
83   for (Int_t i = 0; i < 5; i++) fP[i] = 0;
84   for (Int_t i = 0; i < 15; i++) fC[i] = 0;
85 }
86
87 Double_t AliExternalTrackParam::GetP() const {
88   //---------------------------------------------------------------------
89   // This function returns the track momentum
90   // Results for (nearly) straight tracks are meaningless !
91   //---------------------------------------------------------------------
92   if (TMath::Abs(fP[4])<=kAlmost0) return kVeryBig;
93   return TMath::Sqrt(1.+ fP[3]*fP[3])/TMath::Abs(fP[4]);
94 }
95
96 Double_t AliExternalTrackParam::Get1P() const {
97   //---------------------------------------------------------------------
98   // This function returns the 1/(track momentum)
99   //---------------------------------------------------------------------
100   return TMath::Abs(fP[4])/TMath::Sqrt(1.+ fP[3]*fP[3]);
101 }
102
103 //_______________________________________________________________________
104 Double_t AliExternalTrackParam::GetD(Double_t x,Double_t y,Double_t b) const {
105   //------------------------------------------------------------------
106   // This function calculates the transverse impact parameter
107   // with respect to a point with global coordinates (x,y)
108   // in the magnetic field "b" (kG)
109   //------------------------------------------------------------------
110   if (TMath::Abs(b) < kAlmost0Field) return GetLinearD(x,y);
111   Double_t rp4=kB2C*b*fP[4];
112
113   Double_t xt=fX, yt=fP[0];
114
115   Double_t sn=TMath::Sin(fAlpha), cs=TMath::Cos(fAlpha);
116   Double_t a = x*cs + y*sn;
117   y = -x*sn + y*cs; x=a;
118   xt-=x; yt-=y;
119
120   sn=rp4*xt - fP[2]; cs=rp4*yt + TMath::Sqrt(1.- fP[2]*fP[2]);
121   a=2*(xt*fP[2] - yt*TMath::Sqrt(1.- fP[2]*fP[2]))-rp4*(xt*xt + yt*yt);
122   if (rp4<0) a=-a;
123   return a/(1 + TMath::Sqrt(sn*sn + cs*cs));
124 }
125
126 //_______________________________________________________________________
127 Double_t AliExternalTrackParam::GetLinearD(Double_t xv,Double_t yv) const {
128   //------------------------------------------------------------------
129   // This function calculates the transverse impact parameter
130   // with respect to a point with global coordinates (xv,yv)
131   // neglecting the track curvature.
132   //------------------------------------------------------------------
133   Double_t sn=TMath::Sin(fAlpha), cs=TMath::Cos(fAlpha);
134   Double_t x= xv*cs + yv*sn;
135   Double_t y=-xv*sn + yv*cs;
136
137   Double_t d = (fX-x)*fP[2] - (fP[0]-y)*TMath::Sqrt(1.- fP[2]*fP[2]);
138
139   return d;
140 }
141
142 Bool_t AliExternalTrackParam::
143 CorrectForMaterial(Double_t d,  Double_t x0, Double_t mass) {
144   //------------------------------------------------------------------
145   // This function corrects the track parameters for the crossed material
146   // "d"    - the thickness (fraction of the radiation length)
147   // "x0"   - the radiation length (g/cm^2) 
148   // "mass" - the mass of this particle (GeV/c^2)
149   //------------------------------------------------------------------
150   Double_t &fP2=fP[2];
151   Double_t &fP3=fP[3];
152   Double_t &fP4=fP[4];
153
154   Double_t &fC22=fC[5];
155   Double_t &fC33=fC[9];
156   Double_t &fC43=fC[13];
157   Double_t &fC44=fC[14];
158
159   Double_t p2=(1.+ fP3*fP3)/(fP4*fP4);
160   Double_t beta2=p2/(p2 + mass*mass);
161   d*=TMath::Sqrt((1.+ fP3*fP3)/(1.- fP2*fP2));
162
163   //Multiple scattering******************
164   if (d!=0) {
165      Double_t theta2=14.1*14.1/(beta2*p2*1e6)*TMath::Abs(d);
166      //Double_t theta2=1.0259e-6*14*14/28/(beta2*p2)*TMath::Abs(d)*9.36*2.33;
167      fC22 += theta2*(1.- fP2*fP2)*(1. + fP3*fP3);
168      fC33 += theta2*(1. + fP3*fP3)*(1. + fP3*fP3);
169      fC43 += theta2*fP3*fP4*(1. + fP3*fP3);
170      fC44 += theta2*fP3*fP4*fP3*fP4;
171   }
172
173   //Energy losses************************
174   if (x0!=0.) {
175      d*=x0;
176      Double_t dE=0.153e-3/beta2*(log(5940*beta2/(1-beta2)) - beta2)*d;
177      if (beta2/(1-beta2)>3.5*3.5)
178        dE=0.153e-3/beta2*(log(3.5*5940)+0.5*log(beta2/(1-beta2)) - beta2)*d;
179
180      fP4*=(1.- TMath::Sqrt(p2 + mass*mass)/p2*dE);
181   }
182
183   return kTRUE;
184 }
185
186 Bool_t AliExternalTrackParam::Rotate(Double_t alpha) {
187   //------------------------------------------------------------------
188   // Transform this track to the local coord. system rotated
189   // by angle "alpha" (rad) with respect to the global coord. system. 
190   //------------------------------------------------------------------
191   if      (alpha < -TMath::Pi()) alpha += 2*TMath::Pi();
192   else if (alpha >= TMath::Pi()) alpha -= 2*TMath::Pi();
193
194   Double_t &fP0=fP[0];
195   Double_t &fP2=fP[2];
196   Double_t &fC00=fC[0];
197   Double_t &fC10=fC[1];
198   Double_t &fC20=fC[3];
199   Double_t &fC21=fC[4];
200   Double_t &fC22=fC[5];
201   Double_t &fC30=fC[6];
202   Double_t &fC32=fC[8];
203   Double_t &fC40=fC[10];
204   Double_t &fC42=fC[12];
205
206   Double_t x=fX;
207   Double_t ca=TMath::Cos(alpha-fAlpha), sa=TMath::Sin(alpha-fAlpha);
208   Double_t sf=fP2, cf=TMath::Sqrt(1.- fP2*fP2);
209
210   fAlpha = alpha;
211   fX =  x*ca + fP0*sa;
212   fP0= -x*sa + fP0*ca;
213   fP2=  sf*ca - cf*sa;
214
215   if (TMath::Abs(cf)<kAlmost0) {
216     AliError(Form("Too small cosine value %f",cf)); 
217     cf = kAlmost0;
218   } 
219
220   Double_t rr=(ca+sf/cf*sa);  
221
222   fC00 *= (ca*ca);
223   fC10 *= ca;
224   fC20 *= ca*rr;
225   fC21 *= rr;
226   fC22 *= rr*rr;
227   fC30 *= ca;
228   fC32 *= rr;
229   fC40 *= ca;
230   fC42 *= rr;
231
232   return kTRUE;
233 }
234
235 Bool_t AliExternalTrackParam::PropagateTo(Double_t xk, Double_t b) {
236   //----------------------------------------------------------------
237   // Propagate this track to the plane X=xk (cm) in the field "b" (kG)
238   //----------------------------------------------------------------
239   Double_t dx=xk-fX;
240   if (TMath::Abs(dx)<=0)  return kTRUE;
241
242   Double_t crv=kB2C*b*fP[4];
243   if (TMath::Abs(b) < kAlmost0Field) crv=0.;
244
245   Double_t f1=fP[2], f2=f1 + crv*dx;
246   if (TMath::Abs(f1) >= kAlmost1) return kFALSE;
247   if (TMath::Abs(f2) >= kAlmost1) return kFALSE;
248
249   Double_t &fP0=fP[0], &fP1=fP[1], &fP2=fP[2], &fP3=fP[3], &fP4=fP[4];
250   Double_t 
251   &fC00=fC[0],
252   &fC10=fC[1],   &fC11=fC[2],  
253   &fC20=fC[3],   &fC21=fC[4],   &fC22=fC[5],
254   &fC30=fC[6],   &fC31=fC[7],   &fC32=fC[8],   &fC33=fC[9],  
255   &fC40=fC[10],  &fC41=fC[11],  &fC42=fC[12],  &fC43=fC[13], &fC44=fC[14];
256
257   Double_t r1=TMath::Sqrt(1.- f1*f1), r2=TMath::Sqrt(1.- f2*f2);
258
259   fX=xk;
260   fP0 += dx*(f1+f2)/(r1+r2);
261   fP1 += dx*(r2 + f2*(f1+f2)/(r1+r2))*fP3;  // Many thanks to P.Hristov !
262   fP2 += dx*crv;
263
264   //f = F - 1
265    
266   Double_t f02=    dx/(r1*r1*r1);            Double_t cc=crv/fP4;
267   Double_t f04=0.5*dx*dx/(r1*r1*r1);         f04*=cc;
268   Double_t f12=    dx*fP3*f1/(r1*r1*r1);
269   Double_t f14=0.5*dx*dx*fP3*f1/(r1*r1*r1);  f14*=cc;
270   Double_t f13=    dx/r1;
271   Double_t f24=    dx;                       f24*=cc;
272   
273   //b = C*ft
274   Double_t b00=f02*fC20 + f04*fC40, b01=f12*fC20 + f14*fC40 + f13*fC30;
275   Double_t b02=f24*fC40;
276   Double_t b10=f02*fC21 + f04*fC41, b11=f12*fC21 + f14*fC41 + f13*fC31;
277   Double_t b12=f24*fC41;
278   Double_t b20=f02*fC22 + f04*fC42, b21=f12*fC22 + f14*fC42 + f13*fC32;
279   Double_t b22=f24*fC42;
280   Double_t b40=f02*fC42 + f04*fC44, b41=f12*fC42 + f14*fC44 + f13*fC43;
281   Double_t b42=f24*fC44;
282   Double_t b30=f02*fC32 + f04*fC43, b31=f12*fC32 + f14*fC43 + f13*fC33;
283   Double_t b32=f24*fC43;
284   
285   //a = f*b = f*C*ft
286   Double_t a00=f02*b20+f04*b40,a01=f02*b21+f04*b41,a02=f02*b22+f04*b42;
287   Double_t a11=f12*b21+f14*b41+f13*b31,a12=f12*b22+f14*b42+f13*b32;
288   Double_t a22=f24*b42;
289
290   //F*C*Ft = C + (b + bt + a)
291   fC00 += b00 + b00 + a00;
292   fC10 += b10 + b01 + a01; 
293   fC20 += b20 + b02 + a02;
294   fC30 += b30;
295   fC40 += b40;
296   fC11 += b11 + b11 + a11;
297   fC21 += b21 + b12 + a12;
298   fC31 += b31; 
299   fC41 += b41;
300   fC22 += b22 + b22 + a22;
301   fC32 += b32;
302   fC42 += b42;
303
304   return kTRUE;
305 }
306
307 Double_t 
308 AliExternalTrackParam::GetPredictedChi2(Double_t p[2],Double_t cov[3]) const {
309   //----------------------------------------------------------------
310   // Estimate the chi2 of the space point "p" with the cov. matrix "cov"
311   //----------------------------------------------------------------
312   Double_t sdd = fC[0] + cov[0]; 
313   Double_t sdz = fC[1] + cov[1];
314   Double_t szz = fC[2] + cov[2];
315   Double_t det = sdd*szz - sdz*sdz;
316
317   if (TMath::Abs(det) < kAlmost0) return kVeryBig;
318
319   Double_t d = fP[0] - p[0];
320   Double_t z = fP[1] - p[1];
321
322   return (d*szz*d - 2*d*sdz*z + z*sdd*z)/det;
323 }
324
325 Bool_t AliExternalTrackParam::Update(Double_t p[2], Double_t cov[3]) {
326   //------------------------------------------------------------------
327   // Update the track parameters with the space point "p" having
328   // the covariance matrix "cov"
329   //------------------------------------------------------------------
330   Double_t &fP0=fP[0], &fP1=fP[1], &fP2=fP[2], &fP3=fP[3], &fP4=fP[4];
331   Double_t 
332   &fC00=fC[0],
333   &fC10=fC[1],   &fC11=fC[2],  
334   &fC20=fC[3],   &fC21=fC[4],   &fC22=fC[5],
335   &fC30=fC[6],   &fC31=fC[7],   &fC32=fC[8],   &fC33=fC[9],  
336   &fC40=fC[10],  &fC41=fC[11],  &fC42=fC[12],  &fC43=fC[13], &fC44=fC[14];
337
338   Double_t r00=cov[0], r01=cov[1], r11=cov[2];
339   r00+=fC00; r01+=fC10; r11+=fC11;
340   Double_t det=r00*r11 - r01*r01;
341
342   if (TMath::Abs(det) < kAlmost0) return kFALSE;
343
344
345   Double_t tmp=r00; r00=r11/det; r11=tmp/det; r01=-r01/det;
346  
347   Double_t k00=fC00*r00+fC10*r01, k01=fC00*r01+fC10*r11;
348   Double_t k10=fC10*r00+fC11*r01, k11=fC10*r01+fC11*r11;
349   Double_t k20=fC20*r00+fC21*r01, k21=fC20*r01+fC21*r11;
350   Double_t k30=fC30*r00+fC31*r01, k31=fC30*r01+fC31*r11;
351   Double_t k40=fC40*r00+fC41*r01, k41=fC40*r01+fC41*r11;
352
353   Double_t dy=p[0] - fP0, dz=p[1] - fP1;
354   Double_t sf=fP2 + k20*dy + k21*dz;
355   if (TMath::Abs(sf) > kAlmost1) return kFALSE;  
356   
357   fP0 += k00*dy + k01*dz;
358   fP1 += k10*dy + k11*dz;
359   fP2  = sf;
360   fP3 += k30*dy + k31*dz;
361   fP4 += k40*dy + k41*dz;
362   
363   Double_t c01=fC10, c02=fC20, c03=fC30, c04=fC40;
364   Double_t c12=fC21, c13=fC31, c14=fC41;
365
366   fC00-=k00*fC00+k01*fC10; fC10-=k00*c01+k01*fC11;
367   fC20-=k00*c02+k01*c12;   fC30-=k00*c03+k01*c13;
368   fC40-=k00*c04+k01*c14; 
369
370   fC11-=k10*c01+k11*fC11;
371   fC21-=k10*c02+k11*c12;   fC31-=k10*c03+k11*c13;
372   fC41-=k10*c04+k11*c14; 
373
374   fC22-=k20*c02+k21*c12;   fC32-=k20*c03+k21*c13;
375   fC42-=k20*c04+k21*c14; 
376
377   fC33-=k30*c03+k31*c13;
378   fC43-=k30*c04+k31*c14; 
379
380   fC44-=k40*c04+k41*c14; 
381
382   return kTRUE;
383 }
384
385 void 
386 AliExternalTrackParam::GetHelixParameters(Double_t hlx[6], Double_t b) const {
387   //--------------------------------------------------------------------
388   // External track parameters -> helix parameters 
389   // "b" - magnetic field (kG)
390   //--------------------------------------------------------------------
391   Double_t cs=TMath::Cos(fAlpha), sn=TMath::Sin(fAlpha);
392   
393   hlx[0]=fP[0]; hlx[1]=fP[1]; hlx[2]=fP[2]; hlx[3]=fP[3]; hlx[4]=fP[4];
394
395   hlx[5]=fX*cs - hlx[0]*sn;               // x0
396   hlx[0]=fX*sn + hlx[0]*cs;               // y0
397 //hlx[1]=                                 // z0
398   hlx[2]=TMath::ASin(hlx[2]) + fAlpha;    // phi0
399 //hlx[3]=                                 // tgl
400   hlx[4]=hlx[4]*kB2C*b;                   // C
401 }
402
403
404 static void Evaluate(const Double_t *h, Double_t t,
405                      Double_t r[3],  //radius vector
406                      Double_t g[3],  //first defivatives
407                      Double_t gg[3]) //second derivatives
408 {
409   //--------------------------------------------------------------------
410   // Calculate position of a point on a track and some derivatives
411   //--------------------------------------------------------------------
412   Double_t phase=h[4]*t+h[2];
413   Double_t sn=TMath::Sin(phase), cs=TMath::Cos(phase);
414
415   r[0] = h[5] + (sn - h[6])/h[4];
416   r[1] = h[0] - (cs - h[7])/h[4];  
417   r[2] = h[1] + h[3]*t;
418
419   g[0] = cs; g[1]=sn; g[2]=h[3];
420   
421   gg[0]=-h[4]*sn; gg[1]=h[4]*cs; gg[2]=0.;
422 }
423
424 Double_t AliExternalTrackParam::GetDCA(const AliExternalTrackParam *p, 
425 Double_t b, Double_t &xthis, Double_t &xp) const {
426   //------------------------------------------------------------
427   // Returns the (weighed !) distance of closest approach between 
428   // this track and the track "p".
429   // Other returned values:
430   //   xthis, xt - coordinates of tracks' reference planes at the DCA 
431   //-----------------------------------------------------------
432   Double_t dy2=GetSigmaY2() + p->GetSigmaY2();
433   Double_t dz2=GetSigmaZ2() + p->GetSigmaZ2();
434   Double_t dx2=dy2; 
435
436   //dx2=dy2=dz2=1.;
437
438   Double_t p1[8]; GetHelixParameters(p1,b);
439   p1[6]=TMath::Sin(p1[2]); p1[7]=TMath::Cos(p1[2]);
440   Double_t p2[8]; p->GetHelixParameters(p2,b);
441   p2[6]=TMath::Sin(p2[2]); p2[7]=TMath::Cos(p2[2]);
442
443
444   Double_t r1[3],g1[3],gg1[3]; Double_t t1=0.;
445   Evaluate(p1,t1,r1,g1,gg1);
446   Double_t r2[3],g2[3],gg2[3]; Double_t t2=0.;
447   Evaluate(p2,t2,r2,g2,gg2);
448
449   Double_t dx=r2[0]-r1[0], dy=r2[1]-r1[1], dz=r2[2]-r1[2];
450   Double_t dm=dx*dx/dx2 + dy*dy/dy2 + dz*dz/dz2;
451
452   Int_t max=27;
453   while (max--) {
454      Double_t gt1=-(dx*g1[0]/dx2 + dy*g1[1]/dy2 + dz*g1[2]/dz2);
455      Double_t gt2=+(dx*g2[0]/dx2 + dy*g2[1]/dy2 + dz*g2[2]/dz2);
456      Double_t h11=(g1[0]*g1[0] - dx*gg1[0])/dx2 + 
457                   (g1[1]*g1[1] - dy*gg1[1])/dy2 +
458                   (g1[2]*g1[2] - dz*gg1[2])/dz2;
459      Double_t h22=(g2[0]*g2[0] + dx*gg2[0])/dx2 + 
460                   (g2[1]*g2[1] + dy*gg2[1])/dy2 +
461                   (g2[2]*g2[2] + dz*gg2[2])/dz2;
462      Double_t h12=-(g1[0]*g2[0]/dx2 + g1[1]*g2[1]/dy2 + g1[2]*g2[2]/dz2);
463
464      Double_t det=h11*h22-h12*h12;
465
466      Double_t dt1,dt2;
467      if (TMath::Abs(det)<1.e-33) {
468         //(quasi)singular Hessian
469         dt1=-gt1; dt2=-gt2;
470      } else {
471         dt1=-(gt1*h22 - gt2*h12)/det; 
472         dt2=-(h11*gt2 - h12*gt1)/det;
473      }
474
475      if ((dt1*gt1+dt2*gt2)>0) {dt1=-dt1; dt2=-dt2;}
476
477      //check delta(phase1) ?
478      //check delta(phase2) ?
479
480      if (TMath::Abs(dt1)/(TMath::Abs(t1)+1.e-3) < 1.e-4)
481      if (TMath::Abs(dt2)/(TMath::Abs(t2)+1.e-3) < 1.e-4) {
482         if ((gt1*gt1+gt2*gt2) > 1.e-4/dy2/dy2) 
483           AliWarning(" stopped at not a stationary point !");
484         Double_t lmb=h11+h22; lmb=lmb-TMath::Sqrt(lmb*lmb-4*det);
485         if (lmb < 0.) 
486           AliWarning(" stopped at not a minimum !");
487         break;
488      }
489
490      Double_t dd=dm;
491      for (Int_t div=1 ; ; div*=2) {
492         Evaluate(p1,t1+dt1,r1,g1,gg1);
493         Evaluate(p2,t2+dt2,r2,g2,gg2);
494         dx=r2[0]-r1[0]; dy=r2[1]-r1[1]; dz=r2[2]-r1[2];
495         dd=dx*dx/dx2 + dy*dy/dy2 + dz*dz/dz2;
496         if (dd<dm) break;
497         dt1*=0.5; dt2*=0.5;
498         if (div>512) {
499            AliWarning(" overshoot !"); break;
500         }   
501      }
502      dm=dd;
503
504      t1+=dt1;
505      t2+=dt2;
506
507   }
508
509   if (max<=0) AliWarning(" too many iterations !");
510
511   Double_t cs=TMath::Cos(GetAlpha());
512   Double_t sn=TMath::Sin(GetAlpha());
513   xthis=r1[0]*cs + r1[1]*sn;
514
515   cs=TMath::Cos(p->GetAlpha());
516   sn=TMath::Sin(p->GetAlpha());
517   xp=r2[0]*cs + r2[1]*sn;
518
519   return TMath::Sqrt(dm*TMath::Sqrt(dy2*dz2));
520 }
521  
522 Double_t AliExternalTrackParam::
523 PropagateToDCA(AliExternalTrackParam *p, Double_t b) {
524   //--------------------------------------------------------------
525   // Propagates this track and the argument track to the position of the
526   // distance of closest approach.
527   // Returns the (weighed !) distance of closest approach.
528   //--------------------------------------------------------------
529   Double_t xthis,xp;
530   Double_t dca=GetDCA(p,b,xthis,xp);
531
532   if (!PropagateTo(xthis,b)) {
533     //AliWarning(" propagation failed !");
534     return 1e+33;
535   }
536
537   if (!p->PropagateTo(xp,b)) {
538     //AliWarning(" propagation failed !";
539     return 1e+33;
540   }
541
542   return dca;
543 }
544
545
546
547
548 Bool_t AliExternalTrackParam::PropagateToDCA(const AliESDVertex *vtx, Double_t b, Double_t maxd){
549   //
550   // Try to relate this track to the vertex "vtx", 
551   // if the (rough) transverse impact parameter is not bigger then "maxd". 
552   //            Magnetic field is "b" (kG).
553   //
554   // a) The track gets extapolated to the DCA to the vertex.
555   // b) The impact parameters and their covariance matrix are calculated.
556   //
557   //    In the case of success, the returned value is kTRUE
558   //    (otherwise, it's kFALSE)
559   //  
560   Double_t alpha=GetAlpha();
561   Double_t sn=TMath::Sin(alpha), cs=TMath::Cos(alpha);
562   Double_t x=GetX(), y=GetParameter()[0], snp=GetParameter()[2];
563   Double_t xv= vtx->GetXv()*cs + vtx->GetYv()*sn;
564   Double_t yv=-vtx->GetXv()*sn + vtx->GetYv()*cs;
565   x-=xv; y-=yv;
566
567   //Estimate the impact parameter neglecting the track curvature
568   Double_t d=TMath::Abs(x*snp - y*TMath::Sqrt(1.- snp*snp));
569   if (d > maxd) return kFALSE; 
570
571   //Propagate to the DCA
572   Double_t crv=0.299792458e-3*b*GetParameter()[4];
573   Double_t tgfv=-(crv*x - snp)/(crv*y + TMath::Sqrt(1.-snp*snp));
574   sn=tgfv/TMath::Sqrt(1.+ tgfv*tgfv); cs=TMath::Sqrt(1.- sn*sn);
575
576   x = xv*cs + yv*sn;
577   yv=-xv*sn + yv*cs; xv=x;
578
579   if (!Propagate(alpha+TMath::ASin(sn),xv,b)) return kFALSE;
580   return kTRUE;
581 }
582
583
584
585
586 Bool_t Local2GlobalMomentum(Double_t p[3],Double_t alpha) {
587   //----------------------------------------------------------------
588   // This function performs local->global transformation of the
589   // track momentum.
590   // When called, the arguments are:
591   //    p[0] = 1/pt of the track;
592   //    p[1] = sine of local azim. angle of the track momentum;
593   //    p[2] = tangent of the track momentum dip angle;
594   //   alpha - rotation angle. 
595   // The result is returned as:
596   //    p[0] = px
597   //    p[1] = py
598   //    p[2] = pz
599   // Results for (nearly) straight tracks are meaningless !
600   //----------------------------------------------------------------
601   if (TMath::Abs(p[0])<=0)        return kFALSE;
602   if (TMath::Abs(p[1])> kAlmost1) return kFALSE;
603
604   Double_t pt=1./TMath::Abs(p[0]);
605   Double_t cs=TMath::Cos(alpha), sn=TMath::Sin(alpha);
606   Double_t r=TMath::Sqrt(1 - p[1]*p[1]);
607   p[0]=pt*(r*cs - p[1]*sn); p[1]=pt*(p[1]*cs + r*sn); p[2]=pt*p[2];
608
609   return kTRUE;
610 }
611
612 Bool_t Local2GlobalPosition(Double_t r[3],Double_t alpha) {
613   //----------------------------------------------------------------
614   // This function performs local->global transformation of the
615   // track position.
616   // When called, the arguments are:
617   //    r[0] = local x
618   //    r[1] = local y
619   //    r[2] = local z
620   //   alpha - rotation angle. 
621   // The result is returned as:
622   //    r[0] = global x
623   //    r[1] = global y
624   //    r[2] = global z
625   //----------------------------------------------------------------
626   Double_t cs=TMath::Cos(alpha), sn=TMath::Sin(alpha), x=r[0];
627   r[0]=x*cs - r[1]*sn; r[1]=x*sn + r[1]*cs;
628
629   return kTRUE;
630 }
631
632 Bool_t AliExternalTrackParam::GetPxPyPz(Double_t *p) const {
633   //---------------------------------------------------------------------
634   // This function returns the global track momentum components
635   // Results for (nearly) straight tracks are meaningless !
636   //---------------------------------------------------------------------
637   p[0]=fP[4]; p[1]=fP[2]; p[2]=fP[3];
638   return Local2GlobalMomentum(p,fAlpha);
639 }
640
641 Bool_t AliExternalTrackParam::GetXYZ(Double_t *r) const {
642   //---------------------------------------------------------------------
643   // This function returns the global track position
644   //---------------------------------------------------------------------
645   r[0]=fX; r[1]=fP[0]; r[2]=fP[1];
646   return Local2GlobalPosition(r,fAlpha);
647 }
648
649 Bool_t AliExternalTrackParam::GetCovarianceXYZPxPyPz(Double_t cv[21]) const {
650   //---------------------------------------------------------------------
651   // This function returns the global covariance matrix of the track params
652   // 
653   // Cov(x,x) ... :   cv[0]
654   // Cov(y,x) ... :   cv[1]  cv[2]
655   // Cov(z,x) ... :   cv[3]  cv[4]  cv[5]
656   // Cov(px,x)... :   cv[6]  cv[7]  cv[8]  cv[9]
657   // Cov(py,x)... :   cv[10] cv[11] cv[12] cv[13] cv[14]
658   // Cov(pz,x)... :   cv[15] cv[16] cv[17] cv[18] cv[19] cv[20]
659   //
660   // Results for (nearly) straight tracks are meaningless !
661   //---------------------------------------------------------------------
662   if (TMath::Abs(fP[4])<=0) {
663      for (Int_t i=0; i<21; i++) cv[i]=0.;
664      return kFALSE;
665   }
666   if (TMath::Abs(fP[2]) > kAlmost1) {
667      for (Int_t i=0; i<21; i++) cv[i]=0.;
668      return kFALSE;
669   }
670   Double_t pt=1./TMath::Abs(fP[4]);
671   Double_t cs=TMath::Cos(fAlpha), sn=TMath::Sin(fAlpha);
672   Double_t r=TMath::Sqrt(1-fP[2]*fP[2]);
673
674   Double_t m00=-sn, m10=cs;
675   Double_t m23=-pt*(sn + fP[2]*cs/r), m43=-pt*pt*(r*cs - fP[2]*sn);
676   Double_t m24= pt*(cs - fP[2]*sn/r), m44=-pt*pt*(r*sn + fP[2]*cs);
677   Double_t m35=pt, m45=-pt*pt*fP[3];
678
679   cv[0 ] = fC[0]*m00*m00;
680   cv[1 ] = fC[0]*m00*m10; 
681   cv[2 ] = fC[0]*m10*m10;
682   cv[3 ] = fC[1]*m00; 
683   cv[4 ] = fC[1]*m10; 
684   cv[5 ] = fC[2];
685   cv[6 ] = m00*(fC[3]*m23 + fC[10]*m43); 
686   cv[7 ] = m10*(fC[3]*m23 + fC[10]*m43); 
687   cv[8 ] = fC[4]*m23 + fC[11]*m43; 
688   cv[9 ] = m23*(fC[5]*m23 + fC[12]*m43)  +  m43*(fC[12]*m23 + fC[14]*m43);
689   cv[10] = m00*(fC[3]*m24 + fC[10]*m44); 
690   cv[11] = m10*(fC[3]*m24 + fC[10]*m44); 
691   cv[12] = fC[4]*m24 + fC[11]*m44; 
692   cv[13] = m23*(fC[5]*m24 + fC[12]*m44)  +  m43*(fC[12]*m24 + fC[14]*m44);
693   cv[14] = m24*(fC[5]*m24 + fC[12]*m44)  +  m44*(fC[12]*m24 + fC[14]*m44);
694   cv[15] = m00*(fC[6]*m35 + fC[10]*m45); 
695   cv[16] = m10*(fC[6]*m35 + fC[10]*m45); 
696   cv[17] = fC[7]*m35 + fC[11]*m45; 
697   cv[18] = m23*(fC[8]*m35 + fC[12]*m45)  +  m43*(fC[13]*m35 + fC[14]*m45);
698   cv[19] = m24*(fC[8]*m35 + fC[12]*m45)  +  m44*(fC[13]*m35 + fC[14]*m45); 
699   cv[20] = m35*(fC[9]*m35 + fC[13]*m45)  +  m45*(fC[13]*m35 + fC[14]*m45);
700
701   return kTRUE;
702 }
703
704
705 Bool_t 
706 AliExternalTrackParam::GetPxPyPzAt(Double_t x, Double_t b, Double_t *p) const {
707   //---------------------------------------------------------------------
708   // This function returns the global track momentum extrapolated to
709   // the radial position "x" (cm) in the magnetic field "b" (kG)
710   //---------------------------------------------------------------------
711   p[0]=fP[4]; 
712   p[1]=fP[2]+(x-fX)*fP[4]*b*kB2C; 
713   p[2]=fP[3];
714   return Local2GlobalMomentum(p,fAlpha);
715 }
716
717 Bool_t 
718 AliExternalTrackParam::GetXYZAt(Double_t x, Double_t b, Double_t *r) const {
719   //---------------------------------------------------------------------
720   // This function returns the global track position extrapolated to
721   // the radial position "x" (cm) in the magnetic field "b" (kG)
722   //---------------------------------------------------------------------
723   Double_t dx=x-fX;
724   Double_t f1=fP[2], f2=f1 + dx*fP[4]*b*kB2C;
725
726   if (TMath::Abs(f2) >= kAlmost1) return kFALSE;
727   
728   Double_t r1=TMath::Sqrt(1.- f1*f1), r2=TMath::Sqrt(1.- f2*f2);
729   r[0] = x;
730   r[1] = fP[0] + dx*(f1+f2)/(r1+r2);
731   r[2] = fP[1] + dx*(f1+f2)/(f1*r2 + f2*r1)*fP[3];
732   return Local2GlobalPosition(r,fAlpha);
733 }
734
735 //_____________________________________________________________________________
736 void AliExternalTrackParam::Print(Option_t* /*option*/) const
737 {
738 // print the parameters and the covariance matrix
739
740   printf("AliExternalTrackParam: x = %-12g  alpha = %-12g\n", fX, fAlpha);
741   printf("  parameters: %12g %12g %12g %12g %12g\n",
742          fP[0], fP[1], fP[2], fP[3], fP[4]);
743   printf("  covariance: %12g\n", fC[0]);
744   printf("              %12g %12g\n", fC[1], fC[2]);
745   printf("              %12g %12g %12g\n", fC[3], fC[4], fC[5]);
746   printf("              %12g %12g %12g %12g\n", 
747          fC[6], fC[7], fC[8], fC[9]);
748   printf("              %12g %12g %12g %12g %12g\n", 
749          fC[10], fC[11], fC[12], fC[13], fC[14]);
750 }
751
752
753 Bool_t AliExternalTrackParam::PropagateTo(Double_t xToGo, Double_t b, Double_t mass, Double_t maxStep, Bool_t rotateTo){
754   //----------------------------------------------------------------
755   //
756   // Very expensive function !  Don't abuse it !
757   //
758   // Propagates this track to the plane X=xk (cm) 
759   // in the magnetic field "b" (kG),
760   // the correction for the material is included
761   //
762   //  Requires acces to geomanager
763   //
764   // mass     - mass used in propagation - used for energy loss correction
765   // maxStep  - maximal step for propagation
766   //----------------------------------------------------------------
767   const Double_t kEpsilon = 0.00001;
768   Double_t xpos     = GetX();
769   Double_t dir      = (xpos<xToGo) ? 1.:-1.;
770   //
771   while ( (xToGo-xpos)*dir > kEpsilon){
772     if (TMath::Abs(fP[2]) >= kAlmost1) { return kFALSE;}
773     Double_t step = dir*TMath::Min(TMath::Abs(xToGo-xpos), maxStep);
774     Double_t x    = xpos+step;
775     Double_t xyz0[3],xyz1[3],param[7];
776     GetXYZ(xyz0);   //starting global position
777     if (!GetXYZAt(x,b,xyz1)) return kFALSE;   // no prolongation
778     AliKalmanTrack::MeanMaterialBudget(xyz0,xyz1,param);        
779     if (!PropagateTo(x,b))  return kFALSE;
780
781     Double_t rho=param[0],x0=param[1],distance=param[4];
782     Double_t d=distance*rho/x0;
783
784     if (!CorrectForMaterial(d,x0,mass)) return kFALSE;
785     if (rotateTo){
786       GetXYZ(xyz0);   // global position
787       Double_t alphan = TMath::ATan2(xyz0[1], xyz0[0]);
788       if (!Rotate(alphan)) return kFALSE;
789     }
790     xpos = GetX();
791   }
792   return kTRUE;
793 }
794
795
796
797