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