]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliExternalTrackParam.cxx
Do not change the track parameters if the propagation fails (Yuri)
[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 <TMatrixDSym.h>
29 #include <TPolyMarker3D.h>
30 #include <TVector3.h>
31
32 #include "AliExternalTrackParam.h"
33 #include "AliVVertex.h"
34 #include "AliLog.h"
35
36 ClassImp(AliExternalTrackParam)
37
38 Double32_t AliExternalTrackParam::fgMostProbablePt=kMostProbablePt;
39  
40 //_____________________________________________________________________________
41 AliExternalTrackParam::AliExternalTrackParam() :
42   AliVTrack(),
43   fX(0),
44   fAlpha(0)
45 {
46   //
47   // default constructor
48   //
49   for (Int_t i = 0; i < 5; i++) fP[i] = 0;
50   for (Int_t i = 0; i < 15; i++) fC[i] = 0;
51 }
52
53 //_____________________________________________________________________________
54 AliExternalTrackParam::AliExternalTrackParam(const AliExternalTrackParam &track):
55   AliVTrack(track),
56   fX(track.fX),
57   fAlpha(track.fAlpha)
58 {
59   //
60   // copy constructor
61   //
62   for (Int_t i = 0; i < 5; i++) fP[i] = track.fP[i];
63   for (Int_t i = 0; i < 15; i++) fC[i] = track.fC[i];
64 }
65
66 //_____________________________________________________________________________
67 AliExternalTrackParam& AliExternalTrackParam::operator=(const AliExternalTrackParam &trkPar)
68 {
69   //
70   // assignment operator
71   //
72   
73   if (this!=&trkPar) {
74     AliVTrack::operator=(trkPar);
75     fX = trkPar.fX;
76     fAlpha = trkPar.fAlpha;
77
78     for (Int_t i = 0; i < 5; i++) fP[i] = trkPar.fP[i];
79     for (Int_t i = 0; i < 15; i++) fC[i] = trkPar.fC[i];
80   }
81
82   return *this;
83 }
84
85 //_____________________________________________________________________________
86 AliExternalTrackParam::AliExternalTrackParam(Double_t x, Double_t alpha, 
87                                              const Double_t param[5], 
88                                              const Double_t covar[15]) :
89   AliVTrack(),
90   fX(x),
91   fAlpha(alpha)
92 {
93   //
94   // create external track parameters from given arguments
95   //
96   for (Int_t i = 0; i < 5; i++)  fP[i] = param[i];
97   for (Int_t i = 0; i < 15; i++) fC[i] = covar[i];
98 }
99
100 //_____________________________________________________________________________
101 AliExternalTrackParam::AliExternalTrackParam(const AliVTrack *vTrack) :
102   AliVTrack(),
103   fX(0.),
104   fAlpha(0.)
105 {
106   //
107   // Constructor from virtual track,
108   // This is not a copy contructor !
109   //
110
111   if (vTrack->InheritsFrom("AliExternalTrackParam")) {
112      AliError("This is not a copy constructor. Use AliExternalTrackParam(const AliExternalTrackParam &) !");
113      AliWarning("Calling the default constructor...");
114      AliExternalTrackParam();
115      return;
116   }
117
118   Double_t xyz[3],pxpypz[3],cv[21];
119   vTrack->GetXYZ(xyz);
120   pxpypz[0]=vTrack->Px();
121   pxpypz[1]=vTrack->Py();
122   pxpypz[2]=vTrack->Pz();
123   vTrack->GetCovarianceXYZPxPyPz(cv);
124   Short_t sign = (Short_t)vTrack->Charge();
125
126   Set(xyz,pxpypz,cv,sign);
127 }
128
129 //_____________________________________________________________________________
130 AliExternalTrackParam::AliExternalTrackParam(Double_t xyz[3],Double_t pxpypz[3],
131                                              Double_t cv[21],Short_t sign) :
132   AliVTrack(),
133   fX(0.),
134   fAlpha(0.)
135 {
136   //
137   // constructor from the global parameters
138   //
139
140   Set(xyz,pxpypz,cv,sign);
141 }
142
143 //_____________________________________________________________________________
144 void AliExternalTrackParam::Set(Double_t xyz[3],Double_t pxpypz[3],
145                                 Double_t cv[21],Short_t sign) 
146 {
147   //
148   // create external track parameters from the global parameters
149   // x,y,z,px,py,pz and their 6x6 covariance matrix
150   // A.Dainese 10.10.08
151
152   // Calculate alpha: the rotation angle of the corresponding local system
153   fAlpha = TMath::ATan2(pxpypz[1],pxpypz[0]);
154
155   // Get the vertex of origin and the momentum
156   TVector3 ver(xyz[0],xyz[1],xyz[2]);
157   TVector3 mom(pxpypz[0],pxpypz[1],pxpypz[2]);
158
159   // Rotate to the local coordinate system
160   ver.RotateZ(-fAlpha);
161   mom.RotateZ(-fAlpha);
162
163   // x of the reference plane
164   fX = ver.X();
165
166   Double_t charge = (Double_t)sign;
167
168   fP[0] = ver.Y();
169   fP[1] = ver.Z();
170   fP[2] = TMath::Sin(mom.Phi());
171   fP[3] = mom.Pz()/mom.Pt();
172   fP[4] = TMath::Sign(1/mom.Pt(),charge);
173
174   // Covariance matrix (formulas to be simplified)
175
176   Double_t pt=1./TMath::Abs(fP[4]);
177   Double_t cs=TMath::Cos(fAlpha), sn=TMath::Sin(fAlpha);
178   Double_t r=TMath::Sqrt((1.-fP[2])*(1.+fP[2]));
179
180   Double_t m00=-sn;// m10=cs;
181   Double_t m23=-pt*(sn + fP[2]*cs/r), m43=-pt*pt*(r*cs - fP[2]*sn);
182   Double_t m24= pt*(cs - fP[2]*sn/r), m44=-pt*pt*(r*sn + fP[2]*cs);
183   Double_t m35=pt, m45=-pt*pt*fP[3];
184
185   m43*=GetSign();
186   m44*=GetSign();
187   m45*=GetSign();
188
189   Double_t cv34 = TMath::Sqrt(cv[3 ]*cv[3 ]+cv[4 ]*cv[4 ]);
190   Double_t a1=cv[13]-cv[9]*(m23*m44+m43*m24)/m23/m43;
191   Double_t a2=m23*m24-m23*(m23*m44+m43*m24)/m43;
192   Double_t a3=m43*m44-m43*(m23*m44+m43*m24)/m23;
193   Double_t a4=cv[14]-2.*cv[9]*m24*m44/m23/m43;
194   Double_t a5=m24*m24-2.*m24*m44*m23/m43;
195   Double_t a6=m44*m44-2.*m24*m44*m43/m23;
196
197   fC[0 ] = cv[0 ]+cv[2 ];  
198   fC[1 ] = TMath::Sign(cv34,cv[3 ]/m00); 
199   fC[2 ] = cv[5 ]; 
200   fC[3 ] = (cv[10]/m44-cv[6]/m43)/(m24/m44-m23/m43)/m00; 
201   fC[10] = (cv[6]/m00-fC[3 ]*m23)/m43; 
202   fC[6 ] = (cv[15]/m00-fC[10]*m45)/m35; 
203   fC[4 ] = (cv[12]-cv[8]*m44/m43)/(m24-m23*m44/m43); 
204   fC[11] = (cv[8]-fC[4]*m23)/m43; 
205   fC[7 ] = cv[17]/m35-fC[11]*m45/m35; 
206   fC[5 ] = TMath::Abs((a4-a6*a1/a3)/(a5-a6*a2/a3));
207   fC[14] = TMath::Abs(a1/a3-a2*fC[5]/a3);
208   fC[12] = (cv[9]-fC[5]*m23*m23-fC[14]*m43*m43)/m23/m43;
209   Double_t b1=cv[18]-fC[12]*m23*m45-fC[14]*m43*m45;
210   Double_t b2=m23*m35;
211   Double_t b3=m43*m35;
212   Double_t b4=cv[19]-fC[12]*m24*m45-fC[14]*m44*m45;
213   Double_t b5=m24*m35;
214   Double_t b6=m44*m35;
215   fC[8 ] = (b4-b6*b1/b3)/(b5-b6*b2/b3);
216   fC[13] = b1/b3-b2*fC[8]/b3;
217   fC[9 ] = TMath::Abs((cv[20]-fC[14]*(m45*m45)-fC[13]*2.*m35*m45)/(m35*m35));
218
219   return;
220 }
221
222 //_____________________________________________________________________________
223 void AliExternalTrackParam::Reset() {
224   //
225   // Resets all the parameters to 0 
226   //
227   fX=fAlpha=0.;
228   for (Int_t i = 0; i < 5; i++) fP[i] = 0;
229   for (Int_t i = 0; i < 15; i++) fC[i] = 0;
230 }
231
232 //_____________________________________________________________________________
233 void AliExternalTrackParam::AddCovariance(const Double_t c[15]) {
234   //
235   // Add "something" to the track covarince matrix.
236   // May be needed to account for unknown mis-calibration/mis-alignment
237   //
238     fC[0] +=c[0];
239     fC[1] +=c[1];  fC[2] +=c[2];
240     fC[3] +=c[3];  fC[4] +=c[4];  fC[5] +=c[5];
241     fC[6] +=c[6];  fC[7] +=c[7];  fC[8] +=c[8];  fC[9] +=c[9];
242     fC[10]+=c[10]; fC[11]+=c[11]; fC[12]+=c[12]; fC[13]+=c[13]; fC[14]+=c[14];
243 }
244
245
246 Double_t AliExternalTrackParam::GetP() const {
247   //---------------------------------------------------------------------
248   // This function returns the track momentum
249   // Results for (nearly) straight tracks are meaningless !
250   //---------------------------------------------------------------------
251   if (TMath::Abs(fP[4])<=kAlmost0) return kVeryBig;
252   return TMath::Sqrt(1.+ fP[3]*fP[3])/TMath::Abs(fP[4]);
253 }
254
255 Double_t AliExternalTrackParam::Get1P() const {
256   //---------------------------------------------------------------------
257   // This function returns the 1/(track momentum)
258   //---------------------------------------------------------------------
259   return TMath::Abs(fP[4])/TMath::Sqrt(1.+ fP[3]*fP[3]);
260 }
261
262 //_______________________________________________________________________
263 Double_t AliExternalTrackParam::GetD(Double_t x,Double_t y,Double_t b) const {
264   //------------------------------------------------------------------
265   // This function calculates the transverse impact parameter
266   // with respect to a point with global coordinates (x,y)
267   // in the magnetic field "b" (kG)
268   //------------------------------------------------------------------
269   if (TMath::Abs(b) < kAlmost0Field) return GetLinearD(x,y);
270   Double_t rp4=GetC(b);
271
272   Double_t xt=fX, yt=fP[0];
273
274   Double_t sn=TMath::Sin(fAlpha), cs=TMath::Cos(fAlpha);
275   Double_t a = x*cs + y*sn;
276   y = -x*sn + y*cs; x=a;
277   xt-=x; yt-=y;
278
279   sn=rp4*xt - fP[2]; cs=rp4*yt + TMath::Sqrt(1.- fP[2]*fP[2]);
280   a=2*(xt*fP[2] - yt*TMath::Sqrt(1.- fP[2]*fP[2]))-rp4*(xt*xt + yt*yt);
281   return  -a/(1 + TMath::Sqrt(sn*sn + cs*cs));
282 }
283
284 //_______________________________________________________________________
285 void AliExternalTrackParam::
286 GetDZ(Double_t x, Double_t y, Double_t z, Double_t b, Float_t dz[2]) const {
287   //------------------------------------------------------------------
288   // This function calculates the transverse and longitudinal impact parameters
289   // with respect to a point with global coordinates (x,y)
290   // in the magnetic field "b" (kG)
291   //------------------------------------------------------------------
292   Double_t f1 = fP[2], r1 = TMath::Sqrt(1. - f1*f1);
293   Double_t xt=fX, yt=fP[0];
294   Double_t sn=TMath::Sin(fAlpha), cs=TMath::Cos(fAlpha);
295   Double_t a = x*cs + y*sn;
296   y = -x*sn + y*cs; x=a;
297   xt-=x; yt-=y;
298
299   Double_t rp4=GetC(b);
300   if ((TMath::Abs(b) < kAlmost0Field) || (TMath::Abs(rp4) < kAlmost0)) {
301      dz[0] = -(xt*f1 - yt*r1);
302      dz[1] = fP[1] + (dz[0]*f1 - xt)/r1*fP[3] - z;
303      return;
304   }
305
306   sn=rp4*xt - f1; cs=rp4*yt + r1;
307   a=2*(xt*f1 - yt*r1)-rp4*(xt*xt + yt*yt);
308   Double_t rr=TMath::Sqrt(sn*sn + cs*cs);
309   dz[0] = -a/(1 + rr);
310   Double_t f2 = -sn/rr, r2 = TMath::Sqrt(1. - f2*f2);
311   dz[1] = fP[1] + fP[3]/rp4*TMath::ASin(f2*r1 - f1*r2) - z;
312 }
313
314 //_______________________________________________________________________
315 Double_t AliExternalTrackParam::GetLinearD(Double_t xv,Double_t yv) const {
316   //------------------------------------------------------------------
317   // This function calculates the transverse impact parameter
318   // with respect to a point with global coordinates (xv,yv)
319   // neglecting the track curvature.
320   //------------------------------------------------------------------
321   Double_t sn=TMath::Sin(fAlpha), cs=TMath::Cos(fAlpha);
322   Double_t x= xv*cs + yv*sn;
323   Double_t y=-xv*sn + yv*cs;
324
325   Double_t d = (fX-x)*fP[2] - (fP[0]-y)*TMath::Sqrt(1.- fP[2]*fP[2]);
326
327   return -d;
328 }
329
330 Bool_t AliExternalTrackParam::CorrectForMeanMaterial
331 (Double_t xOverX0,  Double_t xTimesRho, Double_t mass, Bool_t anglecorr, 
332  Double_t (*Bethe)(Double_t)) {
333   //------------------------------------------------------------------
334   // This function corrects the track parameters for the crossed material.
335   // "xOverX0"   - X/X0, the thickness in units of the radiation length.
336   // "xTimesRho" - is the product length*density (g/cm^2). 
337   // "mass" - the mass of this particle (GeV/c^2).
338   //------------------------------------------------------------------
339   Double_t &fP2=fP[2];
340   Double_t &fP3=fP[3];
341   Double_t &fP4=fP[4];
342
343   Double_t &fC22=fC[5];
344   Double_t &fC33=fC[9];
345   Double_t &fC43=fC[13];
346   Double_t &fC44=fC[14];
347
348   //Apply angle correction, if requested
349   if(anglecorr) {
350     Double_t angle=TMath::Sqrt((1.+ fP3*fP3)/(1.- fP2*fP2));
351     xOverX0 *=angle;
352     xTimesRho *=angle;
353   } 
354
355   Double_t p=GetP();
356   Double_t p2=p*p;
357   Double_t beta2=p2/(p2 + mass*mass);
358
359   //Calculating the multiple scattering corrections******************
360   Double_t cC22 = 0.;
361   Double_t cC33 = 0.;
362   Double_t cC43 = 0.;
363   Double_t cC44 = 0.;
364   if (xOverX0 != 0) {
365      Double_t theta2=14.1*14.1/(beta2*p2*1e6)*TMath::Abs(xOverX0);
366      //Double_t theta2=1.0259e-6*14*14/28/(beta2*p2)*TMath::Abs(d)*9.36*2.33;
367      if(theta2>TMath::Pi()*TMath::Pi()) return kFALSE;
368      cC22 = theta2*(1.- fP2*fP2)*(1. + fP3*fP3);
369      cC33 = theta2*(1. + fP3*fP3)*(1. + fP3*fP3);
370      cC43 = theta2*fP3*fP4*(1. + fP3*fP3);
371      cC44 = theta2*fP3*fP4*fP3*fP4;
372   }
373
374   //Calculating the energy loss corrections************************
375   Double_t cP4=1.;
376   if ((xTimesRho != 0.) && (beta2 < 1.)) {
377      Double_t dE=Bethe(p/mass)*xTimesRho;
378      Double_t e=TMath::Sqrt(p2 + mass*mass);
379      if ( TMath::Abs(dE) > 0.3*e ) return kFALSE; //30% energy loss is too much!
380      cP4 = (1.- e/p2*dE);
381      if (TMath::Abs(fP4*cP4)>100.) return kFALSE; //Do not track below 10 MeV/c
382
383
384      // Approximate energy loss fluctuation (M.Ivanov)
385      const Double_t knst=0.07; // To be tuned.  
386      Double_t sigmadE=knst*TMath::Sqrt(TMath::Abs(dE)); 
387      cC44 += ((sigmadE*e/p2*fP4)*(sigmadE*e/p2*fP4)); 
388  
389   }
390
391   //Applying the corrections*****************************
392   fC22 += cC22;
393   fC33 += cC33;
394   fC43 += cC43;
395   fC44 += cC44;
396   fP4  *= cP4;
397
398   return kTRUE;
399 }
400
401
402 Bool_t AliExternalTrackParam::CorrectForMaterial
403 (Double_t d,  Double_t x0, Double_t mass, Double_t (*Bethe)(Double_t)) {
404   //------------------------------------------------------------------
405   //                    Deprecated function !   
406   //       Better use CorrectForMeanMaterial instead of it.
407   //
408   // This function corrects the track parameters for the crossed material
409   // "d"    - the thickness (fraction of the radiation length)
410   // "x0"   - the radiation length (g/cm^2) 
411   // "mass" - the mass of this particle (GeV/c^2)
412   //------------------------------------------------------------------
413   Double_t &fP2=fP[2];
414   Double_t &fP3=fP[3];
415   Double_t &fP4=fP[4];
416
417   Double_t &fC22=fC[5];
418   Double_t &fC33=fC[9];
419   Double_t &fC43=fC[13];
420   Double_t &fC44=fC[14];
421
422   Double_t p=GetP();
423   Double_t p2=p*p;
424   Double_t beta2=p2/(p2 + mass*mass);
425   d*=TMath::Sqrt((1.+ fP3*fP3)/(1.- fP2*fP2));
426
427   //Multiple scattering******************
428   Double_t cC22 = 0.;
429   Double_t cC33 = 0.;
430   Double_t cC43 = 0.;
431   Double_t cC44 = 0.;
432   if (d!=0) {
433      Double_t theta2=14.1*14.1/(beta2*p2*1e6)*TMath::Abs(d);
434      //Double_t theta2=1.0259e-6*14*14/28/(beta2*p2)*TMath::Abs(d)*9.36*2.33;
435      if(theta2>TMath::Pi()*TMath::Pi()) return kFALSE;
436      cC22 = theta2*(1.- fP2*fP2)*(1. + fP3*fP3);
437      cC33 = theta2*(1. + fP3*fP3)*(1. + fP3*fP3);
438      cC43 = theta2*fP3*fP4*(1. + fP3*fP3);
439      cC44 = theta2*fP3*fP4*fP3*fP4;
440   }
441
442   //Energy losses************************
443   Double_t cP4=1.;
444   if (x0!=0. && beta2<1) {
445      d*=x0;
446      Double_t dE=Bethe(p/mass)*d;
447      Double_t e=TMath::Sqrt(p2 + mass*mass);
448      if ( TMath::Abs(dE) > 0.3*e ) return kFALSE; //30% energy loss is too much!
449      cP4 = (1.- e/p2*dE);
450
451      // Approximate energy loss fluctuation (M.Ivanov)
452      const Double_t knst=0.07; // To be tuned.  
453      Double_t sigmadE=knst*TMath::Sqrt(TMath::Abs(dE)); 
454      cC44 += ((sigmadE*e/p2*fP4)*(sigmadE*e/p2*fP4)); 
455  
456   }
457
458   fC22 += cC22;
459   fC33 += cC33;
460   fC43 += cC43;
461   fC44 += cC44;
462   fP4  *= cP4;
463
464   return kTRUE;
465 }
466
467 Double_t AliExternalTrackParam::BetheBlochAleph(Double_t bg,
468          Double_t kp1,
469          Double_t kp2,
470          Double_t kp3,
471          Double_t kp4,
472          Double_t kp5) {
473   //
474   // This is the empirical ALEPH parameterization of the Bethe-Bloch formula.
475   // It is normalized to 1 at the minimum.
476   //
477   // bg - beta*gamma
478   // 
479   // The default values for the kp* parameters are for ALICE TPC.
480   // The returned value is in MIP units
481   //
482
483   Double_t beta = bg/TMath::Sqrt(1.+ bg*bg);
484
485   Double_t aa = TMath::Power(beta,kp4);
486   Double_t bb = TMath::Power(1./bg,kp5);
487
488   bb=TMath::Log(kp3+bb);
489   
490   return (kp2-aa-bb)*kp1/aa;
491 }
492
493 Double_t AliExternalTrackParam::BetheBlochGeant(Double_t bg,
494          Double_t kp0,
495          Double_t kp1,
496          Double_t kp2,
497          Double_t kp3,
498          Double_t kp4) {
499   //
500   // This is the parameterization of the Bethe-Bloch formula inspired by Geant.
501   //
502   // bg  - beta*gamma
503   // kp0 - density [g/cm^3]
504   // kp1 - density effect first junction point
505   // kp2 - density effect second junction point
506   // kp3 - mean excitation energy [GeV]
507   // kp4 - mean Z/A
508   //
509   // The default values for the kp* parameters are for silicon. 
510   // The returned value is in [GeV/(g/cm^2)].
511   // 
512
513   const Double_t mK  = 0.307075e-3; // [GeV*cm^2/g]
514   const Double_t me  = 0.511e-3;    // [GeV/c^2]
515   const Double_t rho = kp0;
516   const Double_t x0  = kp1*2.303;
517   const Double_t x1  = kp2*2.303;
518   const Double_t mI  = kp3;
519   const Double_t mZA = kp4;
520   const Double_t bg2 = bg*bg;
521   const Double_t maxT= 2*me*bg2;    // neglecting the electron mass
522   
523   //*** Density effect
524   Double_t d2=0.; 
525   const Double_t x=TMath::Log(bg);
526   const Double_t lhwI=TMath::Log(28.816*1e-9*TMath::Sqrt(rho*mZA)/mI);
527   if (x > x1) {
528     d2 = lhwI + x - 0.5;
529   } else if (x > x0) {
530     const Double_t r=(x1-x)/(x1-x0);
531     d2 = lhwI + x - 0.5 + (0.5 - lhwI - x0)*r*r*r;
532   }
533
534   return mK*mZA*(1+bg2)/bg2*
535          (0.5*TMath::Log(2*me*bg2*maxT/(mI*mI)) - bg2/(1+bg2) - d2);
536 }
537
538 Double_t AliExternalTrackParam::BetheBlochSolid(Double_t bg) {
539   //------------------------------------------------------------------
540   // This is an approximation of the Bethe-Bloch formula, 
541   // reasonable for solid materials. 
542   // All the parameters are, in fact, for Si.
543   // The returned value is in [GeV]
544   //------------------------------------------------------------------
545
546   return BetheBlochGeant(bg);
547 }
548
549 Double_t AliExternalTrackParam::BetheBlochGas(Double_t bg) {
550   //------------------------------------------------------------------
551   // This is an approximation of the Bethe-Bloch formula, 
552   // reasonable for gas materials.
553   // All the parameters are, in fact, for Ne.
554   // The returned value is in [GeV]
555   //------------------------------------------------------------------
556
557   const Double_t rho = 0.9e-3;
558   const Double_t x0  = 2.;
559   const Double_t x1  = 4.;
560   const Double_t mI  = 140.e-9;
561   const Double_t mZA = 0.49555;
562
563   return BetheBlochGeant(bg,rho,x0,x1,mI,mZA);
564 }
565
566 Bool_t AliExternalTrackParam::Rotate(Double_t alpha) {
567   //------------------------------------------------------------------
568   // Transform this track to the local coord. system rotated
569   // by angle "alpha" (rad) with respect to the global coord. system. 
570   //------------------------------------------------------------------
571   if (TMath::Abs(fP[2]) >= kAlmost1) {
572      AliError(Form("Precondition is not satisfied: |sin(phi)|>1 ! %f",fP[2])); 
573      return kFALSE;
574   }
575  
576   if      (alpha < -TMath::Pi()) alpha += 2*TMath::Pi();
577   else if (alpha >= TMath::Pi()) alpha -= 2*TMath::Pi();
578
579   Double_t &fP0=fP[0];
580   Double_t &fP2=fP[2];
581   Double_t &fC00=fC[0];
582   Double_t &fC10=fC[1];
583   Double_t &fC20=fC[3];
584   Double_t &fC21=fC[4];
585   Double_t &fC22=fC[5];
586   Double_t &fC30=fC[6];
587   Double_t &fC32=fC[8];
588   Double_t &fC40=fC[10];
589   Double_t &fC42=fC[12];
590
591   Double_t x=fX;
592   Double_t ca=TMath::Cos(alpha-fAlpha), sa=TMath::Sin(alpha-fAlpha);
593   Double_t sf=fP2, cf=TMath::Sqrt(1.- fP2*fP2);
594
595   Double_t tmp=sf*ca - cf*sa;
596   if (TMath::Abs(tmp) >= kAlmost1) return kFALSE;
597
598   fAlpha = alpha;
599   fX =  x*ca + fP0*sa;
600   fP0= -x*sa + fP0*ca;
601   fP2=  tmp;
602
603   if (TMath::Abs(cf)<kAlmost0) {
604     AliError(Form("Too small cosine value %f",cf)); 
605     cf = kAlmost0;
606   } 
607
608   Double_t rr=(ca+sf/cf*sa);  
609
610   fC00 *= (ca*ca);
611   fC10 *= ca;
612   fC20 *= ca*rr;
613   fC21 *= rr;
614   fC22 *= rr*rr;
615   fC30 *= ca;
616   fC32 *= rr;
617   fC40 *= ca;
618   fC42 *= rr;
619
620   return kTRUE;
621 }
622
623 Bool_t AliExternalTrackParam::PropagateTo(Double_t xk, Double_t b) {
624   //----------------------------------------------------------------
625   // Propagate this track to the plane X=xk (cm) in the field "b" (kG)
626   //----------------------------------------------------------------
627   Double_t dx=xk-fX;
628   if (TMath::Abs(dx)<=kAlmost0)  return kTRUE;
629
630   Double_t crv=GetC(b);
631   if (TMath::Abs(b) < kAlmost0Field) crv=0.;
632
633   Double_t f1=fP[2], f2=f1 + crv*dx;
634   if (TMath::Abs(f1) >= kAlmost1) return kFALSE;
635   if (TMath::Abs(f2) >= kAlmost1) return kFALSE;
636
637   Double_t &fP0=fP[0], &fP1=fP[1], &fP2=fP[2], &fP3=fP[3], &fP4=fP[4];
638   Double_t 
639   &fC00=fC[0],
640   &fC10=fC[1],   &fC11=fC[2],  
641   &fC20=fC[3],   &fC21=fC[4],   &fC22=fC[5],
642   &fC30=fC[6],   &fC31=fC[7],   &fC32=fC[8],   &fC33=fC[9],  
643   &fC40=fC[10],  &fC41=fC[11],  &fC42=fC[12],  &fC43=fC[13], &fC44=fC[14];
644
645   Double_t r1=TMath::Sqrt(1.- f1*f1), r2=TMath::Sqrt(1.- f2*f2);
646
647   fX=xk;
648   fP0 += dx*(f1+f2)/(r1+r2);
649   fP1 += dx*(r2 + f2*(f1+f2)/(r1+r2))*fP3;  // Many thanks to P.Hristov !
650   fP2 += dx*crv;
651
652   //f = F - 1
653    
654   Double_t f02=    dx/(r1*r1*r1);            Double_t cc=crv/fP4;
655   Double_t f04=0.5*dx*dx/(r1*r1*r1);         f04*=cc;
656   Double_t f12=    dx*fP3*f1/(r1*r1*r1);
657   Double_t f14=0.5*dx*dx*fP3*f1/(r1*r1*r1);  f14*=cc;
658   Double_t f13=    dx/r1;
659   Double_t f24=    dx;                       f24*=cc;
660   
661   //b = C*ft
662   Double_t b00=f02*fC20 + f04*fC40, b01=f12*fC20 + f14*fC40 + f13*fC30;
663   Double_t b02=f24*fC40;
664   Double_t b10=f02*fC21 + f04*fC41, b11=f12*fC21 + f14*fC41 + f13*fC31;
665   Double_t b12=f24*fC41;
666   Double_t b20=f02*fC22 + f04*fC42, b21=f12*fC22 + f14*fC42 + f13*fC32;
667   Double_t b22=f24*fC42;
668   Double_t b40=f02*fC42 + f04*fC44, b41=f12*fC42 + f14*fC44 + f13*fC43;
669   Double_t b42=f24*fC44;
670   Double_t b30=f02*fC32 + f04*fC43, b31=f12*fC32 + f14*fC43 + f13*fC33;
671   Double_t b32=f24*fC43;
672   
673   //a = f*b = f*C*ft
674   Double_t a00=f02*b20+f04*b40,a01=f02*b21+f04*b41,a02=f02*b22+f04*b42;
675   Double_t a11=f12*b21+f14*b41+f13*b31,a12=f12*b22+f14*b42+f13*b32;
676   Double_t a22=f24*b42;
677
678   //F*C*Ft = C + (b + bt + a)
679   fC00 += b00 + b00 + a00;
680   fC10 += b10 + b01 + a01; 
681   fC20 += b20 + b02 + a02;
682   fC30 += b30;
683   fC40 += b40;
684   fC11 += b11 + b11 + a11;
685   fC21 += b21 + b12 + a12;
686   fC31 += b31; 
687   fC41 += b41;
688   fC22 += b22 + b22 + a22;
689   fC32 += b32;
690   fC42 += b42;
691
692   return kTRUE;
693 }
694
695 Bool_t 
696 AliExternalTrackParam::Propagate(Double_t alpha, Double_t x, Double_t b) {
697   //------------------------------------------------------------------
698   // Transform this track to the local coord. system rotated
699   // by angle "alpha" (rad) with respect to the global coord. system, 
700   // and propagate this track to the plane X=xk (cm) in the field "b" (kG)
701   //------------------------------------------------------------------
702   
703   //Save the parameters
704   Double_t as=fAlpha;
705   Double_t xs=fX;
706   Double_t ps[5], cs[15];
707   for (Int_t i=0; i<5;  i++) ps[i]=fP[i]; 
708   for (Int_t i=0; i<15; i++) cs[i]=fC[i]; 
709
710   if (Rotate(alpha))
711      if (PropagateTo(x,b)) return kTRUE;
712
713   //Restore the parameters, if the operation failed
714   fAlpha=as;
715   fX=xs;
716   for (Int_t i=0; i<5;  i++) fP[i]=ps[i]; 
717   for (Int_t i=0; i<15; i++) fC[i]=cs[i]; 
718   return kFALSE;
719 }
720
721
722 void AliExternalTrackParam::Propagate(Double_t len, Double_t x[3],
723 Double_t p[3], Double_t bz) const {
724   //+++++++++++++++++++++++++++++++++++++++++    
725   // Origin: K. Shileev (Kirill.Shileev@cern.ch)
726   // Extrapolate track along simple helix in magnetic field
727   // Arguments: len -distance alogn helix, [cm]
728   //            bz  - mag field, [kGaus]   
729   // Returns: x and p contain extrapolated positon and momentum  
730   // The momentum returned for straight-line tracks is meaningless !
731   //+++++++++++++++++++++++++++++++++++++++++    
732   GetXYZ(x);
733     
734   if (OneOverPt() < kAlmost0 || TMath::Abs(bz) < kAlmost0Field || GetC(bz) < kAlmost0){ //straight-line tracks
735      Double_t unit[3]; GetDirection(unit);
736      x[0]+=unit[0]*len;   
737      x[1]+=unit[1]*len;   
738      x[2]+=unit[2]*len;
739
740      p[0]=unit[0]/kAlmost0;   
741      p[1]=unit[1]/kAlmost0;   
742      p[2]=unit[2]/kAlmost0;   
743   } else {
744      GetPxPyPz(p);
745      Double_t pp=GetP();
746      Double_t a = -kB2C*bz*GetSign();
747      Double_t rho = a/pp;
748      x[0] += p[0]*TMath::Sin(rho*len)/a - p[1]*(1-TMath::Cos(rho*len))/a;
749      x[1] += p[1]*TMath::Sin(rho*len)/a + p[0]*(1-TMath::Cos(rho*len))/a;
750      x[2] += p[2]*len/pp;
751
752      Double_t p0=p[0];
753      p[0] = p0  *TMath::Cos(rho*len) - p[1]*TMath::Sin(rho*len);
754      p[1] = p[1]*TMath::Cos(rho*len) + p0  *TMath::Sin(rho*len);
755   }
756 }
757
758 Bool_t AliExternalTrackParam::Intersect(Double_t pnt[3], Double_t norm[3],
759 Double_t bz) const {
760   //+++++++++++++++++++++++++++++++++++++++++    
761   // Origin: K. Shileev (Kirill.Shileev@cern.ch)
762   // Finds point of intersection (if exists) of the helix with the plane. 
763   // Stores result in fX and fP.   
764   // Arguments: planePoint,planeNorm - the plane defined by any plane's point 
765   // and vector, normal to the plane
766   // Returns: kTrue if helix intersects the plane, kFALSE otherwise.
767   //+++++++++++++++++++++++++++++++++++++++++    
768   Double_t x0[3]; GetXYZ(x0); //get track position in MARS
769   
770   //estimates initial helix length up to plane
771   Double_t s=
772     (pnt[0]-x0[0])*norm[0] + (pnt[1]-x0[1])*norm[1] + (pnt[2]-x0[2])*norm[2];
773   Double_t dist=99999,distPrev=dist;
774   Double_t x[3],p[3]; 
775   while(TMath::Abs(dist)>0.00001){
776     //calculates helix at the distance s from x0 ALONG the helix
777     Propagate(s,x,p,bz);
778
779     //distance between current helix position and plane
780     dist=(x[0]-pnt[0])*norm[0]+(x[1]-pnt[1])*norm[1]+(x[2]-pnt[2])*norm[2];
781
782     if(TMath::Abs(dist) >= TMath::Abs(distPrev)) {return kFALSE;}
783     distPrev=dist;
784     s-=dist;
785   }
786   //on exit pnt is intersection point,norm is track vector at that point, 
787   //all in MARS
788   for (Int_t i=0; i<3; i++) {pnt[i]=x[i]; norm[i]=p[i];}
789   return kTRUE;
790 }
791
792 Double_t 
793 AliExternalTrackParam::GetPredictedChi2(Double_t p[2],Double_t cov[3]) const {
794   //----------------------------------------------------------------
795   // Estimate the chi2 of the space point "p" with the cov. matrix "cov"
796   //----------------------------------------------------------------
797   Double_t sdd = fC[0] + cov[0]; 
798   Double_t sdz = fC[1] + cov[1];
799   Double_t szz = fC[2] + cov[2];
800   Double_t det = sdd*szz - sdz*sdz;
801
802   if (TMath::Abs(det) < kAlmost0) return kVeryBig;
803
804   Double_t d = fP[0] - p[0];
805   Double_t z = fP[1] - p[1];
806
807   return (d*szz*d - 2*d*sdz*z + z*sdd*z)/det;
808 }
809
810 Double_t AliExternalTrackParam::
811 GetPredictedChi2(Double_t p[3],Double_t covyz[3],Double_t covxyz[3]) const {
812   //----------------------------------------------------------------
813   // Estimate the chi2 of the 3D space point "p" and
814   // the full covariance matrix "covyz" and "covxyz"
815   //
816   // Cov(x,x) ... :   covxyz[0]
817   // Cov(y,x) ... :   covxyz[1]  covyz[0]
818   // Cov(z,x) ... :   covxyz[2]  covyz[1]  covyz[2]
819   //----------------------------------------------------------------
820
821   Double_t res[3] = {
822     GetX() - p[0],
823     GetY() - p[1],
824     GetZ() - p[2]
825   };
826
827   Double_t f=GetSnp();
828   if (TMath::Abs(f) >= kAlmost1) return kVeryBig;
829   Double_t r=TMath::Sqrt(1.- f*f);
830   Double_t a=f/r, b=GetTgl()/r;
831
832   Double_t s2=333.*333.;  //something reasonably big (cm^2)
833  
834   TMatrixDSym v(3);
835   v(0,0)=  s2;  v(0,1)=  a*s2;                 v(0,2)=  b*s2;;
836   v(1,0)=a*s2;  v(1,1)=a*a*s2 + GetSigmaY2();  v(1,2)=a*b*s2 + GetSigmaZY();
837   v(2,0)=b*s2;  v(2,1)=a*b*s2 + GetSigmaZY();  v(2,2)=b*b*s2 + GetSigmaZ2();
838
839   v(0,0)+=covxyz[0]; v(0,1)+=covxyz[1]; v(0,2)+=covxyz[2];
840   v(1,0)+=covxyz[1]; v(1,1)+=covyz[0];  v(1,2)+=covyz[1];
841   v(2,0)+=covxyz[2]; v(2,1)+=covyz[1];  v(2,2)+=covyz[2];
842
843   v.Invert();
844   if (!v.IsValid()) return kVeryBig;
845
846   Double_t chi2=0.;
847   for (Int_t i = 0; i < 3; i++)
848     for (Int_t j = 0; j < 3; j++) chi2 += res[i]*res[j]*v(i,j);
849
850   return chi2;  
851
852
853 }
854
855 Bool_t AliExternalTrackParam::
856 PropagateTo(Double_t p[3],Double_t covyz[3],Double_t covxyz[3],Double_t bz) {
857   //----------------------------------------------------------------
858   // Propagate this track to the plane 
859   // the 3D space point "p" (with the covariance matrix "covyz" and "covxyz")
860   // belongs to.
861   // The magnetic field is "bz" (kG)
862   //
863   // The track curvature and the change of the covariance matrix
864   // of the track parameters are negleted !
865   // (So the "step" should be small compared with 1/curvature)
866   //----------------------------------------------------------------
867
868   Double_t f=GetSnp();
869   if (TMath::Abs(f) >= kAlmost1) return kFALSE;
870   Double_t r=TMath::Sqrt(1.- f*f);
871   Double_t a=f/r, b=GetTgl()/r;
872
873   Double_t s2=333.*333.;  //something reasonably big (cm^2)
874  
875   TMatrixDSym tV(3);
876   tV(0,0)=  s2;  tV(0,1)=  a*s2;  tV(0,2)=  b*s2;
877   tV(1,0)=a*s2;  tV(1,1)=a*a*s2;  tV(1,2)=a*b*s2;
878   tV(2,0)=b*s2;  tV(2,1)=a*b*s2;  tV(2,2)=b*b*s2;
879
880   TMatrixDSym pV(3);
881   pV(0,0)=covxyz[0]; pV(0,1)=covxyz[1]; pV(0,2)=covxyz[2];
882   pV(1,0)=covxyz[1]; pV(1,1)=covyz[0];  pV(1,2)=covyz[1];
883   pV(2,0)=covxyz[2]; pV(2,1)=covyz[1];  pV(2,2)=covyz[2];
884
885   TMatrixDSym tpV(tV);
886   tpV+=pV;
887   tpV.Invert();
888   if (!tpV.IsValid()) return kFALSE;
889
890   TMatrixDSym pW(3),tW(3);
891   for (Int_t i=0; i<3; i++)
892     for (Int_t j=0; j<3; j++) {
893       pW(i,j)=tW(i,j)=0.;
894       for (Int_t k=0; k<3; k++) {
895         pW(i,j) += tV(i,k)*tpV(k,j);
896         tW(i,j) += pV(i,k)*tpV(k,j);
897       }
898     }
899
900   Double_t t[3] = {GetX(), GetY(), GetZ()};
901
902   Double_t x=0.;
903   for (Int_t i=0; i<3; i++) x += (tW(0,i)*t[i] + pW(0,i)*p[i]);  
904   Double_t crv=GetC(bz);
905   if (TMath::Abs(b) < kAlmost0Field) crv=0.;
906   f += crv*(x-fX);
907   if (TMath::Abs(f) >= kAlmost1) return kFALSE;
908   fX=x;  
909
910   fP[0]=0.;
911   for (Int_t i=0; i<3; i++) fP[0] += (tW(1,i)*t[i] + pW(1,i)*p[i]);  
912   fP[1]=0.;
913   for (Int_t i=0; i<3; i++) fP[1] += (tW(2,i)*t[i] + pW(2,i)*p[i]);  
914
915   return kTRUE;  
916 }
917
918 Double_t *AliExternalTrackParam::GetResiduals(
919 Double_t *p,Double_t *cov,Bool_t updated) const {
920   //------------------------------------------------------------------
921   // Returns the track residuals with the space point "p" having
922   // the covariance matrix "cov".
923   // If "updated" is kTRUE, the track parameters expected to be updated,
924   // otherwise they must be predicted.  
925   //------------------------------------------------------------------
926   static Double_t res[2];
927
928   Double_t r00=cov[0], r01=cov[1], r11=cov[2];
929   if (updated) {
930      r00-=fC[0]; r01-=fC[1]; r11-=fC[2];
931   } else {
932      r00+=fC[0]; r01+=fC[1]; r11+=fC[2];
933   }
934   Double_t det=r00*r11 - r01*r01;
935
936   if (TMath::Abs(det) < kAlmost0) return 0;
937
938   Double_t tmp=r00; r00=r11/det; r11=tmp/det;
939
940   if (r00 < 0.) return 0;
941   if (r11 < 0.) return 0;
942
943   Double_t dy = fP[0] - p[0];
944   Double_t dz = fP[1] - p[1];
945
946   res[0]=dy*TMath::Sqrt(r00);
947   res[1]=dz*TMath::Sqrt(r11);
948
949   return res;
950 }
951
952 Bool_t AliExternalTrackParam::Update(Double_t p[2], Double_t cov[3]) {
953   //------------------------------------------------------------------
954   // Update the track parameters with the space point "p" having
955   // the covariance matrix "cov"
956   //------------------------------------------------------------------
957   Double_t &fP0=fP[0], &fP1=fP[1], &fP2=fP[2], &fP3=fP[3], &fP4=fP[4];
958   Double_t 
959   &fC00=fC[0],
960   &fC10=fC[1],   &fC11=fC[2],  
961   &fC20=fC[3],   &fC21=fC[4],   &fC22=fC[5],
962   &fC30=fC[6],   &fC31=fC[7],   &fC32=fC[8],   &fC33=fC[9],  
963   &fC40=fC[10],  &fC41=fC[11],  &fC42=fC[12],  &fC43=fC[13], &fC44=fC[14];
964
965   Double_t r00=cov[0], r01=cov[1], r11=cov[2];
966   r00+=fC00; r01+=fC10; r11+=fC11;
967   Double_t det=r00*r11 - r01*r01;
968
969   if (TMath::Abs(det) < kAlmost0) return kFALSE;
970
971
972   Double_t tmp=r00; r00=r11/det; r11=tmp/det; r01=-r01/det;
973  
974   Double_t k00=fC00*r00+fC10*r01, k01=fC00*r01+fC10*r11;
975   Double_t k10=fC10*r00+fC11*r01, k11=fC10*r01+fC11*r11;
976   Double_t k20=fC20*r00+fC21*r01, k21=fC20*r01+fC21*r11;
977   Double_t k30=fC30*r00+fC31*r01, k31=fC30*r01+fC31*r11;
978   Double_t k40=fC40*r00+fC41*r01, k41=fC40*r01+fC41*r11;
979
980   Double_t dy=p[0] - fP0, dz=p[1] - fP1;
981   Double_t sf=fP2 + k20*dy + k21*dz;
982   if (TMath::Abs(sf) > kAlmost1) return kFALSE;  
983   
984   fP0 += k00*dy + k01*dz;
985   fP1 += k10*dy + k11*dz;
986   fP2  = sf;
987   fP3 += k30*dy + k31*dz;
988   fP4 += k40*dy + k41*dz;
989   
990   Double_t c01=fC10, c02=fC20, c03=fC30, c04=fC40;
991   Double_t c12=fC21, c13=fC31, c14=fC41;
992
993   fC00-=k00*fC00+k01*fC10; fC10-=k00*c01+k01*fC11;
994   fC20-=k00*c02+k01*c12;   fC30-=k00*c03+k01*c13;
995   fC40-=k00*c04+k01*c14; 
996
997   fC11-=k10*c01+k11*fC11;
998   fC21-=k10*c02+k11*c12;   fC31-=k10*c03+k11*c13;
999   fC41-=k10*c04+k11*c14; 
1000
1001   fC22-=k20*c02+k21*c12;   fC32-=k20*c03+k21*c13;
1002   fC42-=k20*c04+k21*c14; 
1003
1004   fC33-=k30*c03+k31*c13;
1005   fC43-=k30*c04+k31*c14; 
1006
1007   fC44-=k40*c04+k41*c14; 
1008
1009   return kTRUE;
1010 }
1011
1012 void 
1013 AliExternalTrackParam::GetHelixParameters(Double_t hlx[6], Double_t b) const {
1014   //--------------------------------------------------------------------
1015   // External track parameters -> helix parameters 
1016   // "b" - magnetic field (kG)
1017   //--------------------------------------------------------------------
1018   Double_t cs=TMath::Cos(fAlpha), sn=TMath::Sin(fAlpha);
1019   
1020   hlx[0]=fP[0]; hlx[1]=fP[1]; hlx[2]=fP[2]; hlx[3]=fP[3];
1021
1022   hlx[5]=fX*cs - hlx[0]*sn;               // x0
1023   hlx[0]=fX*sn + hlx[0]*cs;               // y0
1024 //hlx[1]=                                 // z0
1025   hlx[2]=TMath::ASin(hlx[2]) + fAlpha;    // phi0
1026 //hlx[3]=                                 // tgl
1027   hlx[4]=GetC(b);                         // C
1028 }
1029
1030
1031 static void Evaluate(const Double_t *h, Double_t t,
1032                      Double_t r[3],  //radius vector
1033                      Double_t g[3],  //first defivatives
1034                      Double_t gg[3]) //second derivatives
1035 {
1036   //--------------------------------------------------------------------
1037   // Calculate position of a point on a track and some derivatives
1038   //--------------------------------------------------------------------
1039   Double_t phase=h[4]*t+h[2];
1040   Double_t sn=TMath::Sin(phase), cs=TMath::Cos(phase);
1041
1042   r[0] = h[5] + (sn - h[6])/h[4];
1043   r[1] = h[0] - (cs - h[7])/h[4];  
1044   r[2] = h[1] + h[3]*t;
1045
1046   g[0] = cs; g[1]=sn; g[2]=h[3];
1047   
1048   gg[0]=-h[4]*sn; gg[1]=h[4]*cs; gg[2]=0.;
1049 }
1050
1051 Double_t AliExternalTrackParam::GetDCA(const AliExternalTrackParam *p, 
1052 Double_t b, Double_t &xthis, Double_t &xp) const {
1053   //------------------------------------------------------------
1054   // Returns the (weighed !) distance of closest approach between 
1055   // this track and the track "p".
1056   // Other returned values:
1057   //   xthis, xt - coordinates of tracks' reference planes at the DCA 
1058   //-----------------------------------------------------------
1059   Double_t dy2=GetSigmaY2() + p->GetSigmaY2();
1060   Double_t dz2=GetSigmaZ2() + p->GetSigmaZ2();
1061   Double_t dx2=dy2; 
1062
1063   Double_t p1[8]; GetHelixParameters(p1,b);
1064   p1[6]=TMath::Sin(p1[2]); p1[7]=TMath::Cos(p1[2]);
1065   Double_t p2[8]; p->GetHelixParameters(p2,b);
1066   p2[6]=TMath::Sin(p2[2]); p2[7]=TMath::Cos(p2[2]);
1067
1068
1069   Double_t r1[3],g1[3],gg1[3]; Double_t t1=0.;
1070   Evaluate(p1,t1,r1,g1,gg1);
1071   Double_t r2[3],g2[3],gg2[3]; Double_t t2=0.;
1072   Evaluate(p2,t2,r2,g2,gg2);
1073
1074   Double_t dx=r2[0]-r1[0], dy=r2[1]-r1[1], dz=r2[2]-r1[2];
1075   Double_t dm=dx*dx/dx2 + dy*dy/dy2 + dz*dz/dz2;
1076
1077   Int_t max=27;
1078   while (max--) {
1079      Double_t gt1=-(dx*g1[0]/dx2 + dy*g1[1]/dy2 + dz*g1[2]/dz2);
1080      Double_t gt2=+(dx*g2[0]/dx2 + dy*g2[1]/dy2 + dz*g2[2]/dz2);
1081      Double_t h11=(g1[0]*g1[0] - dx*gg1[0])/dx2 + 
1082                   (g1[1]*g1[1] - dy*gg1[1])/dy2 +
1083                   (g1[2]*g1[2] - dz*gg1[2])/dz2;
1084      Double_t h22=(g2[0]*g2[0] + dx*gg2[0])/dx2 + 
1085                   (g2[1]*g2[1] + dy*gg2[1])/dy2 +
1086                   (g2[2]*g2[2] + dz*gg2[2])/dz2;
1087      Double_t h12=-(g1[0]*g2[0]/dx2 + g1[1]*g2[1]/dy2 + g1[2]*g2[2]/dz2);
1088
1089      Double_t det=h11*h22-h12*h12;
1090
1091      Double_t dt1,dt2;
1092      if (TMath::Abs(det)<1.e-33) {
1093         //(quasi)singular Hessian
1094         dt1=-gt1; dt2=-gt2;
1095      } else {
1096         dt1=-(gt1*h22 - gt2*h12)/det; 
1097         dt2=-(h11*gt2 - h12*gt1)/det;
1098      }
1099
1100      if ((dt1*gt1+dt2*gt2)>0) {dt1=-dt1; dt2=-dt2;}
1101
1102      //check delta(phase1) ?
1103      //check delta(phase2) ?
1104
1105      if (TMath::Abs(dt1)/(TMath::Abs(t1)+1.e-3) < 1.e-4)
1106      if (TMath::Abs(dt2)/(TMath::Abs(t2)+1.e-3) < 1.e-4) {
1107         if ((gt1*gt1+gt2*gt2) > 1.e-4/dy2/dy2) 
1108           AliDebug(1," stopped at not a stationary point !");
1109         Double_t lmb=h11+h22; lmb=lmb-TMath::Sqrt(lmb*lmb-4*det);
1110         if (lmb < 0.) 
1111           AliDebug(1," stopped at not a minimum !");
1112         break;
1113      }
1114
1115      Double_t dd=dm;
1116      for (Int_t div=1 ; ; div*=2) {
1117         Evaluate(p1,t1+dt1,r1,g1,gg1);
1118         Evaluate(p2,t2+dt2,r2,g2,gg2);
1119         dx=r2[0]-r1[0]; dy=r2[1]-r1[1]; dz=r2[2]-r1[2];
1120         dd=dx*dx/dx2 + dy*dy/dy2 + dz*dz/dz2;
1121         if (dd<dm) break;
1122         dt1*=0.5; dt2*=0.5;
1123         if (div>512) {
1124           AliDebug(1," overshoot !"); break;
1125         }   
1126      }
1127      dm=dd;
1128
1129      t1+=dt1;
1130      t2+=dt2;
1131
1132   }
1133
1134   if (max<=0) AliDebug(1," too many iterations !");
1135
1136   Double_t cs=TMath::Cos(GetAlpha());
1137   Double_t sn=TMath::Sin(GetAlpha());
1138   xthis=r1[0]*cs + r1[1]*sn;
1139
1140   cs=TMath::Cos(p->GetAlpha());
1141   sn=TMath::Sin(p->GetAlpha());
1142   xp=r2[0]*cs + r2[1]*sn;
1143
1144   return TMath::Sqrt(dm*TMath::Sqrt(dy2*dz2));
1145 }
1146  
1147 Double_t AliExternalTrackParam::
1148 PropagateToDCA(AliExternalTrackParam *p, Double_t b) {
1149   //--------------------------------------------------------------
1150   // Propagates this track and the argument track to the position of the
1151   // distance of closest approach.
1152   // Returns the (weighed !) distance of closest approach.
1153   //--------------------------------------------------------------
1154   Double_t xthis,xp;
1155   Double_t dca=GetDCA(p,b,xthis,xp);
1156
1157   if (!PropagateTo(xthis,b)) {
1158     //AliWarning(" propagation failed !");
1159     return 1e+33;
1160   }
1161
1162   if (!p->PropagateTo(xp,b)) {
1163     //AliWarning(" propagation failed !";
1164     return 1e+33;
1165   }
1166
1167   return dca;
1168 }
1169
1170
1171 Bool_t AliExternalTrackParam::PropagateToDCA(const AliVVertex *vtx, 
1172 Double_t b, Double_t maxd, Double_t dz[2], Double_t covar[3]) {
1173   //
1174   // Propagate this track to the DCA to vertex "vtx", 
1175   // if the (rough) transverse impact parameter is not bigger then "maxd". 
1176   //            Magnetic field is "b" (kG).
1177   //
1178   // a) The track gets extapolated to the DCA to the vertex.
1179   // b) The impact parameters and their covariance matrix are calculated.
1180   //
1181   //    In the case of success, the returned value is kTRUE
1182   //    (otherwise, it's kFALSE)
1183   //  
1184   Double_t alpha=GetAlpha();
1185   Double_t sn=TMath::Sin(alpha), cs=TMath::Cos(alpha);
1186   Double_t x=GetX(), y=GetParameter()[0], snp=GetParameter()[2];
1187   Double_t xv= vtx->GetX()*cs + vtx->GetY()*sn;
1188   Double_t yv=-vtx->GetX()*sn + vtx->GetY()*cs, zv=vtx->GetZ();
1189   x-=xv; y-=yv;
1190
1191   //Estimate the impact parameter neglecting the track curvature
1192   Double_t d=TMath::Abs(x*snp - y*TMath::Sqrt(1.- snp*snp));
1193   if (d > maxd) return kFALSE; 
1194
1195   //Propagate to the DCA
1196   Double_t crv=GetC(b);
1197   if (TMath::Abs(b) < kAlmost0Field) crv=0.;
1198
1199   Double_t tgfv=-(crv*x - snp)/(crv*y + TMath::Sqrt(1.-snp*snp));
1200   sn=tgfv/TMath::Sqrt(1.+ tgfv*tgfv); cs=TMath::Sqrt(1.- sn*sn);
1201   if (TMath::Abs(tgfv)>0.) cs = sn/tgfv;
1202   else cs=1.;
1203
1204   x = xv*cs + yv*sn;
1205   yv=-xv*sn + yv*cs; xv=x;
1206
1207   if (!Propagate(alpha+TMath::ASin(sn),xv,b)) return kFALSE;
1208
1209   if (dz==0) return kTRUE;
1210   dz[0] = GetParameter()[0] - yv;
1211   dz[1] = GetParameter()[1] - zv;
1212   
1213   if (covar==0) return kTRUE;
1214   Double_t cov[6]; vtx->GetCovarianceMatrix(cov);
1215
1216   //***** Improvements by A.Dainese
1217   alpha=GetAlpha(); sn=TMath::Sin(alpha); cs=TMath::Cos(alpha);
1218   Double_t s2ylocvtx = cov[0]*sn*sn + cov[2]*cs*cs - 2.*cov[1]*cs*sn;
1219   covar[0] = GetCovariance()[0] + s2ylocvtx;   // neglecting correlations
1220   covar[1] = GetCovariance()[1];               // between (x,y) and z
1221   covar[2] = GetCovariance()[2] + cov[5];      // in vertex's covariance matrix
1222   //*****
1223
1224   return kTRUE;
1225 }
1226
1227
1228 void AliExternalTrackParam::GetDirection(Double_t d[3]) const {
1229   //----------------------------------------------------------------
1230   // This function returns a unit vector along the track direction
1231   // in the global coordinate system.
1232   //----------------------------------------------------------------
1233   Double_t cs=TMath::Cos(fAlpha), sn=TMath::Sin(fAlpha);
1234   Double_t snp=fP[2];
1235   Double_t csp =TMath::Sqrt((1.- snp)*(1.+snp));
1236   Double_t norm=TMath::Sqrt(1.+ fP[3]*fP[3]);
1237   d[0]=(csp*cs - snp*sn)/norm; 
1238   d[1]=(snp*cs + csp*sn)/norm; 
1239   d[2]=fP[3]/norm;
1240 }
1241
1242 Bool_t AliExternalTrackParam::GetPxPyPz(Double_t p[3]) const {
1243   //---------------------------------------------------------------------
1244   // This function returns the global track momentum components
1245   // Results for (nearly) straight tracks are meaningless !
1246   //---------------------------------------------------------------------
1247   p[0]=fP[4]; p[1]=fP[2]; p[2]=fP[3];
1248   return Local2GlobalMomentum(p,fAlpha);
1249 }
1250
1251 Double_t AliExternalTrackParam::Px() const {
1252   //---------------------------------------------------------------------
1253   // Returns x-component of momentum
1254   // Result for (nearly) straight tracks is meaningless !
1255   //---------------------------------------------------------------------
1256
1257   Double_t p[3]={kVeryBig,kVeryBig,kVeryBig};
1258   GetPxPyPz(p);
1259
1260   return p[0];
1261 }
1262
1263 Double_t AliExternalTrackParam::Py() const {
1264   //---------------------------------------------------------------------
1265   // Returns y-component of momentum
1266   // Result for (nearly) straight tracks is meaningless !
1267   //---------------------------------------------------------------------
1268
1269   Double_t p[3]={kVeryBig,kVeryBig,kVeryBig};
1270   GetPxPyPz(p);
1271
1272   return p[1];
1273 }
1274
1275 Double_t AliExternalTrackParam::Pz() const {
1276   //---------------------------------------------------------------------
1277   // Returns z-component of momentum
1278   // Result for (nearly) straight tracks is meaningless !
1279   //---------------------------------------------------------------------
1280
1281   Double_t p[3]={kVeryBig,kVeryBig,kVeryBig};
1282   GetPxPyPz(p);
1283
1284   return p[2];
1285 }
1286
1287 Double_t AliExternalTrackParam::Xv() const {
1288   //---------------------------------------------------------------------
1289   // Returns x-component of first track point
1290   //---------------------------------------------------------------------
1291
1292   Double_t r[3]={0.,0.,0.};
1293   GetXYZ(r);
1294
1295   return r[0];
1296 }
1297
1298 Double_t AliExternalTrackParam::Yv() const {
1299   //---------------------------------------------------------------------
1300   // Returns y-component of first track point
1301   //---------------------------------------------------------------------
1302
1303   Double_t r[3]={0.,0.,0.};
1304   GetXYZ(r);
1305
1306   return r[1];
1307 }
1308
1309 Double_t AliExternalTrackParam::Zv() const {
1310   //---------------------------------------------------------------------
1311   // Returns z-component of first track point
1312   //---------------------------------------------------------------------
1313
1314   Double_t r[3]={0.,0.,0.};
1315   GetXYZ(r);
1316
1317   return r[2];
1318 }
1319
1320 Double_t AliExternalTrackParam::Theta() const {
1321   // return theta angle of momentum
1322
1323   return 0.5*TMath::Pi() - TMath::ATan(fP[3]);
1324 }
1325
1326 Double_t AliExternalTrackParam::Phi() const {
1327   //---------------------------------------------------------------------
1328   // Returns the azimuthal angle of momentum
1329   // 0 <= phi < 2*pi
1330   //---------------------------------------------------------------------
1331
1332   Double_t phi=TMath::ASin(fP[2]) + fAlpha;
1333   if (phi<0.) phi+=2.*TMath::Pi();
1334   else if (phi>=2.*TMath::Pi()) phi-=2.*TMath::Pi();
1335  
1336   return phi;
1337 }
1338
1339 Double_t AliExternalTrackParam::M() const {
1340   // return particle mass
1341
1342   // No mass information available so far.
1343   // Redifine in derived class!
1344
1345   return -999.;
1346 }
1347
1348 Double_t AliExternalTrackParam::E() const {
1349   // return particle energy
1350
1351   // No PID information available so far.
1352   // Redifine in derived class!
1353
1354   return -999.;
1355 }
1356
1357 Double_t AliExternalTrackParam::Eta() const { 
1358   // return pseudorapidity
1359
1360   return -TMath::Log(TMath::Tan(0.5 * Theta())); 
1361 }
1362
1363 Double_t AliExternalTrackParam::Y() const {
1364   // return rapidity
1365
1366   // No PID information available so far.
1367   // Redifine in derived class!
1368
1369   return -999.;
1370 }
1371
1372 Bool_t AliExternalTrackParam::GetXYZ(Double_t *r) const {
1373   //---------------------------------------------------------------------
1374   // This function returns the global track position
1375   //---------------------------------------------------------------------
1376   r[0]=fX; r[1]=fP[0]; r[2]=fP[1];
1377   return Local2GlobalPosition(r,fAlpha);
1378 }
1379
1380 Bool_t AliExternalTrackParam::GetCovarianceXYZPxPyPz(Double_t cv[21]) const {
1381   //---------------------------------------------------------------------
1382   // This function returns the global covariance matrix of the track params
1383   // 
1384   // Cov(x,x) ... :   cv[0]
1385   // Cov(y,x) ... :   cv[1]  cv[2]
1386   // Cov(z,x) ... :   cv[3]  cv[4]  cv[5]
1387   // Cov(px,x)... :   cv[6]  cv[7]  cv[8]  cv[9]
1388   // Cov(py,x)... :   cv[10] cv[11] cv[12] cv[13] cv[14]
1389   // Cov(pz,x)... :   cv[15] cv[16] cv[17] cv[18] cv[19] cv[20]
1390   //
1391   // Results for (nearly) straight tracks are meaningless !
1392   //---------------------------------------------------------------------
1393   if (TMath::Abs(fP[4])<=kAlmost0) {
1394      for (Int_t i=0; i<21; i++) cv[i]=0.;
1395      return kFALSE;
1396   }
1397   if (TMath::Abs(fP[2]) > kAlmost1) {
1398      for (Int_t i=0; i<21; i++) cv[i]=0.;
1399      return kFALSE;
1400   }
1401   Double_t pt=1./TMath::Abs(fP[4]);
1402   Double_t cs=TMath::Cos(fAlpha), sn=TMath::Sin(fAlpha);
1403   Double_t r=TMath::Sqrt((1.-fP[2])*(1.+fP[2]));
1404
1405   Double_t m00=-sn, m10=cs;
1406   Double_t m23=-pt*(sn + fP[2]*cs/r), m43=-pt*pt*(r*cs - fP[2]*sn);
1407   Double_t m24= pt*(cs - fP[2]*sn/r), m44=-pt*pt*(r*sn + fP[2]*cs);
1408   Double_t m35=pt, m45=-pt*pt*fP[3];
1409
1410   m43*=GetSign();
1411   m44*=GetSign();
1412   m45*=GetSign();
1413
1414   cv[0 ] = fC[0]*m00*m00;
1415   cv[1 ] = fC[0]*m00*m10; 
1416   cv[2 ] = fC[0]*m10*m10;
1417   cv[3 ] = fC[1]*m00; 
1418   cv[4 ] = fC[1]*m10; 
1419   cv[5 ] = fC[2];
1420   cv[6 ] = m00*(fC[3]*m23 + fC[10]*m43); 
1421   cv[7 ] = m10*(fC[3]*m23 + fC[10]*m43); 
1422   cv[8 ] = fC[4]*m23 + fC[11]*m43; 
1423   cv[9 ] = m23*(fC[5]*m23 + fC[12]*m43)  +  m43*(fC[12]*m23 + fC[14]*m43);
1424   cv[10] = m00*(fC[3]*m24 + fC[10]*m44); 
1425   cv[11] = m10*(fC[3]*m24 + fC[10]*m44); 
1426   cv[12] = fC[4]*m24 + fC[11]*m44; 
1427   cv[13] = m23*(fC[5]*m24 + fC[12]*m44)  +  m43*(fC[12]*m24 + fC[14]*m44);
1428   cv[14] = m24*(fC[5]*m24 + fC[12]*m44)  +  m44*(fC[12]*m24 + fC[14]*m44);
1429   cv[15] = m00*(fC[6]*m35 + fC[10]*m45); 
1430   cv[16] = m10*(fC[6]*m35 + fC[10]*m45); 
1431   cv[17] = fC[7]*m35 + fC[11]*m45; 
1432   cv[18] = m23*(fC[8]*m35 + fC[12]*m45)  +  m43*(fC[13]*m35 + fC[14]*m45);
1433   cv[19] = m24*(fC[8]*m35 + fC[12]*m45)  +  m44*(fC[13]*m35 + fC[14]*m45); 
1434   cv[20] = m35*(fC[9]*m35 + fC[13]*m45)  +  m45*(fC[13]*m35 + fC[14]*m45);
1435
1436   return kTRUE;
1437 }
1438
1439
1440 Bool_t 
1441 AliExternalTrackParam::GetPxPyPzAt(Double_t x, Double_t b, Double_t *p) const {
1442   //---------------------------------------------------------------------
1443   // This function returns the global track momentum extrapolated to
1444   // the radial position "x" (cm) in the magnetic field "b" (kG)
1445   //---------------------------------------------------------------------
1446   p[0]=fP[4]; 
1447   p[1]=fP[2]+(x-fX)*GetC(b); 
1448   p[2]=fP[3];
1449   return Local2GlobalMomentum(p,fAlpha);
1450 }
1451
1452 Bool_t 
1453 AliExternalTrackParam::GetYAt(Double_t x, Double_t b, Double_t &y) const {
1454   //---------------------------------------------------------------------
1455   // This function returns the local Y-coordinate of the intersection 
1456   // point between this track and the reference plane "x" (cm). 
1457   // Magnetic field "b" (kG)
1458   //---------------------------------------------------------------------
1459   Double_t dx=x-fX;
1460   if(TMath::Abs(dx)<=kAlmost0) {y=fP[0]; return kTRUE;}
1461
1462   Double_t f1=fP[2], f2=f1 + dx*GetC(b);
1463
1464   if (TMath::Abs(f1) >= kAlmost1) return kFALSE;
1465   if (TMath::Abs(f2) >= kAlmost1) return kFALSE;
1466   
1467   Double_t r1=TMath::Sqrt(1.- f1*f1), r2=TMath::Sqrt(1.- f2*f2);
1468   y = fP[0] + dx*(f1+f2)/(r1+r2);
1469   return kTRUE;
1470 }
1471
1472 Bool_t 
1473 AliExternalTrackParam::GetZAt(Double_t x, Double_t b, Double_t &z) const {
1474   //---------------------------------------------------------------------
1475   // This function returns the local Z-coordinate of the intersection 
1476   // point between this track and the reference plane "x" (cm). 
1477   // Magnetic field "b" (kG)
1478   //---------------------------------------------------------------------
1479   Double_t dx=x-fX;
1480   if(TMath::Abs(dx)<=kAlmost0) {z=fP[1]; return kTRUE;}
1481
1482   Double_t f1=fP[2], f2=f1 + dx*GetC(b);
1483
1484   if (TMath::Abs(f1) >= kAlmost1) return kFALSE;
1485   if (TMath::Abs(f2) >= kAlmost1) return kFALSE;
1486   
1487   Double_t r1=sqrt(1.- f1*f1), r2=sqrt(1.- f2*f2);
1488   z = fP[1] + dx*(r2 + f2*(f1+f2)/(r1+r2))*fP[3]; // Many thanks to P.Hristov !
1489   return kTRUE;
1490 }
1491
1492 Bool_t 
1493 AliExternalTrackParam::GetXYZAt(Double_t x, Double_t b, Double_t *r) const {
1494   //---------------------------------------------------------------------
1495   // This function returns the global track position extrapolated to
1496   // the radial position "x" (cm) in the magnetic field "b" (kG)
1497   //---------------------------------------------------------------------
1498   Double_t dx=x-fX;
1499   if(TMath::Abs(dx)<=kAlmost0) return GetXYZ(r);
1500
1501   Double_t f1=fP[2], f2=f1 + dx*GetC(b);
1502
1503   if (TMath::Abs(f1) >= kAlmost1) return kFALSE;
1504   if (TMath::Abs(f2) >= kAlmost1) return kFALSE;
1505   
1506   Double_t r1=TMath::Sqrt(1.- f1*f1), r2=TMath::Sqrt(1.- f2*f2);
1507   r[0] = x;
1508   r[1] = fP[0] + dx*(f1+f2)/(r1+r2);
1509   r[2] = fP[1] + dx*(r2 + f2*(f1+f2)/(r1+r2))*fP[3];//Thanks to Andrea & Peter
1510
1511   return Local2GlobalPosition(r,fAlpha);
1512 }
1513
1514 //_____________________________________________________________________________
1515 void AliExternalTrackParam::Print(Option_t* /*option*/) const
1516 {
1517 // print the parameters and the covariance matrix
1518
1519   printf("AliExternalTrackParam: x = %-12g  alpha = %-12g\n", fX, fAlpha);
1520   printf("  parameters: %12g %12g %12g %12g %12g\n",
1521          fP[0], fP[1], fP[2], fP[3], fP[4]);
1522   printf("  covariance: %12g\n", fC[0]);
1523   printf("              %12g %12g\n", fC[1], fC[2]);
1524   printf("              %12g %12g %12g\n", fC[3], fC[4], fC[5]);
1525   printf("              %12g %12g %12g %12g\n", 
1526          fC[6], fC[7], fC[8], fC[9]);
1527   printf("              %12g %12g %12g %12g %12g\n", 
1528          fC[10], fC[11], fC[12], fC[13], fC[14]);
1529 }
1530
1531 Double_t AliExternalTrackParam::GetSnpAt(Double_t x,Double_t b) const {
1532   //
1533   // Get sinus at given x
1534   //
1535   Double_t crv=GetC(b);
1536   if (TMath::Abs(b) < kAlmost0Field) crv=0.;
1537   Double_t dx = x-fX;
1538   Double_t res = fP[2]+dx*crv;
1539   return res;
1540 }
1541
1542 Bool_t AliExternalTrackParam::GetDistance(AliExternalTrackParam *param2, Double_t x, Double_t dist[3], Double_t bz){
1543   //------------------------------------------------------------------------
1544   // Get the distance between two tracks at the local position x 
1545   // working in the local frame of this track.
1546   // Origin :   Marian.Ivanov@cern.ch
1547   //-----------------------------------------------------------------------
1548   Double_t xyz[3];
1549   Double_t xyz2[3];
1550   xyz[0]=x;
1551   if (!GetYAt(x,bz,xyz[1])) return kFALSE;
1552   if (!GetZAt(x,bz,xyz[2])) return kFALSE;
1553   //  
1554   //
1555   if (TMath::Abs(GetAlpha()-param2->GetAlpha())<kAlmost0){
1556     xyz2[0]=x;
1557     if (!param2->GetYAt(x,bz,xyz2[1])) return kFALSE;
1558     if (!param2->GetZAt(x,bz,xyz2[2])) return kFALSE;
1559   }else{
1560     //
1561     Double_t xyz1[3];
1562     Double_t dfi = param2->GetAlpha()-GetAlpha();
1563     Double_t ca = TMath::Cos(dfi), sa = TMath::Sin(dfi);
1564     xyz2[0] =  xyz[0]*ca+xyz[1]*sa;
1565     xyz2[1] = -xyz[0]*sa+xyz[1]*ca;
1566     //
1567     xyz1[0]=xyz2[0];
1568     if (!param2->GetYAt(xyz2[0],bz,xyz1[1])) return kFALSE;
1569     if (!param2->GetZAt(xyz2[0],bz,xyz1[2])) return kFALSE;
1570     //
1571     xyz2[0] =  xyz1[0]*ca-xyz1[1]*sa;
1572     xyz2[1] = +xyz1[0]*sa+xyz1[1]*ca;
1573     xyz2[2] = xyz1[2];
1574   }
1575   dist[0] = xyz[0]-xyz2[0];
1576   dist[1] = xyz[1]-xyz2[1];
1577   dist[2] = xyz[2]-xyz2[2];
1578
1579   return kTRUE;
1580 }
1581
1582
1583 //
1584 // Draw functionality.
1585 // Origin: Marian Ivanov, Marian.Ivanov@cern.ch
1586 //
1587
1588 void  AliExternalTrackParam::DrawTrack(Float_t magf, Float_t minR, Float_t maxR, Float_t stepR){
1589   //
1590   // Draw track line
1591   //
1592   if (minR>maxR) return ;
1593   if (stepR<=0) return ;
1594   Int_t npoints = TMath::Nint((maxR-minR)/stepR)+1;
1595   if (npoints<1) return;
1596   TPolyMarker3D *polymarker = new TPolyMarker3D(npoints);
1597   FillPolymarker(polymarker, magf,minR,maxR,stepR);
1598   polymarker->Draw();
1599 }
1600
1601 //
1602 void AliExternalTrackParam::FillPolymarker(TPolyMarker3D *pol, Float_t magF, Float_t minR, Float_t maxR, Float_t stepR){
1603   //
1604   // Fill points in the polymarker
1605   //
1606   Int_t counter=0;
1607   for (Double_t r=minR; r<maxR; r+=stepR){
1608     Double_t point[3];
1609     GetXYZAt(r,magF,point);
1610     pol->SetPoint(counter,point[0],point[1], point[2]);
1611     printf("xyz\t%f\t%f\t%f\n",point[0], point[1],point[2]);
1612     counter++;
1613   }
1614 }
1615
1616 Int_t AliExternalTrackParam::GetIndex(Int_t i, Int_t j) const {
1617   //
1618   Int_t min = TMath::Min(i,j);
1619   Int_t max = TMath::Max(i,j);
1620
1621   return min+(max+1)*max/2;
1622 }