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