]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONTrackExtrap.cxx
e7762acbff43e0f3e0e75dd3297a308ecaf9475b
[u/mrichter/AliRoot.git] / MUON / AliMUONTrackExtrap.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 // Class AliMUONTrackExtrap
20 // ------------------------
21 // Tools for track extrapolation in ALICE dimuon spectrometer
22 // Author: Philippe Pillot
23 //-----------------------------------------------------------------------------
24
25 #include "AliMUONTrackExtrap.h" 
26 #include "AliMUONTrackParam.h"
27 #include "AliMUONConstants.h"
28 #include "AliMUONReconstructor.h"
29 #include "AliMUONRecoParam.h"
30
31 #include "AliMagF.h" 
32
33 #include <TMath.h>
34 #include <TGeoManager.h>
35
36 #include <Riostream.h>
37
38 /// \cond CLASSIMP
39 ClassImp(AliMUONTrackExtrap) // Class implementation in ROOT context
40 /// \endcond
41
42 const AliMagF* AliMUONTrackExtrap::fgkField = 0x0;
43 const Double_t AliMUONTrackExtrap::fgkSimpleBPosition = 0.5 * (AliMUONConstants::CoilZ() + AliMUONConstants::YokeZ());
44 const Double_t AliMUONTrackExtrap::fgkSimpleBLength = 0.5 * (AliMUONConstants::CoilL() + AliMUONConstants::YokeL());
45       Double_t AliMUONTrackExtrap::fgSimpleBValue = 0.;
46       Bool_t   AliMUONTrackExtrap::fgFieldON = kFALSE;
47 const Bool_t   AliMUONTrackExtrap::fgkUseHelix = kFALSE;
48 const Int_t    AliMUONTrackExtrap::fgkMaxStepNumber = 5000;
49 const Double_t AliMUONTrackExtrap::fgkHelixStepLength = 6.;
50 const Double_t AliMUONTrackExtrap::fgkRungeKuttaMaxResidue = 0.002;
51
52 //__________________________________________________________________________
53 void AliMUONTrackExtrap::SetField(const AliMagF* magField)
54 {
55   /// set magnetic field
56   
57   // set field map
58   fgkField = magField;
59   if (!fgkField) {
60     cout<<"E-AliMUONTrackExtrap::SetField: fgkField = 0x0"<<endl;
61     return;
62   }
63   
64   // set field on/off flag
65   fgFieldON = (fgkField->Factor() == 0.) ? kFALSE : kTRUE;
66   
67   // set field at the centre of the dipole
68   if (fgFieldON) {
69     Float_t b[3] = {0.,0.,0.}, x[3] = {50.,50.,(Float_t) fgkSimpleBPosition};
70     fgkField->Field(x,b);
71     fgSimpleBValue = (Double_t) b[0];
72   } else fgSimpleBValue = 0.;
73   
74 }
75
76 //__________________________________________________________________________
77 Double_t AliMUONTrackExtrap::GetImpactParamFromBendingMomentum(Double_t bendingMomentum)
78 {
79   /// Returns impact parameter at vertex in bending plane (cm),
80   /// from the signed bending momentum "BendingMomentum" in bending plane (GeV/c),
81   /// using simple values for dipole magnetic field.
82   /// The sign of "BendingMomentum" is the sign of the charge.
83   
84   if (bendingMomentum == 0.) return 1.e10;
85   
86   if (!fgkField) {
87     cout<<"F-AliMUONTrackExtrap::GetField: fgkField = 0x0"<<endl;
88     exit(-1);
89   }
90   
91   const Double_t kCorrectionFactor = 0.9; // impact parameter is 10% overestimated
92   
93   return kCorrectionFactor * (-0.0003 * fgSimpleBValue * fgkSimpleBLength * fgkSimpleBPosition / bendingMomentum);
94 }
95
96 //__________________________________________________________________________
97 Double_t 
98 AliMUONTrackExtrap::GetBendingMomentumFromImpactParam(Double_t impactParam)
99 {
100   /// Returns signed bending momentum in bending plane (GeV/c),
101   /// the sign being the sign of the charge for particles moving forward in Z,
102   /// from the impact parameter "ImpactParam" at vertex in bending plane (cm),
103   /// using simple values for dipole magnetic field.
104   
105   if (impactParam == 0.) return 1.e10;
106   
107   if (!fgkField) {
108     cout<<"F-AliMUONTrackExtrap::GetField: fgkField = 0x0"<<endl;
109     exit(-1);
110   }
111   
112   const Double_t kCorrectionFactor = 1.1; // bending momentum is 10% underestimated
113   
114   if (fgFieldON) 
115   {
116     return kCorrectionFactor * (-0.0003 * fgSimpleBValue * fgkSimpleBLength * fgkSimpleBPosition / impactParam);
117   }
118   else 
119   {
120     return AliMUONConstants::GetMostProbBendingMomentum();
121   }
122 }
123
124 //__________________________________________________________________________
125 void AliMUONTrackExtrap::LinearExtrapToZ(AliMUONTrackParam* trackParam, Double_t zEnd, Bool_t updatePropagator)
126 {
127   /// Track parameters (and their covariances if any) linearly extrapolated to the plane at "zEnd".
128   /// On return, results from the extrapolation are updated in trackParam.
129   
130   if (trackParam->GetZ() == zEnd) return; // nothing to be done if same z
131   
132   // Compute track parameters
133   Double_t dZ = zEnd - trackParam->GetZ();
134   trackParam->SetNonBendingCoor(trackParam->GetNonBendingCoor() + trackParam->GetNonBendingSlope() * dZ);
135   trackParam->SetBendingCoor(trackParam->GetBendingCoor() + trackParam->GetBendingSlope() * dZ);
136   trackParam->SetZ(zEnd);
137   
138   // Update track parameters covariances if any
139   if (trackParam->CovariancesExist()) {
140     TMatrixD paramCov(trackParam->GetCovariances());
141     paramCov(0,0) += dZ * dZ * paramCov(1,1) + 2. * dZ * paramCov(0,1);
142     paramCov(0,1) += dZ * paramCov(1,1);
143     paramCov(1,0) = paramCov(0,1);
144     paramCov(2,2) += dZ * dZ * paramCov(3,3) + 2. * dZ * paramCov(2,3);
145     paramCov(2,3) += dZ * paramCov(3,3);
146     paramCov(3,2) = paramCov(2,3);
147     trackParam->SetCovariances(paramCov);
148     
149     // Update the propagator if required
150     if (updatePropagator) {
151       TMatrixD jacob(5,5);
152       jacob.UnitMatrix();
153       jacob(0,1) = dZ;
154       jacob(2,3) = dZ;
155       trackParam->UpdatePropagator(jacob);
156     }
157     
158   }
159   
160 }
161
162 //__________________________________________________________________________
163 void AliMUONTrackExtrap::ExtrapToZ(AliMUONTrackParam* trackParam, Double_t zEnd)
164 {
165   /// Interface to track parameter extrapolation to the plane at "Z" using Helix or Rungekutta algorithm.
166   /// On return, the track parameters resulting from the extrapolation are updated in trackParam.
167   if (!fgFieldON) AliMUONTrackExtrap::LinearExtrapToZ(trackParam,zEnd);
168   else if (fgkUseHelix) AliMUONTrackExtrap::ExtrapToZHelix(trackParam,zEnd);
169   else AliMUONTrackExtrap::ExtrapToZRungekutta(trackParam,zEnd);
170 }
171
172 //__________________________________________________________________________
173 void AliMUONTrackExtrap::ExtrapToZHelix(AliMUONTrackParam* trackParam, Double_t zEnd)
174 {
175   /// Track parameter extrapolation to the plane at "Z" using Helix algorithm.
176   /// On return, the track parameters resulting from the extrapolation are updated in trackParam.
177   if (trackParam->GetZ() == zEnd) return; // nothing to be done if same Z
178   Double_t forwardBackward; // +1 if forward, -1 if backward
179   if (zEnd < trackParam->GetZ()) forwardBackward = 1.0; // spectro. z<0 
180   else forwardBackward = -1.0;
181   Double_t v3[7], v3New[7]; // 7 in parameter ????
182   Int_t i3, stepNumber;
183   // For safety: return kTRUE or kFALSE ????
184   // Parameter vector for calling EXTRAP_ONESTEP
185   ConvertTrackParamForExtrap(trackParam, forwardBackward, v3);
186   // sign of charge (sign of fInverseBendingMomentum if forward motion)
187   // must be changed if backward extrapolation
188   Double_t chargeExtrap = forwardBackward * TMath::Sign(Double_t(1.0), trackParam->GetInverseBendingMomentum());
189   // Extrapolation loop
190   stepNumber = 0;
191   while (((-forwardBackward * (v3[2] - zEnd)) <= 0.0) && (stepNumber < fgkMaxStepNumber)) { // spectro. z<0
192     stepNumber++;
193     ExtrapOneStepHelix(chargeExtrap, fgkHelixStepLength, v3, v3New);
194     if ((-forwardBackward * (v3New[2] - zEnd)) > 0.0) break; // one is beyond Z spectro. z<0
195                                                              // better use TArray ????
196     for (i3 = 0; i3 < 7; i3++) {v3[i3] = v3New[i3];}
197   }
198   // check fgkMaxStepNumber ????
199   // Interpolation back to exact Z (2nd order)
200   // should be in function ???? using TArray ????
201   Double_t dZ12 = v3New[2] - v3[2]; // 1->2
202   if (TMath::Abs(dZ12) > 0) {
203     Double_t dZ1i = zEnd - v3[2]; // 1-i
204     Double_t dZi2 = v3New[2] - zEnd; // i->2
205     Double_t xPrime = (v3New[0] - v3[0]) / dZ12;
206     Double_t xSecond = ((v3New[3] / v3New[5]) - (v3[3] / v3[5])) / dZ12;
207     Double_t yPrime = (v3New[1] - v3[1]) / dZ12;
208     Double_t ySecond = ((v3New[4] / v3New[5]) - (v3[4] / v3[5])) / dZ12;
209     v3[0] = v3[0] + xPrime * dZ1i - 0.5 * xSecond * dZ1i * dZi2; // X
210     v3[1] = v3[1] + yPrime * dZ1i - 0.5 * ySecond * dZ1i * dZi2; // Y
211     v3[2] = zEnd; // Z
212     Double_t xPrimeI = xPrime - 0.5 * xSecond * (dZi2 - dZ1i);
213     Double_t yPrimeI = yPrime - 0.5 * ySecond * (dZi2 - dZ1i);
214     // (PX, PY, PZ)/PTOT assuming forward motion
215     v3[5] = 1.0 / TMath::Sqrt(1.0 + xPrimeI * xPrimeI + yPrimeI * yPrimeI); // PZ/PTOT
216     v3[3] = xPrimeI * v3[5]; // PX/PTOT
217     v3[4] = yPrimeI * v3[5]; // PY/PTOT
218   } else {
219     cout<<"W-AliMUONTrackExtrap::ExtrapToZHelix: Extrap. to Z not reached, Z = "<<zEnd<<endl;
220   }
221   // Recover track parameters (charge back for forward motion)
222   RecoverTrackParam(v3, chargeExtrap * forwardBackward, trackParam);
223 }
224
225 //__________________________________________________________________________
226 void AliMUONTrackExtrap::ExtrapToZRungekutta(AliMUONTrackParam* trackParam, Double_t zEnd)
227 {
228   /// Track parameter extrapolation to the plane at "Z" using Rungekutta algorithm.
229   /// On return, the track parameters resulting from the extrapolation are updated in trackParam.
230   if (trackParam->GetZ() == zEnd) return; // nothing to be done if same Z
231   Double_t forwardBackward; // +1 if forward, -1 if backward
232   if (zEnd < trackParam->GetZ()) forwardBackward = 1.0; // spectro. z<0 
233   else forwardBackward = -1.0;
234   // sign of charge (sign of fInverseBendingMomentum if forward motion)
235   // must be changed if backward extrapolation
236   Double_t chargeExtrap = forwardBackward * TMath::Sign(Double_t(1.0), trackParam->GetInverseBendingMomentum());
237   Double_t v3[7], v3New[7];
238   Double_t dZ, step;
239   Int_t stepNumber = 0;
240   
241   // Extrapolation loop (until within tolerance)
242   Double_t residue = zEnd - trackParam->GetZ();
243   while (TMath::Abs(residue) > fgkRungeKuttaMaxResidue && stepNumber <= fgkMaxStepNumber) {
244     dZ = zEnd - trackParam->GetZ();
245     // step lenght assuming linear trajectory
246     step = dZ * TMath::Sqrt(1.0 + trackParam->GetBendingSlope()*trackParam->GetBendingSlope() +
247                             trackParam->GetNonBendingSlope()*trackParam->GetNonBendingSlope());
248     ConvertTrackParamForExtrap(trackParam, forwardBackward, v3);
249     do { // reduce step lenght while zEnd oversteped
250       if (stepNumber > fgkMaxStepNumber) {
251         cout<<"W-AliMUONTrackExtrap::ExtrapToZRungekutta: Too many trials: "<<stepNumber<<endl;
252         break;
253       }
254       stepNumber ++;
255       step = TMath::Abs(step);
256       AliMUONTrackExtrap::ExtrapOneStepRungekutta(chargeExtrap,step,v3,v3New);
257       residue = zEnd - v3New[2];
258       step *= dZ/(v3New[2]-trackParam->GetZ());
259     } while (residue*dZ < 0 && TMath::Abs(residue) > fgkRungeKuttaMaxResidue);
260     RecoverTrackParam(v3New, chargeExtrap * forwardBackward, trackParam);
261   }
262   
263   // terminate the extropolation with a straight line up to the exact "zEnd" value
264   trackParam->SetNonBendingCoor(trackParam->GetNonBendingCoor() + residue * trackParam->GetNonBendingSlope());
265   trackParam->SetBendingCoor(trackParam->GetBendingCoor() + residue * trackParam->GetBendingSlope());
266   trackParam->SetZ(zEnd);
267 }
268
269 //__________________________________________________________________________
270 void AliMUONTrackExtrap::ConvertTrackParamForExtrap(AliMUONTrackParam* trackParam, Double_t forwardBackward, Double_t *v3)
271 {
272   /// Set vector of Geant3 parameters pointed to by "v3" from track parameters in trackParam.
273   /// Since AliMUONTrackParam is only geometry, one uses "forwardBackward"
274   /// to know whether the particle is going forward (+1) or backward (-1).
275   v3[0] = trackParam->GetNonBendingCoor(); // X
276   v3[1] = trackParam->GetBendingCoor(); // Y
277   v3[2] = trackParam->GetZ(); // Z
278   Double_t pYZ = TMath::Abs(1.0 / trackParam->GetInverseBendingMomentum());
279   Double_t pZ = pYZ / TMath::Sqrt(1.0 + trackParam->GetBendingSlope() * trackParam->GetBendingSlope());
280   v3[6] = TMath::Sqrt(pYZ * pYZ + pZ * pZ * trackParam->GetNonBendingSlope() * trackParam->GetNonBendingSlope()); // PTOT
281   v3[5] = -forwardBackward * pZ / v3[6]; // PZ/PTOT spectro. z<0
282   v3[3] = trackParam->GetNonBendingSlope() * v3[5]; // PX/PTOT
283   v3[4] = trackParam->GetBendingSlope() * v3[5]; // PY/PTOT
284 }
285
286 //__________________________________________________________________________
287 void AliMUONTrackExtrap::RecoverTrackParam(Double_t *v3, Double_t charge, AliMUONTrackParam* trackParam)
288 {
289   /// Set track parameters in trackParam from Geant3 parameters pointed to by "v3",
290   /// assumed to be calculated for forward motion in Z.
291   /// "InverseBendingMomentum" is signed with "charge".
292   trackParam->SetNonBendingCoor(v3[0]); // X
293   trackParam->SetBendingCoor(v3[1]); // Y
294   trackParam->SetZ(v3[2]); // Z
295   Double_t pYZ = v3[6] * TMath::Sqrt(1.0 - v3[3] * v3[3]);
296   trackParam->SetInverseBendingMomentum(charge/pYZ);
297   trackParam->SetBendingSlope(v3[4]/v3[5]);
298   trackParam->SetNonBendingSlope(v3[3]/v3[5]);
299 }
300
301 //__________________________________________________________________________
302 void AliMUONTrackExtrap::ExtrapToZCov(AliMUONTrackParam* trackParam, Double_t zEnd, Bool_t updatePropagator)
303 {
304   /// Track parameters and their covariances extrapolated to the plane at "zEnd".
305   /// On return, results from the extrapolation are updated in trackParam.
306   
307   if (trackParam->GetZ() == zEnd) return; // nothing to be done if same z
308   
309   if (!fgFieldON) { // linear extrapolation if no magnetic field
310     AliMUONTrackExtrap::LinearExtrapToZ(trackParam,zEnd,updatePropagator);
311     return;
312   }
313   
314   // No need to propagate the covariance matrix if it does not exist
315   if (!trackParam->CovariancesExist()) {
316     cout<<"W-AliMUONTrackExtrap::ExtrapToZCov: Covariance matrix does not exist"<<endl;
317     // Extrapolate track parameters to "zEnd"
318     ExtrapToZ(trackParam,zEnd);
319     return;
320   }
321   
322   // Save the actual track parameters
323   AliMUONTrackParam trackParamSave(*trackParam);
324   TMatrixD paramSave(trackParamSave.GetParameters());
325   Double_t zBegin = trackParamSave.GetZ();
326   
327   // Get reference to the parameter covariance matrix
328   const TMatrixD& kParamCov = trackParam->GetCovariances();
329         
330   // Extrapolate track parameters to "zEnd"
331   ExtrapToZ(trackParam,zEnd);
332   
333   // Get reference to the extrapolated parameters
334   const TMatrixD& extrapParam = trackParam->GetParameters();
335   
336   // Calculate the jacobian related to the track parameters extrapolation to "zEnd"
337   TMatrixD jacob(5,5);
338   jacob.Zero();
339   TMatrixD dParam(5,1);
340   for (Int_t i=0; i<5; i++) {
341     // Skip jacobian calculation for parameters with no associated error
342     if (kParamCov(i,i) <= 0.) continue;
343     
344     // Small variation of parameter i only
345     for (Int_t j=0; j<5; j++) {
346       if (j==i) {
347         dParam(j,0) = TMath::Sqrt(kParamCov(i,i));
348         if (j == 4) dParam(j,0) *= TMath::Sign(1.,-paramSave(4,0)); // variation always in the same direction
349       } else dParam(j,0) = 0.;
350     }
351     
352     // Set new parameters
353     trackParamSave.SetParameters(paramSave);
354     trackParamSave.AddParameters(dParam);
355     trackParamSave.SetZ(zBegin);
356     
357     // Extrapolate new track parameters to "zEnd"
358     ExtrapToZ(&trackParamSave,zEnd);
359     
360     // Calculate the jacobian
361     TMatrixD jacobji(trackParamSave.GetParameters(),TMatrixD::kMinus,extrapParam);
362     jacobji *= 1. / dParam(i,0);
363     jacob.SetSub(0,i,jacobji);
364   }
365   
366   // Extrapolate track parameter covariances to "zEnd"
367   TMatrixD tmp(kParamCov,TMatrixD::kMultTranspose,jacob);
368   TMatrixD tmp2(jacob,TMatrixD::kMult,tmp);
369   trackParam->SetCovariances(tmp2);
370   
371   // Update the propagator if required
372   if (updatePropagator) trackParam->UpdatePropagator(jacob);
373 }
374
375 //__________________________________________________________________________
376 void AliMUONTrackExtrap::AddMCSEffectInAbsorber(AliMUONTrackParam* param, Double_t pathLength, Double_t f0, Double_t f1, Double_t f2)
377 {
378   /// Add to the track parameter covariances the effects of multiple Coulomb scattering
379   /// The absorber correction parameters are supposed to be calculated at the current track z-position
380   
381   // absorber related covariance parameters
382   Double_t bendingSlope = param->GetBendingSlope();
383   Double_t nonBendingSlope = param->GetNonBendingSlope();
384   Double_t inverseBendingMomentum = param->GetInverseBendingMomentum();
385   Double_t alpha2 = 0.0136 * 0.0136 * inverseBendingMomentum * inverseBendingMomentum * (1.0 + bendingSlope * bendingSlope) /
386                     (1.0 + bendingSlope *bendingSlope + nonBendingSlope * nonBendingSlope); // velocity = 1
387   Double_t varCoor = alpha2 * (pathLength * pathLength * f0 - 2. * pathLength * f1 + f2);
388   Double_t covCorrSlope = alpha2 * (pathLength * f0 - f1);
389   Double_t varSlop = alpha2 * f0;
390   
391   // compute derivative d(q/Pxy) / dSlopeX and d(q/Pxy) / dSlopeX
392   Double_t dqPxydSlopeX = inverseBendingMomentum * nonBendingSlope / (1. + nonBendingSlope*nonBendingSlope + bendingSlope*bendingSlope);
393   Double_t dqPxydSlopeY = - inverseBendingMomentum * nonBendingSlope*nonBendingSlope * bendingSlope /
394                             (1. + bendingSlope*bendingSlope) / (1. + nonBendingSlope*nonBendingSlope + bendingSlope*bendingSlope);
395   
396   // Set MCS covariance matrix
397   TMatrixD newParamCov(param->GetCovariances());
398   // Non bending plane
399   newParamCov(0,0) += varCoor;       newParamCov(0,1) += covCorrSlope;
400   newParamCov(1,0) += covCorrSlope;  newParamCov(1,1) += varSlop;
401   // Bending plane
402   newParamCov(2,2) += varCoor;       newParamCov(2,3) += covCorrSlope;
403   newParamCov(3,2) += covCorrSlope;  newParamCov(3,3) += varSlop;
404   // Inverse bending momentum (due to dependences with bending and non bending slopes)
405   newParamCov(4,0) += dqPxydSlopeX * covCorrSlope; newParamCov(0,4) += dqPxydSlopeX * covCorrSlope;
406   newParamCov(4,1) += dqPxydSlopeX * varSlop;      newParamCov(1,4) += dqPxydSlopeX * varSlop;
407   newParamCov(4,2) += dqPxydSlopeY * covCorrSlope; newParamCov(2,4) += dqPxydSlopeY * covCorrSlope;
408   newParamCov(4,3) += dqPxydSlopeY * varSlop;      newParamCov(3,4) += dqPxydSlopeY * varSlop;
409   newParamCov(4,4) += (dqPxydSlopeX*dqPxydSlopeX + dqPxydSlopeY*dqPxydSlopeY) * varSlop;
410   
411   // Set new covariances
412   param->SetCovariances(newParamCov);
413 }
414
415 //__________________________________________________________________________
416 void AliMUONTrackExtrap::CorrectMCSEffectInAbsorber(AliMUONTrackParam* param,
417                                                     Double_t xVtx, Double_t yVtx, Double_t zVtx,
418                                                     Double_t errXVtx, Double_t errYVtx,
419                                                     Double_t absZBeg, Double_t pathLength, Double_t f0, Double_t f1, Double_t f2)
420 {
421   /// Correct parameters and corresponding covariances using Branson correction
422   /// - input param are parameters and covariances at the end of absorber
423   /// - output param are parameters and covariances at vertex
424   /// Absorber correction parameters are supposed to be calculated at the current track z-position
425   
426   // Position of the Branson plane (spectro. (z<0))
427   Double_t zB = (f1>0.) ? absZBeg - f2/f1 : 0.;
428   
429   // Add MCS effects to current parameter covariances
430   AddMCSEffectInAbsorber(param, pathLength, f0, f1, f2);
431   
432   // Get track parameters and covariances in the Branson plane corrected for magnetic field effect
433   ExtrapToZCov(param,zVtx);
434   LinearExtrapToZ(param,zB);
435   
436   // compute track parameters at vertex
437   TMatrixD newParam(5,1);
438   newParam(0,0) = xVtx;
439   newParam(1,0) = (param->GetNonBendingCoor() - xVtx) / (zB - zVtx);
440   newParam(2,0) = yVtx;
441   newParam(3,0) = (param->GetBendingCoor() - yVtx) / (zB - zVtx);
442   newParam(4,0) = param->GetCharge() / param->P() *
443                   TMath::Sqrt(1.0 + newParam(1,0)*newParam(1,0) + newParam(3,0)*newParam(3,0)) /
444                   TMath::Sqrt(1.0 + newParam(3,0)*newParam(3,0));
445   
446   // Get covariances in (X, SlopeX, Y, SlopeY, q*PTot) coordinate system
447   TMatrixD paramCovP(param->GetCovariances());
448   Cov2CovP(param->GetParameters(),paramCovP);
449   
450   // Get the covariance matrix in the (XVtx, X, YVtx, Y, q*PTot) coordinate system
451   TMatrixD paramCovVtx(5,5);
452   paramCovVtx.Zero();
453   paramCovVtx(0,0) = errXVtx * errXVtx;
454   paramCovVtx(1,1) = paramCovP(0,0);
455   paramCovVtx(2,2) = errYVtx * errYVtx;
456   paramCovVtx(3,3) = paramCovP(2,2);
457   paramCovVtx(4,4) = paramCovP(4,4);
458   paramCovVtx(1,3) = paramCovP(0,2);
459   paramCovVtx(3,1) = paramCovP(2,0);
460   paramCovVtx(1,4) = paramCovP(0,4);
461   paramCovVtx(4,1) = paramCovP(4,0);
462   paramCovVtx(3,4) = paramCovP(2,4);
463   paramCovVtx(4,3) = paramCovP(4,2);
464   
465   // Jacobian of the transformation (XVtx, X, YVtx, Y, q*PTot) -> (XVtx, SlopeXVtx, YVtx, SlopeYVtx, q*PTotVtx)
466   TMatrixD jacob(5,5);
467   jacob.UnitMatrix();
468   jacob(1,0) = - 1. / (zB - zVtx);
469   jacob(1,1) = 1. / (zB - zVtx);
470   jacob(3,2) = - 1. / (zB - zVtx);
471   jacob(3,3) = 1. / (zB - zVtx);
472   
473   // Compute covariances at vertex in the (XVtx, SlopeXVtx, YVtx, SlopeYVtx, q*PTotVtx) coordinate system
474   TMatrixD tmp(paramCovVtx,TMatrixD::kMultTranspose,jacob);
475   TMatrixD newParamCov(jacob,TMatrixD::kMult,tmp);
476   
477   // Compute covariances at vertex in the (XVtx, SlopeXVtx, YVtx, SlopeYVtx, q/PyzVtx) coordinate system
478   CovP2Cov(newParam,newParamCov);
479   
480   // Set parameters and covariances at vertex
481   param->SetParameters(newParam);
482   param->SetZ(zVtx);
483   param->SetCovariances(newParamCov);
484 }
485
486 //__________________________________________________________________________
487 void AliMUONTrackExtrap::CorrectELossEffectInAbsorber(AliMUONTrackParam* param, Double_t eLoss, Double_t sigmaELoss2)
488 {
489   /// Correct parameters for energy loss and add energy loss fluctuation effect to covariances
490   
491   // Get parameter covariances in (X, SlopeX, Y, SlopeY, q*PTot) coordinate system
492   TMatrixD newParamCov(param->GetCovariances());
493   Cov2CovP(param->GetParameters(),newParamCov);
494   
495   // Add effects of energy loss fluctuation to covariances
496   newParamCov(4,4) += sigmaELoss2;
497   
498   // Compute new parameters corrected for energy loss
499   Double_t nonBendingSlope = param->GetNonBendingSlope();
500   Double_t bendingSlope = param->GetBendingSlope();
501   param->SetInverseBendingMomentum(param->GetCharge() / (param->P() + eLoss) *
502                                    TMath::Sqrt(1.0 + nonBendingSlope*nonBendingSlope + bendingSlope*bendingSlope) /
503                                    TMath::Sqrt(1.0 + bendingSlope*bendingSlope));
504   
505   // Get new parameter covariances in (X, SlopeX, Y, SlopeY, q/Pyz) coordinate system
506   CovP2Cov(param->GetParameters(),newParamCov);
507   
508   // Set new parameter covariances
509   param->SetCovariances(newParamCov);
510 }
511
512 //__________________________________________________________________________
513 Bool_t AliMUONTrackExtrap::GetAbsorberCorrectionParam(Double_t trackXYZIn[3], Double_t trackXYZOut[3], Double_t pTotal,
514                                                       Double_t &pathLength, Double_t &f0, Double_t &f1, Double_t &f2,
515                                                       Double_t &meanRho, Double_t &totalELoss, Double_t &sigmaELoss2)
516 {
517   /// Parameters used to correct for Multiple Coulomb Scattering and energy loss in absorber
518   /// Calculated assuming a linear propagation from trackXYZIn to trackXYZOut (order is important)
519   // pathLength: path length between trackXYZIn and trackXYZOut (cm)
520   // f0:         0th moment of z calculated with the inverse radiation-length distribution
521   // f1:         1st moment of z calculated with the inverse radiation-length distribution
522   // f2:         2nd moment of z calculated with the inverse radiation-length distribution
523   // meanRho:    average density of crossed material (g/cm3)
524   // totalELoss: total energy loss in absorber
525   
526   // Reset absorber's parameters
527   pathLength = 0.;
528   f0 = 0.;
529   f1 = 0.;
530   f2 = 0.;
531   meanRho = 0.;
532   totalELoss = 0.;
533   sigmaELoss2 = 0.;
534   
535   // Check whether the geometry is available
536   if (!gGeoManager) {
537     cout<<"E-AliMUONTrackExtrap::GetAbsorberCorrectionParam: no TGeo"<<endl;
538     return kFALSE;
539   }
540   
541   // Initialize starting point and direction
542   pathLength = TMath::Sqrt((trackXYZOut[0] - trackXYZIn[0])*(trackXYZOut[0] - trackXYZIn[0])+
543                            (trackXYZOut[1] - trackXYZIn[1])*(trackXYZOut[1] - trackXYZIn[1])+
544                            (trackXYZOut[2] - trackXYZIn[2])*(trackXYZOut[2] - trackXYZIn[2]));
545   if (pathLength < TGeoShape::Tolerance()) return kFALSE;
546   Double_t b[3];
547   b[0] = (trackXYZOut[0] - trackXYZIn[0]) / pathLength;
548   b[1] = (trackXYZOut[1] - trackXYZIn[1]) / pathLength;
549   b[2] = (trackXYZOut[2] - trackXYZIn[2]) / pathLength;
550   TGeoNode *currentnode = gGeoManager->InitTrack(trackXYZIn, b);
551   if (!currentnode) {
552     cout<<"E-AliMUONTrackExtrap::GetAbsorberCorrectionParam: start point out of geometry"<<endl;
553     return kFALSE;
554   }
555   
556   // loop over absorber slices and calculate absorber's parameters
557   Double_t rho = 0.; // material density (g/cm3)
558   Double_t x0 = 0.;  // radiation-length (cm-1)
559   Double_t atomicA = 0.; // A of material
560   Double_t atomicZ = 0.; // Z of material
561   Double_t localPathLength = 0;
562   Double_t remainingPathLength = pathLength;
563   Double_t zB = trackXYZIn[2];
564   Double_t zE, dzB, dzE;
565   do {
566     // Get material properties
567     TGeoMaterial *material = currentnode->GetVolume()->GetMedium()->GetMaterial();
568     rho = material->GetDensity();
569     x0 = material->GetRadLen();
570     if (!material->IsMixture()) x0 /= rho; // different normalization in the modeler for mixture
571     atomicA = material->GetA();
572     atomicZ = material->GetZ();
573     
574     // Get path length within this material
575     gGeoManager->FindNextBoundary(remainingPathLength);
576     localPathLength = gGeoManager->GetStep() + 1.e-6;
577     // Check if boundary within remaining path length. If so, make sure to cross the boundary to prepare the next step
578     if (localPathLength >= remainingPathLength) localPathLength = remainingPathLength;
579     else {
580       currentnode = gGeoManager->Step();
581       if (!currentnode) {
582         cout<<"E-AliMUONTrackExtrap::GetAbsorberCorrectionParam: navigation failed"<<endl;
583         f0 = f1 = f2 = meanRho = totalELoss = sigmaELoss2 = 0.;
584         return kFALSE;
585       }
586       if (!gGeoManager->IsEntering()) {
587         // make another small step to try to enter in new absorber slice
588         gGeoManager->SetStep(0.001);
589         currentnode = gGeoManager->Step();
590         if (!gGeoManager->IsEntering() || !currentnode) {
591           cout<<"E-AliMUONTrackExtrap::GetAbsorberCorrectionParam: navigation failed"<<endl;
592           f0 = f1 = f2 = meanRho = totalELoss = sigmaELoss2 = 0.;
593           return kFALSE;
594         }
595         localPathLength += 0.001;
596       }
597     }
598     
599     // calculate absorber's parameters
600     zE = b[2] * localPathLength + zB;
601     dzB = zB - trackXYZIn[2];
602     dzE = zE - trackXYZIn[2];
603     f0 += localPathLength / x0;
604     f1 += (dzE*dzE - dzB*dzB) / b[2] / b[2] / x0 / 2.;
605     f2 += (dzE*dzE*dzE - dzB*dzB*dzB) / b[2] / b[2] / b[2] / x0 / 3.;
606     meanRho += localPathLength * rho;
607     totalELoss += BetheBloch(pTotal, localPathLength, rho, atomicA, atomicZ);
608     sigmaELoss2 += EnergyLossFluctuation2(pTotal, localPathLength, rho, atomicA, atomicZ);
609     
610     // prepare next step
611     zB = zE;
612     remainingPathLength -= localPathLength;
613   } while (remainingPathLength > TGeoShape::Tolerance());
614   
615   meanRho /= pathLength;
616   
617   return kTRUE;
618 }
619
620 //__________________________________________________________________________
621 Double_t AliMUONTrackExtrap::GetMCSAngle2(const AliMUONTrackParam& param, Double_t dZ, Double_t x0)
622 {
623   /// Return the angular dispersion square due to multiple Coulomb scattering
624   /// through a material of thickness "dZ" and of radiation length "x0"
625   /// assuming linear propagation and using the small angle approximation.
626   
627   Double_t bendingSlope = param.GetBendingSlope();
628   Double_t nonBendingSlope = param.GetNonBendingSlope();
629   Double_t inverseTotalMomentum2 = param.GetInverseBendingMomentum() * param.GetInverseBendingMomentum() *
630                                    (1.0 + bendingSlope * bendingSlope) /
631                                    (1.0 + bendingSlope *bendingSlope + nonBendingSlope * nonBendingSlope); 
632   // Path length in the material
633   Double_t pathLength = TMath::Abs(dZ) * TMath::Sqrt(1.0 + bendingSlope*bendingSlope + nonBendingSlope*nonBendingSlope);
634   // relativistic velocity
635   Double_t velo = 1.;
636   // Angular dispersion square of the track (variance) in a plane perpendicular to the trajectory
637   Double_t theta02 = 0.0136 / velo * (1 + 0.038 * TMath::Log(pathLength/x0));
638   
639   return theta02 * theta02 * inverseTotalMomentum2 * pathLength / x0;
640 }
641
642 //__________________________________________________________________________
643 void AliMUONTrackExtrap::AddMCSEffect(AliMUONTrackParam *param, Double_t dZ, Double_t x0)
644 {
645   /// Add to the track parameter covariances the effects of multiple Coulomb scattering
646   /// through a material of thickness "dZ" and of radiation length "x0"
647   /// assuming linear propagation and using the small angle approximation.
648   
649   Double_t bendingSlope = param->GetBendingSlope();
650   Double_t nonBendingSlope = param->GetNonBendingSlope();
651   Double_t inverseBendingMomentum = param->GetInverseBendingMomentum();
652   Double_t inverseTotalMomentum2 = inverseBendingMomentum * inverseBendingMomentum *
653                                    (1.0 + bendingSlope * bendingSlope) /
654                                    (1.0 + bendingSlope *bendingSlope + nonBendingSlope * nonBendingSlope); 
655   // Path length in the material
656   Double_t pathLength = TMath::Abs(dZ) * TMath::Sqrt(1.0 + bendingSlope*bendingSlope + nonBendingSlope*nonBendingSlope);
657   Double_t pathLength2 = pathLength * pathLength;
658   // relativistic velocity
659   Double_t velo = 1.;
660   // Angular dispersion square of the track (variance) in a plane perpendicular to the trajectory
661   Double_t theta02 = 0.0136 / velo * (1 + 0.038 * TMath::Log(pathLength/x0));
662   theta02 *= theta02 * inverseTotalMomentum2 * pathLength / x0;
663   
664   Double_t varCoor      = pathLength2 * theta02 / 3.;
665   Double_t varSlop      = theta02;
666   Double_t covCorrSlope = pathLength * theta02 / 2.;
667   
668   // compute derivative d(q/Pxy) / dSlopeX and d(q/Pxy) / dSlopeX
669   Double_t dqPxydSlopeX = inverseBendingMomentum * nonBendingSlope / (1. + nonBendingSlope*nonBendingSlope + bendingSlope*bendingSlope);
670   Double_t dqPxydSlopeY = - inverseBendingMomentum * nonBendingSlope*nonBendingSlope * bendingSlope /
671                             (1. + bendingSlope*bendingSlope) / (1. + nonBendingSlope*nonBendingSlope + bendingSlope*bendingSlope);
672   
673   // Set MCS covariance matrix
674   TMatrixD newParamCov(param->GetCovariances());
675   // Non bending plane
676   newParamCov(0,0) += varCoor;       newParamCov(0,1) += covCorrSlope;
677   newParamCov(1,0) += covCorrSlope;  newParamCov(1,1) += varSlop;
678   // Bending plane
679   newParamCov(2,2) += varCoor;       newParamCov(2,3) += covCorrSlope;
680   newParamCov(3,2) += covCorrSlope;  newParamCov(3,3) += varSlop;
681   // Inverse bending momentum (due to dependences with bending and non bending slopes)
682   newParamCov(4,0) += dqPxydSlopeX * covCorrSlope; newParamCov(0,4) += dqPxydSlopeX * covCorrSlope;
683   newParamCov(4,1) += dqPxydSlopeX * varSlop;      newParamCov(1,4) += dqPxydSlopeX * varSlop;
684   newParamCov(4,2) += dqPxydSlopeY * covCorrSlope; newParamCov(2,4) += dqPxydSlopeY * covCorrSlope;
685   newParamCov(4,3) += dqPxydSlopeY * varSlop;      newParamCov(3,4) += dqPxydSlopeY * varSlop;
686   newParamCov(4,4) += (dqPxydSlopeX*dqPxydSlopeX + dqPxydSlopeY*dqPxydSlopeY) * varSlop;
687   
688   // Set new covariances
689   param->SetCovariances(newParamCov);
690 }
691
692 //__________________________________________________________________________
693 void AliMUONTrackExtrap::ExtrapToVertex(AliMUONTrackParam* trackParam,
694                                         Double_t xVtx, Double_t yVtx, Double_t zVtx,
695                                         Double_t errXVtx, Double_t errYVtx,
696                                         Bool_t correctForMCS, Bool_t correctForEnergyLoss)
697 {
698   /// Main method for extrapolation to the vertex:
699   /// Returns the track parameters and covariances resulting from the extrapolation of the current trackParam
700   /// Changes parameters and covariances according to multiple scattering and energy loss corrections:
701   /// if correctForMCS=kTRUE:  compute parameters using Branson correction and add correction resolution to covariances
702   /// if correctForMCS=kFALSE: add parameter dispersion due to MCS in parameter covariances
703   /// if correctForEnergyLoss=kTRUE:  correct parameters for energy loss and add energy loss fluctuation to covariances
704   /// if correctForEnergyLoss=kFALSE: do nothing about energy loss
705   
706   if (trackParam->GetZ() == zVtx) return; // nothing to be done if already at vertex
707   
708   if (trackParam->GetZ() > zVtx) { // spectro. (z<0)
709     cout<<"E-AliMUONTrackExtrap::ExtrapToVertex: Starting Z ("<<trackParam->GetZ()
710         <<") upstream the vertex (zVtx = "<<zVtx<<")"<<endl;
711     return;
712   }
713   
714   // Check the vertex position relatively to the absorber
715   if (zVtx < AliMUONConstants::AbsZBeg() && zVtx > AliMUONConstants::AbsZEnd()) { // spectro. (z<0)
716     cout<<"W-AliMUONTrackExtrap::ExtrapToVertex: Ending Z ("<<zVtx
717         <<") inside the front absorber ("<<AliMUONConstants::AbsZBeg()<<","<<AliMUONConstants::AbsZEnd()<<")"<<endl;
718   } else if (zVtx < AliMUONConstants::AbsZEnd() ) { // spectro. (z<0)
719     cout<<"W-AliMUONTrackExtrap::ExtrapToVertex: Ending Z ("<<zVtx
720         <<") downstream the front absorber (zAbsorberEnd = "<<AliMUONConstants::AbsZEnd()<<")"<<endl;
721     if (trackParam->CovariancesExist()) ExtrapToZCov(trackParam,zVtx);
722     else ExtrapToZ(trackParam,zVtx);
723     return;
724   }
725   
726   // Check the track position relatively to the absorber and extrapolate track parameters to the end of the absorber if needed
727   if (trackParam->GetZ() > AliMUONConstants::AbsZBeg()) { // spectro. (z<0)
728     cout<<"W-AliMUONTrackExtrap::ExtrapToVertex: Starting Z ("<<trackParam->GetZ()
729         <<") upstream the front absorber (zAbsorberBegin = "<<AliMUONConstants::AbsZBeg()<<")"<<endl;
730     if (trackParam->CovariancesExist()) ExtrapToZCov(trackParam,zVtx);
731     else ExtrapToZ(trackParam,zVtx);
732     return;
733   } else if (trackParam->GetZ() > AliMUONConstants::AbsZEnd()) { // spectro. (z<0)
734     cout<<"W-AliMUONTrackExtrap::ExtrapToVertex: Starting Z ("<<trackParam->GetZ()
735         <<") inside the front absorber ("<<AliMUONConstants::AbsZBeg()<<","<<AliMUONConstants::AbsZEnd()<<")"<<endl;
736   } else {
737     if (trackParam->CovariancesExist()) ExtrapToZCov(trackParam,AliMUONConstants::AbsZEnd());
738     else ExtrapToZ(trackParam,AliMUONConstants::AbsZEnd());
739   }
740   
741   // Get absorber correction parameters assuming linear propagation in absorber
742   Double_t trackXYZOut[3];
743   trackXYZOut[0] = trackParam->GetNonBendingCoor();
744   trackXYZOut[1] = trackParam->GetBendingCoor();
745   trackXYZOut[2] = trackParam->GetZ();
746   Double_t trackXYZIn[3];
747   if (correctForMCS) { // assume linear propagation until the vertex
748     trackXYZIn[2] = TMath::Min(zVtx, AliMUONConstants::AbsZBeg()); // spectro. (z<0)
749     trackXYZIn[0] = trackXYZOut[0] + (xVtx - trackXYZOut[0]) / (zVtx - trackXYZOut[2]) * (trackXYZIn[2] - trackXYZOut[2]);
750     trackXYZIn[1] = trackXYZOut[1] + (yVtx - trackXYZOut[1]) / (zVtx - trackXYZOut[2]) * (trackXYZIn[2] - trackXYZOut[2]);
751   } else {
752     AliMUONTrackParam trackParamIn(*trackParam);
753     ExtrapToZ(&trackParamIn, TMath::Min(zVtx, AliMUONConstants::AbsZBeg()));
754     trackXYZIn[0] = trackParamIn.GetNonBendingCoor();
755     trackXYZIn[1] = trackParamIn.GetBendingCoor();
756     trackXYZIn[2] = trackParamIn.GetZ();
757   }
758   Double_t pTot = trackParam->P();
759   Double_t pathLength, f0, f1, f2, meanRho, deltaP, sigmaDeltaP2;
760   if (!GetAbsorberCorrectionParam(trackXYZIn,trackXYZOut,pTot,pathLength,f0,f1,f2,meanRho,deltaP,sigmaDeltaP2)) {
761     cout<<"E-AliMUONTrackExtrap::ExtrapToVertex: Unable to take into account the absorber effects"<<endl;
762     if (trackParam->CovariancesExist()) ExtrapToZCov(trackParam,zVtx);
763     else ExtrapToZ(trackParam,zVtx);
764     return;
765   }
766   
767   // Compute track parameters and covariances at vertex according to correctForMCS and correctForEnergyLoss flags
768   if (correctForMCS) {
769     
770     if (correctForEnergyLoss) {
771       
772       // Correct for multiple scattering and energy loss
773       CorrectELossEffectInAbsorber(trackParam, 0.5*deltaP, 0.5*sigmaDeltaP2);
774       CorrectMCSEffectInAbsorber(trackParam, xVtx, yVtx, zVtx, errXVtx, errYVtx,
775                                  trackXYZIn[2], pathLength, f0, f1, f2);
776       CorrectELossEffectInAbsorber(trackParam, 0.5*deltaP, 0.5*sigmaDeltaP2);
777       
778     } else {
779       
780       // Correct for multiple scattering
781       CorrectMCSEffectInAbsorber(trackParam, xVtx, yVtx, zVtx, errXVtx, errYVtx,
782                                  trackXYZIn[2], pathLength, f0, f1, f2);
783     }
784     
785   } else {
786     
787     if (correctForEnergyLoss) {
788       
789       // Correct for energy loss add multiple scattering dispersion in covariance matrix
790       CorrectELossEffectInAbsorber(trackParam, 0.5*deltaP, 0.5*sigmaDeltaP2);
791       AddMCSEffectInAbsorber(trackParam, pathLength, f0, f1, f2);
792       ExtrapToZCov(trackParam, trackXYZIn[2]);
793       CorrectELossEffectInAbsorber(trackParam, 0.5*deltaP, 0.5*sigmaDeltaP2);
794       ExtrapToZCov(trackParam, zVtx);
795       
796     } else {
797       
798       // add multiple scattering dispersion in covariance matrix
799       AddMCSEffectInAbsorber(trackParam, pathLength, f0, f1, f2);
800       ExtrapToZCov(trackParam, zVtx);
801       
802     }
803     
804   }
805   
806 }
807
808 //__________________________________________________________________________
809 void AliMUONTrackExtrap::ExtrapToVertex(AliMUONTrackParam* trackParam,
810                                         Double_t xVtx, Double_t yVtx, Double_t zVtx,
811                                         Double_t errXVtx, Double_t errYVtx)
812 {
813   /// Extrapolate track parameters to vertex, corrected for multiple scattering and energy loss effects
814   /// Add branson correction resolution and energy loss fluctuation to parameter covariances
815   ExtrapToVertex(trackParam, xVtx, yVtx, zVtx, errXVtx, errYVtx, kTRUE, kTRUE);
816 }
817
818 //__________________________________________________________________________
819 void AliMUONTrackExtrap::ExtrapToVertexWithoutELoss(AliMUONTrackParam* trackParam,
820                                                     Double_t xVtx, Double_t yVtx, Double_t zVtx,
821                                                     Double_t errXVtx, Double_t errYVtx)
822 {
823   /// Extrapolate track parameters to vertex, corrected for multiple scattering effects only
824   /// Add branson correction resolution to parameter covariances
825   ExtrapToVertex(trackParam, xVtx, yVtx, zVtx, errXVtx, errYVtx, kTRUE, kFALSE);
826 }
827
828 //__________________________________________________________________________
829 void AliMUONTrackExtrap::ExtrapToVertexWithoutBranson(AliMUONTrackParam* trackParam, Double_t zVtx)
830 {
831   /// Extrapolate track parameters to vertex, corrected for energy loss effects only
832   /// Add dispersion due to multiple scattering and energy loss fluctuation to parameter covariances
833   ExtrapToVertex(trackParam, 0., 0., zVtx, 0., 0., kFALSE, kTRUE);
834 }
835
836 //__________________________________________________________________________
837 void AliMUONTrackExtrap::ExtrapToVertexUncorrected(AliMUONTrackParam* trackParam, Double_t zVtx)
838 {
839   /// Extrapolate track parameters to vertex without multiple scattering and energy loss corrections
840   /// Add dispersion due to multiple scattering to parameter covariances
841   ExtrapToVertex(trackParam, 0., 0., zVtx, 0., 0., kFALSE, kFALSE);
842 }
843
844 //__________________________________________________________________________
845 Double_t AliMUONTrackExtrap::TotalMomentumEnergyLoss(AliMUONTrackParam* trackParam, Double_t xVtx, Double_t yVtx, Double_t zVtx)
846 {
847   /// Calculate the total momentum energy loss in-between the track position and the vertex assuming a linear propagation
848   
849   if (trackParam->GetZ() == zVtx) return 0.; // nothing to be done if already at vertex
850   
851   // Check whether the geometry is available
852   if (!gGeoManager) {
853     cout<<"E-AliMUONTrackExtrap::TotalMomentumEnergyLoss: no TGeo"<<endl;
854     return 0.;
855   }
856   
857   // Get encountered material correction parameters assuming linear propagation from vertex to the track position
858   Double_t trackXYZOut[3];
859   trackXYZOut[0] = trackParam->GetNonBendingCoor();
860   trackXYZOut[1] = trackParam->GetBendingCoor();
861   trackXYZOut[2] = trackParam->GetZ();
862   Double_t trackXYZIn[3];
863   trackXYZIn[0] = xVtx;
864   trackXYZIn[1] = yVtx;
865   trackXYZIn[2] = zVtx;
866   Double_t pTot = trackParam->P();
867   Double_t pathLength, f0, f1, f2, meanRho, totalELoss, sigmaELoss2;
868   GetAbsorberCorrectionParam(trackXYZIn,trackXYZOut,pTot,pathLength,f0,f1,f2,meanRho,totalELoss,sigmaELoss2);
869   
870   return totalELoss;
871 }
872
873 //__________________________________________________________________________
874 Double_t AliMUONTrackExtrap::BetheBloch(Double_t pTotal, Double_t pathLength, Double_t rho, Double_t atomicA, Double_t atomicZ)
875 {
876   /// Returns the mean total momentum energy loss of muon with total momentum='pTotal'
877   /// in the absorber layer of lenght='pathLength', density='rho', A='atomicA' and Z='atomicZ'
878   Double_t muMass = 0.105658369; // GeV
879   Double_t eMass = 0.510998918e-3; // GeV
880   Double_t k = 0.307075e-3; // GeV.g^-1.cm^2
881   Double_t i = 9.5e-9; // mean exitation energy per atomic Z (GeV)
882   Double_t p2=pTotal*pTotal;
883   Double_t beta2=p2/(p2 + muMass*muMass);
884   
885   Double_t w = k * rho * pathLength * atomicZ / atomicA / beta2;
886   
887   if (beta2/(1-beta2)>3.5*3.5)
888     return w * (log(2.*eMass*3.5/(i*atomicZ)) + 0.5*log(beta2/(1-beta2)) - beta2);
889   
890   return w * (log(2.*eMass*beta2/(1-beta2)/(i*atomicZ)) - beta2);
891 }
892
893 //__________________________________________________________________________
894 Double_t AliMUONTrackExtrap::EnergyLossFluctuation2(Double_t pTotal, Double_t pathLength, Double_t rho, Double_t atomicA, Double_t atomicZ)
895 {
896   /// Returns the total momentum energy loss fluctuation of muon with total momentum='pTotal'
897   /// in the absorber layer of lenght='pathLength', density='rho', A='atomicA' and Z='atomicZ'
898   Double_t muMass = 0.105658369; // GeV
899   //Double_t eMass = 0.510998918e-3; // GeV
900   Double_t k = 0.307075e-3; // GeV.g^-1.cm^2
901   Double_t p2=pTotal*pTotal;
902   Double_t beta2=p2/(p2 + muMass*muMass);
903   
904   Double_t fwhm = 2. * k * rho * pathLength * atomicZ / atomicA / beta2; // FWHM of the energy loss Landau distribution
905   Double_t sigma2 = fwhm * fwhm / (8.*log(2.)); // gaussian: fwmh = 2 * srqt(2*ln(2)) * sigma (i.e. fwmh = 2.35 * sigma)
906   
907   //sigma2 = k * rho * pathLength * atomicZ / atomicA * eMass; // sigma2 of the energy loss gaussian distribution
908   
909   return sigma2;
910 }
911
912 //__________________________________________________________________________
913 void AliMUONTrackExtrap::Cov2CovP(const TMatrixD &param, TMatrixD &cov)
914 {
915   /// change coordinate system: (X, SlopeX, Y, SlopeY, q/Pyz) -> (X, SlopeX, Y, SlopeY, q*PTot)
916   /// parameters (param) are given in the (X, SlopeX, Y, SlopeY, q/Pyz) coordinate system
917   
918   // charge * total momentum
919   Double_t qPTot = TMath::Sqrt(1. + param(1,0)*param(1,0) + param(3,0)*param(3,0)) /
920                    TMath::Sqrt(1. + param(3,0)*param(3,0)) / param(4,0);
921   
922   // Jacobian of the opposite transformation
923   TMatrixD jacob(5,5);
924   jacob.UnitMatrix();
925   jacob(4,1) = qPTot * param(1,0) / (1. + param(1,0)*param(1,0) + param(3,0)*param(3,0));
926   jacob(4,3) = - qPTot * param(1,0) * param(1,0) * param(3,0) /
927                  (1. + param(3,0)*param(3,0)) / (1. + param(1,0)*param(1,0) + param(3,0)*param(3,0));
928   jacob(4,4) = - qPTot / param(4,0);
929   
930   // compute covariances in new coordinate system
931   TMatrixD tmp(cov,TMatrixD::kMultTranspose,jacob);
932   cov.Mult(jacob,tmp);
933 }
934
935 //__________________________________________________________________________
936 void AliMUONTrackExtrap::CovP2Cov(const TMatrixD &param, TMatrixD &covP)
937 {
938   /// change coordinate system: (X, SlopeX, Y, SlopeY, q*PTot) -> (X, SlopeX, Y, SlopeY, q/Pyz)
939   /// parameters (param) are given in the (X, SlopeX, Y, SlopeY, q/Pyz) coordinate system
940   
941   // charge * total momentum
942   Double_t qPTot = TMath::Sqrt(1. + param(1,0)*param(1,0) + param(3,0)*param(3,0)) /
943                    TMath::Sqrt(1. + param(3,0)*param(3,0)) / param(4,0);
944   
945   // Jacobian of the transformation
946   TMatrixD jacob(5,5);
947   jacob.UnitMatrix();
948   jacob(4,1) = param(4,0) * param(1,0) / (1. + param(1,0)*param(1,0) + param(3,0)*param(3,0));
949   jacob(4,3) = - param(4,0) * param(1,0) * param(1,0) * param(3,0) /
950                  (1. + param(3,0)*param(3,0)) / (1. + param(1,0)*param(1,0) + param(3,0)*param(3,0));
951   jacob(4,4) = - param(4,0) / qPTot;
952   
953   // compute covariances in new coordinate system
954   TMatrixD tmp(covP,TMatrixD::kMultTranspose,jacob);
955   covP.Mult(jacob,tmp);
956 }
957
958  //__________________________________________________________________________
959 void AliMUONTrackExtrap::ExtrapOneStepHelix(Double_t charge, Double_t step, Double_t *vect, Double_t *vout)
960 {
961 /// <pre>
962 ///    ******************************************************************
963 ///    *                                                                *
964 ///    *  Performs the tracking of one step in a magnetic field         *
965 ///    *  The trajectory is assumed to be a helix in a constant field   *
966 ///    *  taken at the mid point of the step.                           *
967 ///    *  Parameters:                                                   *
968 ///    *   input                                                        *
969 ///    *     STEP =arc length of the step asked                         *
970 ///    *     VECT =input vector (position,direction cos and momentum)   *
971 ///    *     CHARGE=  electric charge of the particle                   *
972 ///    *   output                                                       *
973 ///    *     VOUT = same as VECT after completion of the step           *
974 ///    *                                                                *
975 ///    *    ==>Called by : USER, GUSWIM                               *
976 ///    *       Author    m.hansroul  *********                          *
977 ///    *       modified  s.egli, s.v.levonian                           *
978 ///    *       modified  v.perevoztchikov
979 ///    *                                                                *
980 ///    ******************************************************************
981 /// </pre>
982
983 // modif: everything in double precision
984
985     Double_t xyz[3], h[4], hxp[3];
986     Double_t h2xy, hp, rho, tet;
987     Double_t sint, sintt, tsint, cos1t;
988     Double_t f1, f2, f3, f4, f5, f6;
989
990     const Int_t kix  = 0;
991     const Int_t kiy  = 1;
992     const Int_t kiz  = 2;
993     const Int_t kipx = 3;
994     const Int_t kipy = 4;
995     const Int_t kipz = 5;
996     const Int_t kipp = 6;
997
998     const Double_t kec = 2.9979251e-4;
999     //
1000     //    ------------------------------------------------------------------
1001     //
1002     //       units are kgauss,centimeters,gev/c
1003     //
1004     vout[kipp] = vect[kipp];
1005     if (TMath::Abs(charge) < 0.00001) {
1006       for (Int_t i = 0; i < 3; i++) {
1007         vout[i] = vect[i] + step * vect[i+3];
1008         vout[i+3] = vect[i+3];
1009       }
1010       return;
1011     }
1012     xyz[0]    = vect[kix] + 0.5 * step * vect[kipx];
1013     xyz[1]    = vect[kiy] + 0.5 * step * vect[kipy];
1014     xyz[2]    = vect[kiz] + 0.5 * step * vect[kipz];
1015
1016     //cmodif: call gufld (xyz, h) changed into:
1017     GetField (xyz, h);
1018  
1019     h2xy = h[0]*h[0] + h[1]*h[1];
1020     h[3] = h[2]*h[2]+ h2xy;
1021     if (h[3] < 1.e-12) {
1022       for (Int_t i = 0; i < 3; i++) {
1023         vout[i] = vect[i] + step * vect[i+3];
1024         vout[i+3] = vect[i+3];
1025       }
1026       return;
1027     }
1028     if (h2xy < 1.e-12*h[3]) {
1029       ExtrapOneStepHelix3(charge*h[2], step, vect, vout);
1030       return;
1031     }
1032     h[3] = TMath::Sqrt(h[3]);
1033     h[0] /= h[3];
1034     h[1] /= h[3];
1035     h[2] /= h[3];
1036     h[3] *= kec;
1037
1038     hxp[0] = h[1]*vect[kipz] - h[2]*vect[kipy];
1039     hxp[1] = h[2]*vect[kipx] - h[0]*vect[kipz];
1040     hxp[2] = h[0]*vect[kipy] - h[1]*vect[kipx];
1041  
1042     hp = h[0]*vect[kipx] + h[1]*vect[kipy] + h[2]*vect[kipz];
1043
1044     rho = -charge*h[3]/vect[kipp];
1045     tet = rho * step;
1046
1047     if (TMath::Abs(tet) > 0.15) {
1048       sint = TMath::Sin(tet);
1049       sintt = (sint/tet);
1050       tsint = (tet-sint)/tet;
1051       cos1t = 2.*(TMath::Sin(0.5*tet))*(TMath::Sin(0.5*tet))/tet;
1052     } else {
1053       tsint = tet*tet/36.;
1054       sintt = (1. - tsint);
1055       sint = tet*sintt;
1056       cos1t = 0.5*tet;
1057     }
1058
1059     f1 = step * sintt;
1060     f2 = step * cos1t;
1061     f3 = step * tsint * hp;
1062     f4 = -tet*cos1t;
1063     f5 = sint;
1064     f6 = tet * cos1t * hp;
1065  
1066     vout[kix] = vect[kix] + f1*vect[kipx] + f2*hxp[0] + f3*h[0];
1067     vout[kiy] = vect[kiy] + f1*vect[kipy] + f2*hxp[1] + f3*h[1];
1068     vout[kiz] = vect[kiz] + f1*vect[kipz] + f2*hxp[2] + f3*h[2];
1069  
1070     vout[kipx] = vect[kipx] + f4*vect[kipx] + f5*hxp[0] + f6*h[0];
1071     vout[kipy] = vect[kipy] + f4*vect[kipy] + f5*hxp[1] + f6*h[1];
1072     vout[kipz] = vect[kipz] + f4*vect[kipz] + f5*hxp[2] + f6*h[2];
1073  
1074     return;
1075 }
1076
1077  //__________________________________________________________________________
1078 void AliMUONTrackExtrap::ExtrapOneStepHelix3(Double_t field, Double_t step, Double_t *vect, Double_t *vout)
1079 {
1080 /// <pre>
1081 ///     ******************************************************************
1082 ///     *                                                                *
1083 ///     *       Tracking routine in a constant field oriented            *
1084 ///     *       along axis 3                                             *
1085 ///     *       Tracking is performed with a conventional                *
1086 ///     *       helix step method                                        *
1087 ///     *                                                                *
1088 ///     *    ==>Called by : USER, GUSWIM                                 *
1089 ///     *       Authors    R.Brun, M.Hansroul  *********                 *
1090 ///     *       Rewritten  V.Perevoztchikov
1091 ///     *                                                                *
1092 ///     ******************************************************************
1093 /// </pre>
1094
1095     Double_t hxp[3];
1096     Double_t h4, hp, rho, tet;
1097     Double_t sint, sintt, tsint, cos1t;
1098     Double_t f1, f2, f3, f4, f5, f6;
1099
1100     const Int_t kix  = 0;
1101     const Int_t kiy  = 1;
1102     const Int_t kiz  = 2;
1103     const Int_t kipx = 3;
1104     const Int_t kipy = 4;
1105     const Int_t kipz = 5;
1106     const Int_t kipp = 6;
1107
1108     const Double_t kec = 2.9979251e-4;
1109
1110 // 
1111 //     ------------------------------------------------------------------
1112 // 
1113 //       units are kgauss,centimeters,gev/c
1114 // 
1115     vout[kipp] = vect[kipp];
1116     h4 = field * kec;
1117
1118     hxp[0] = - vect[kipy];
1119     hxp[1] = + vect[kipx];
1120  
1121     hp = vect[kipz];
1122
1123     rho = -h4/vect[kipp];
1124     tet = rho * step;
1125     if (TMath::Abs(tet) > 0.15) {
1126       sint = TMath::Sin(tet);
1127       sintt = (sint/tet);
1128       tsint = (tet-sint)/tet;
1129       cos1t = 2.* TMath::Sin(0.5*tet) * TMath::Sin(0.5*tet)/tet;
1130     } else {
1131       tsint = tet*tet/36.;
1132       sintt = (1. - tsint);
1133       sint = tet*sintt;
1134       cos1t = 0.5*tet;
1135     }
1136
1137     f1 = step * sintt;
1138     f2 = step * cos1t;
1139     f3 = step * tsint * hp;
1140     f4 = -tet*cos1t;
1141     f5 = sint;
1142     f6 = tet * cos1t * hp;
1143  
1144     vout[kix] = vect[kix] + f1*vect[kipx] + f2*hxp[0];
1145     vout[kiy] = vect[kiy] + f1*vect[kipy] + f2*hxp[1];
1146     vout[kiz] = vect[kiz] + f1*vect[kipz] + f3;
1147  
1148     vout[kipx] = vect[kipx] + f4*vect[kipx] + f5*hxp[0];
1149     vout[kipy] = vect[kipy] + f4*vect[kipy] + f5*hxp[1];
1150     vout[kipz] = vect[kipz] + f4*vect[kipz] + f6;
1151
1152     return;
1153 }
1154
1155  //__________________________________________________________________________
1156 void AliMUONTrackExtrap::ExtrapOneStepRungekutta(Double_t charge, Double_t step, Double_t* vect, Double_t* vout)
1157 {
1158 /// <pre>
1159 ///     ******************************************************************
1160 ///     *                                                                *
1161 ///     *  Runge-Kutta method for tracking a particle through a magnetic *
1162 ///     *  field. Uses Nystroem algorithm (See Handbook Nat. Bur. of     *
1163 ///     *  Standards, procedure 25.5.20)                                 *
1164 ///     *                                                                *
1165 ///     *  Input parameters                                              *
1166 ///     *       CHARGE    Particle charge                                *
1167 ///     *       STEP      Step size                                      *
1168 ///     *       VECT      Initial co-ords,direction cosines,momentum     *
1169 ///     *  Output parameters                                             *
1170 ///     *       VOUT      Output co-ords,direction cosines,momentum      *
1171 ///     *  User routine called                                           *
1172 ///     *       CALL GUFLD(X,F)                                          *
1173 ///     *                                                                *
1174 ///     *    ==>Called by : USER, GUSWIM                                 *
1175 ///     *       Authors    R.Brun, M.Hansroul  *********                 *
1176 ///     *                  V.Perevoztchikov (CUT STEP implementation)    *
1177 ///     *                                                                *
1178 ///     *                                                                *
1179 ///     ******************************************************************
1180 /// </pre>
1181
1182     Double_t h2, h4, f[4];
1183     Double_t xyzt[3], a, b, c, ph,ph2;
1184     Double_t secxs[4],secys[4],seczs[4],hxp[3];
1185     Double_t g1, g2, g3, g4, g5, g6, ang2, dxt, dyt, dzt;
1186     Double_t est, at, bt, ct, cba;
1187     Double_t f1, f2, f3, f4, rho, tet, hnorm, hp, rho1, sint, cost;
1188     
1189     Double_t x;
1190     Double_t y;
1191     Double_t z;
1192     
1193     Double_t xt;
1194     Double_t yt;
1195     Double_t zt;
1196
1197     Double_t maxit = 1992;
1198     Double_t maxcut = 11;
1199
1200     const Double_t kdlt   = 1e-4;
1201     const Double_t kdlt32 = kdlt/32.;
1202     const Double_t kthird = 1./3.;
1203     const Double_t khalf  = 0.5;
1204     const Double_t kec = 2.9979251e-4;
1205
1206     const Double_t kpisqua = 9.86960440109;
1207     const Int_t kix  = 0;
1208     const Int_t kiy  = 1;
1209     const Int_t kiz  = 2;
1210     const Int_t kipx = 3;
1211     const Int_t kipy = 4;
1212     const Int_t kipz = 5;
1213   
1214     // *.
1215     // *.    ------------------------------------------------------------------
1216     // *.
1217     // *             this constant is for units cm,gev/c and kgauss
1218     // *
1219     Int_t iter = 0;
1220     Int_t ncut = 0;
1221     for(Int_t j = 0; j < 7; j++)
1222       vout[j] = vect[j];
1223
1224     Double_t  pinv   = kec * charge / vect[6];
1225     Double_t tl = 0.;
1226     Double_t h = step;
1227     Double_t rest;
1228
1229  
1230     do {
1231       rest  = step - tl;
1232       if (TMath::Abs(h) > TMath::Abs(rest)) h = rest;
1233       //cmodif: call gufld(vout,f) changed into:
1234
1235       GetField(vout,f);
1236
1237       // *
1238       // *             start of integration
1239       // *
1240       x      = vout[0];
1241       y      = vout[1];
1242       z      = vout[2];
1243       a      = vout[3];
1244       b      = vout[4];
1245       c      = vout[5];
1246
1247       h2     = khalf * h;
1248       h4     = khalf * h2;
1249       ph     = pinv * h;
1250       ph2    = khalf * ph;
1251       secxs[0] = (b * f[2] - c * f[1]) * ph2;
1252       secys[0] = (c * f[0] - a * f[2]) * ph2;
1253       seczs[0] = (a * f[1] - b * f[0]) * ph2;
1254       ang2 = (secxs[0]*secxs[0] + secys[0]*secys[0] + seczs[0]*seczs[0]);
1255       if (ang2 > kpisqua) break;
1256
1257       dxt    = h2 * a + h4 * secxs[0];
1258       dyt    = h2 * b + h4 * secys[0];
1259       dzt    = h2 * c + h4 * seczs[0];
1260       xt     = x + dxt;
1261       yt     = y + dyt;
1262       zt     = z + dzt;
1263       // *
1264       // *              second intermediate point
1265       // *
1266
1267       est = TMath::Abs(dxt) + TMath::Abs(dyt) + TMath::Abs(dzt);
1268       if (est > h) {
1269         if (ncut++ > maxcut) break;
1270         h *= khalf;
1271         continue;
1272       }
1273  
1274       xyzt[0] = xt;
1275       xyzt[1] = yt;
1276       xyzt[2] = zt;
1277
1278       //cmodif: call gufld(xyzt,f) changed into:
1279       GetField(xyzt,f);
1280
1281       at     = a + secxs[0];
1282       bt     = b + secys[0];
1283       ct     = c + seczs[0];
1284
1285       secxs[1] = (bt * f[2] - ct * f[1]) * ph2;
1286       secys[1] = (ct * f[0] - at * f[2]) * ph2;
1287       seczs[1] = (at * f[1] - bt * f[0]) * ph2;
1288       at     = a + secxs[1];
1289       bt     = b + secys[1];
1290       ct     = c + seczs[1];
1291       secxs[2] = (bt * f[2] - ct * f[1]) * ph2;
1292       secys[2] = (ct * f[0] - at * f[2]) * ph2;
1293       seczs[2] = (at * f[1] - bt * f[0]) * ph2;
1294       dxt    = h * (a + secxs[2]);
1295       dyt    = h * (b + secys[2]);
1296       dzt    = h * (c + seczs[2]);
1297       xt     = x + dxt;
1298       yt     = y + dyt;
1299       zt     = z + dzt;
1300       at     = a + 2.*secxs[2];
1301       bt     = b + 2.*secys[2];
1302       ct     = c + 2.*seczs[2];
1303
1304       est = TMath::Abs(dxt)+TMath::Abs(dyt)+TMath::Abs(dzt);
1305       if (est > 2.*TMath::Abs(h)) {
1306         if (ncut++ > maxcut) break;
1307         h *= khalf;
1308         continue;
1309       }
1310  
1311       xyzt[0] = xt;
1312       xyzt[1] = yt;
1313       xyzt[2] = zt;
1314
1315       //cmodif: call gufld(xyzt,f) changed into:
1316       GetField(xyzt,f);
1317
1318       z      = z + (c + (seczs[0] + seczs[1] + seczs[2]) * kthird) * h;
1319       y      = y + (b + (secys[0] + secys[1] + secys[2]) * kthird) * h;
1320       x      = x + (a + (secxs[0] + secxs[1] + secxs[2]) * kthird) * h;
1321
1322       secxs[3] = (bt*f[2] - ct*f[1])* ph2;
1323       secys[3] = (ct*f[0] - at*f[2])* ph2;
1324       seczs[3] = (at*f[1] - bt*f[0])* ph2;
1325       a      = a+(secxs[0]+secxs[3]+2. * (secxs[1]+secxs[2])) * kthird;
1326       b      = b+(secys[0]+secys[3]+2. * (secys[1]+secys[2])) * kthird;
1327       c      = c+(seczs[0]+seczs[3]+2. * (seczs[1]+seczs[2])) * kthird;
1328
1329       est    = TMath::Abs(secxs[0]+secxs[3] - (secxs[1]+secxs[2]))
1330         + TMath::Abs(secys[0]+secys[3] - (secys[1]+secys[2]))
1331         + TMath::Abs(seczs[0]+seczs[3] - (seczs[1]+seczs[2]));
1332
1333       if (est > kdlt && TMath::Abs(h) > 1.e-4) {
1334         if (ncut++ > maxcut) break;
1335         h *= khalf;
1336         continue;
1337       }
1338
1339       ncut = 0;
1340       // *               if too many iterations, go to helix
1341       if (iter++ > maxit) break;
1342
1343       tl += h;
1344       if (est < kdlt32) 
1345         h *= 2.;
1346       cba    = 1./ TMath::Sqrt(a*a + b*b + c*c);
1347       vout[0] = x;
1348       vout[1] = y;
1349       vout[2] = z;
1350       vout[3] = cba*a;
1351       vout[4] = cba*b;
1352       vout[5] = cba*c;
1353       rest = step - tl;
1354       if (step < 0.) rest = -rest;
1355       if (rest < 1.e-5*TMath::Abs(step)) return;
1356
1357     } while(1);
1358
1359     // angle too big, use helix
1360
1361     f1  = f[0];
1362     f2  = f[1];
1363     f3  = f[2];
1364     f4  = TMath::Sqrt(f1*f1+f2*f2+f3*f3);
1365     rho = -f4*pinv;
1366     tet = rho * step;
1367  
1368     hnorm = 1./f4;
1369     f1 = f1*hnorm;
1370     f2 = f2*hnorm;
1371     f3 = f3*hnorm;
1372
1373     hxp[0] = f2*vect[kipz] - f3*vect[kipy];
1374     hxp[1] = f3*vect[kipx] - f1*vect[kipz];
1375     hxp[2] = f1*vect[kipy] - f2*vect[kipx];
1376  
1377     hp = f1*vect[kipx] + f2*vect[kipy] + f3*vect[kipz];
1378
1379     rho1 = 1./rho;
1380     sint = TMath::Sin(tet);
1381     cost = 2.*TMath::Sin(khalf*tet)*TMath::Sin(khalf*tet);
1382
1383     g1 = sint*rho1;
1384     g2 = cost*rho1;
1385     g3 = (tet-sint) * hp*rho1;
1386     g4 = -cost;
1387     g5 = sint;
1388     g6 = cost * hp;
1389  
1390     vout[kix] = vect[kix] + g1*vect[kipx] + g2*hxp[0] + g3*f1;
1391     vout[kiy] = vect[kiy] + g1*vect[kipy] + g2*hxp[1] + g3*f2;
1392     vout[kiz] = vect[kiz] + g1*vect[kipz] + g2*hxp[2] + g3*f3;
1393  
1394     vout[kipx] = vect[kipx] + g4*vect[kipx] + g5*hxp[0] + g6*f1;
1395     vout[kipy] = vect[kipy] + g4*vect[kipy] + g5*hxp[1] + g6*f2;
1396     vout[kipz] = vect[kipz] + g4*vect[kipz] + g5*hxp[2] + g6*f3;
1397
1398     return;
1399 }
1400
1401 //___________________________________________________________
1402 void  AliMUONTrackExtrap::GetField(Double_t *Position, Double_t *Field)
1403 {
1404   /// interface for arguments in double precision (Why ? ChF)
1405   Float_t x[3], b[3];
1406   
1407   x[0] = Position[0]; x[1] = Position[1]; x[2] = Position[2];
1408   
1409   if (fgkField) fgkField->Field(x,b);
1410   else {
1411     cout<<"F-AliMUONTrackExtrap::GetField: fgkField = 0x0"<<endl;
1412     exit(-1);
1413   }
1414   
1415   Field[0] = b[0]; Field[1] = b[1]; Field[2] = b[2];
1416   
1417   return;
1418 }