]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONTrackParam.cxx
Logging of Debug, Info and Error Messages follwing AliRoot Standard http://aliweb...
[u/mrichter/AliRoot.git] / MUON / AliMUONTrackParam.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 ///////////////////////////////////////////////////
19 //
20 // Track parameters
21 // in
22 // ALICE
23 // dimuon
24 // spectrometer
25 //
26 ///////////////////////////////////////////////////
27
28 #include <Riostream.h>
29
30 #include "AliCallf77.h" 
31 #include "AliMUON.h"
32 #include "AliMUONTrackParam.h" 
33 #include "AliMUONChamber.h"
34 #include "AliRun.h" 
35 #include "AliMagF.h" 
36 #include "AliLog.h" 
37
38 ClassImp(AliMUONTrackParam) // Class implementation in ROOT context
39
40   // A few calls in Fortran or from Fortran (extrap.F).
41   // Needed, instead of calls to Geant subroutines,
42   // because double precision is necessary for track fit converging with Minuit.
43   // The "extrap" functions should be translated into C++ ????
44 #ifndef WIN32 
45 # define extrap_onestep_helix extrap_onestep_helix_
46 # define extrap_onestep_helix3 extrap_onestep_helix3_
47 # define extrap_onestep_rungekutta extrap_onestep_rungekutta_
48 # define gufld_double gufld_double_
49 #else 
50 # define extrap_onestep_helix EXTRAP_ONESTEP_HELIX
51 # define extrap_onestep_helix3 EXTRAP_ONESTEP_HELIX3
52 # define extrap_onestep_rungekutta EXTRAP_ONESTEP_RUNGEKUTTA
53 # define gufld_double GUFLD_DOUBLE
54 #endif 
55
56 extern "C" {
57   void type_of_call extrap_onestep_helix
58   (Double_t &Charge, Double_t &StepLength, Double_t *VGeant3, Double_t *VGeant3New);
59
60   void type_of_call extrap_onestep_helix3
61   (Double_t &Field, Double_t &StepLength, Double_t *VGeant3, Double_t *VGeant3New);
62
63   void type_of_call extrap_onestep_rungekutta
64   (Double_t &Charge, Double_t &StepLength, Double_t *VGeant3, Double_t *VGeant3New);
65
66   void type_of_call gufld_double(Double_t *Position, Double_t *Field) {
67     // interface to "gAlice->Field()->Field" for arguments in double precision
68     Float_t x[3], b[3];
69     x[0] = Position[0]; x[1] = Position[1]; x[2] = Position[2];
70     gAlice->Field()->Field(x, b);
71     Field[0] = b[0]; Field[1] = b[1]; Field[2] = b[2];
72   }
73 }
74
75   //_________________________________________________________________________
76 AliMUONTrackParam::AliMUONTrackParam()
77   : TObject()
78 {
79 // Constructor
80
81   fInverseBendingMomentum = 0;
82   fBendingSlope = 0;
83   fNonBendingSlope = 0;
84   fZ = 0;
85   fBendingCoor = 0;
86   fNonBendingCoor = 0;
87 }
88
89   //_________________________________________________________________________
90 AliMUONTrackParam& 
91 AliMUONTrackParam::operator=(const AliMUONTrackParam& theMUONTrackParam)
92 {
93   if (this == &theMUONTrackParam)
94     return *this;
95
96   // base class assignement
97   TObject::operator=(theMUONTrackParam);
98
99   fInverseBendingMomentum =  theMUONTrackParam.fInverseBendingMomentum; 
100   fBendingSlope           =  theMUONTrackParam.fBendingSlope; 
101   fNonBendingSlope        =  theMUONTrackParam.fNonBendingSlope; 
102   fZ                      =  theMUONTrackParam.fZ; 
103   fBendingCoor            =  theMUONTrackParam.fBendingCoor; 
104   fNonBendingCoor         =  theMUONTrackParam.fNonBendingCoor;
105
106   return *this;
107 }
108   //_________________________________________________________________________
109 AliMUONTrackParam::AliMUONTrackParam(const AliMUONTrackParam& theMUONTrackParam)
110   : TObject(theMUONTrackParam)
111 {
112   fInverseBendingMomentum =  theMUONTrackParam.fInverseBendingMomentum; 
113   fBendingSlope           =  theMUONTrackParam.fBendingSlope; 
114   fNonBendingSlope        =  theMUONTrackParam.fNonBendingSlope; 
115   fZ                      =  theMUONTrackParam.fZ; 
116   fBendingCoor            =  theMUONTrackParam.fBendingCoor; 
117   fNonBendingCoor         =  theMUONTrackParam.fNonBendingCoor;
118 }
119
120   //__________________________________________________________________________
121 void AliMUONTrackParam::ExtrapToZ(Double_t Z)
122 {
123   // Track parameter extrapolation to the plane at "Z".
124   // On return, the track parameters resulting from the extrapolation
125   // replace the current track parameters.
126   if (this->fZ == Z) return; // nothing to be done if same Z
127   Double_t forwardBackward; // +1 if forward, -1 if backward
128   if (Z < this->fZ) forwardBackward = 1.0; // spectro. z<0 
129   else forwardBackward = -1.0;
130   Double_t vGeant3[7], vGeant3New[7]; // 7 in parameter ????
131   Int_t iGeant3, stepNumber;
132   Int_t maxStepNumber = 5000; // in parameter ????
133   // For safety: return kTRUE or kFALSE ????
134   // Parameter vector for calling EXTRAP_ONESTEP
135   SetGeant3Parameters(vGeant3, forwardBackward);
136   // sign of charge (sign of fInverseBendingMomentum if forward motion)
137   // must be changed if backward extrapolation
138   Double_t chargeExtrap = forwardBackward *
139     TMath::Sign(Double_t(1.0), this->fInverseBendingMomentum);
140   Double_t stepLength = 6.0; // in parameter ????
141   // Extrapolation loop
142   stepNumber = 0;
143   while (((-forwardBackward * (vGeant3[2] - Z)) <= 0.0) &&  // spectro. z<0
144          (stepNumber < maxStepNumber)) {
145     stepNumber++;
146     // Option for switching between helix and Runge-Kutta ???? 
147     // extrap_onestep_rungekutta(chargeExtrap, stepLength, vGeant3, vGeant3New);
148     extrap_onestep_helix(chargeExtrap, stepLength, vGeant3, vGeant3New);
149     if ((-forwardBackward * (vGeant3New[2] - Z)) > 0.0) break; // one is beyond Z spectro. z<0
150     // better use TArray ????
151     for (iGeant3 = 0; iGeant3 < 7; iGeant3++)
152       {vGeant3[iGeant3] = vGeant3New[iGeant3];}
153   }
154   // check maxStepNumber ????
155   // Interpolation back to exact Z (2nd order)
156   // should be in function ???? using TArray ????
157   Double_t dZ12 = vGeant3New[2] - vGeant3[2]; // 1->2
158   Double_t dZ1i = Z - vGeant3[2]; // 1-i
159   Double_t dZi2 = vGeant3New[2] - Z; // i->2
160   Double_t xPrime = (vGeant3New[0] - vGeant3[0]) / dZ12;
161   Double_t xSecond =
162     ((vGeant3New[3] / vGeant3New[5]) - (vGeant3[3] / vGeant3[5])) / dZ12;
163   Double_t yPrime = (vGeant3New[1] - vGeant3[1]) / dZ12;
164   Double_t ySecond =
165     ((vGeant3New[4] / vGeant3New[5]) - (vGeant3[4] / vGeant3[5])) / dZ12;
166   vGeant3[0] = vGeant3[0] + xPrime * dZ1i - 0.5 * xSecond * dZ1i * dZi2; // X
167   vGeant3[1] = vGeant3[1] + yPrime * dZ1i - 0.5 * ySecond * dZ1i * dZi2; // Y
168   vGeant3[2] = Z; // Z
169   Double_t xPrimeI = xPrime - 0.5 * xSecond * (dZi2 - dZ1i);
170   Double_t yPrimeI = yPrime - 0.5 * ySecond * (dZi2 - dZ1i);
171   // (PX, PY, PZ)/PTOT assuming forward motion
172   vGeant3[5] =
173     1.0 / TMath::Sqrt(1.0 + xPrimeI * xPrimeI + yPrimeI * yPrimeI); // PZ/PTOT
174   vGeant3[3] = xPrimeI * vGeant3[5]; // PX/PTOT
175   vGeant3[4] = yPrimeI * vGeant3[5]; // PY/PTOT
176   // Track parameters from Geant3 parameters,
177   // with charge back for forward motion
178   GetFromGeant3Parameters(vGeant3, chargeExtrap * forwardBackward);
179 }
180
181   //__________________________________________________________________________
182 void AliMUONTrackParam::SetGeant3Parameters(Double_t *VGeant3, Double_t ForwardBackward)
183 {
184   // Set vector of Geant3 parameters pointed to by "VGeant3"
185   // from track parameters in current AliMUONTrackParam.
186   // Since AliMUONTrackParam is only geometry, one uses "ForwardBackward"
187   // to know whether the particle is going forward (+1) or backward (-1).
188   VGeant3[0] = this->fNonBendingCoor; // X
189   VGeant3[1] = this->fBendingCoor; // Y
190   VGeant3[2] = this->fZ; // Z
191   Double_t pYZ = TMath::Abs(1.0 / this->fInverseBendingMomentum);
192   Double_t pZ =
193     pYZ / TMath::Sqrt(1.0 + this->fBendingSlope * this->fBendingSlope);
194   VGeant3[6] =
195     TMath::Sqrt(pYZ * pYZ +
196                 pZ * pZ * this->fNonBendingSlope * this->fNonBendingSlope); // PTOT
197   VGeant3[5] = -ForwardBackward * pZ / VGeant3[6]; // PZ/PTOT spectro. z<0
198   VGeant3[3] = this->fNonBendingSlope * VGeant3[5]; // PX/PTOT
199   VGeant3[4] = this->fBendingSlope * VGeant3[5]; // PY/PTOT
200 }
201
202   //__________________________________________________________________________
203 void AliMUONTrackParam::GetFromGeant3Parameters(Double_t *VGeant3, Double_t Charge)
204 {
205   // Get track parameters in current AliMUONTrackParam
206   // from Geant3 parameters pointed to by "VGeant3",
207   // assumed to be calculated for forward motion in Z.
208   // "InverseBendingMomentum" is signed with "Charge".
209   this->fNonBendingCoor = VGeant3[0]; // X
210   this->fBendingCoor = VGeant3[1]; // Y
211   this->fZ = VGeant3[2]; // Z
212   Double_t pYZ = VGeant3[6] * TMath::Sqrt(1.0 - VGeant3[3] * VGeant3[3]);
213   this->fInverseBendingMomentum = Charge / pYZ;
214   this->fBendingSlope = VGeant3[4] / VGeant3[5];
215   this->fNonBendingSlope = VGeant3[3] / VGeant3[5];
216 }
217
218   //__________________________________________________________________________
219 void AliMUONTrackParam::ExtrapToStation(Int_t Station, AliMUONTrackParam *TrackParam)
220 {
221   // Track parameters extrapolated from current track parameters ("this")
222   // to both chambers of the station(0..) "Station"
223   // are returned in the array (dimension 2) of track parameters
224   // pointed to by "TrackParam" (index 0 and 1 for first and second chambers).
225   Double_t extZ[2], z1, z2;
226   Int_t i1 = -1, i2 = -1; // = -1 to avoid compilation warnings
227   AliMUON *pMUON = (AliMUON*) gAlice->GetModule("MUON"); // necessary ????
228   // range of Station to be checked ????
229   z1 = (&(pMUON->Chamber(2 * Station)))->Z(); // Z of first chamber
230   z2 = (&(pMUON->Chamber(2 * Station + 1)))->Z(); // Z of second chamber
231   // First and second Z to extrapolate at
232   if ((z1 > this->fZ) && (z2 > this->fZ)) {i1 = 0; i2 = 1;}
233   else if ((z1 < this->fZ) && (z2 < this->fZ)) {i1 = 1; i2 = 0;}
234   else {
235         AliError(Form("Starting Z (%f) in between z1 (%f) and z2 (%f) of station(0..)%d",this->fZ,z1,z2,Station));
236 //     cout << "ERROR in AliMUONTrackParam::CreateExtrapSegmentInStation" << endl;
237 //     cout << "Starting Z (" << this->fZ << ") in between z1 (" << z1 <<
238 //       ") and z2 (" << z2 << ") of station(0..) " << Station << endl;
239   }
240   extZ[i1] = z1;
241   extZ[i2] = z2;
242   // copy of track parameters
243   TrackParam[i1] = *this;
244   // first extrapolation
245   (&(TrackParam[i1]))->ExtrapToZ(extZ[0]);
246   TrackParam[i2] = TrackParam[i1];
247   // second extrapolation
248   (&(TrackParam[i2]))->ExtrapToZ(extZ[1]);
249   return;
250 }
251
252   //__________________________________________________________________________
253 void AliMUONTrackParam::ExtrapToVertex()
254 {
255   // Extrapolation to the vertex.
256   // Returns the track parameters resulting from the extrapolation,
257   // in the current TrackParam.
258   // Changes parameters according to Branson correction through the absorber 
259   
260   Double_t zAbsorber = -503.0; // to be coherent with the Geant absorber geometry !!!!
261                                // spectro. (z<0) 
262   // Extrapolates track parameters upstream to the "Z" end of the front absorber
263   ExtrapToZ(zAbsorber); // !!!
264   // Makes Branson correction (multiple scattering + energy loss)
265   BransonCorrection();
266   // Makes a simple magnetic field correction through the absorber
267   FieldCorrection(zAbsorber);
268 }
269
270
271 //  Keep this version for future developments
272   //__________________________________________________________________________
273 // void AliMUONTrackParam::BransonCorrection()
274 // {
275 //   // Branson correction of track parameters
276 //   // the entry parameters have to be calculated at the end of the absorber
277 //   Double_t zEndAbsorber, zBP, xBP, yBP;
278 //   Double_t  pYZ, pX, pY, pZ, pTotal, xEndAbsorber, yEndAbsorber, radiusEndAbsorber2, pT, theta;
279 //   Int_t sign;
280 //   // Would it be possible to calculate all that from Geant configuration ????
281 //   // and to get the Branson parameters from a function in ABSO module ????
282 //   // with an eventual contribution from other detectors like START ????
283 //   // Radiation lengths outer part theta > 3 degres
284 //   static Double_t x01[9] = { 18.8,    // C (cm)
285 //                           10.397,   // Concrete (cm)
286 //                           0.56,    // Plomb (cm)
287 //                           47.26,   // Polyethylene (cm)
288 //                           0.56,   // Plomb (cm)
289 //                           47.26,   // Polyethylene (cm)
290 //                           0.56,   // Plomb (cm)
291 //                           47.26,   // Polyethylene (cm)
292 //                           0.56 };   // Plomb (cm)
293 //   // inner part theta < 3 degres
294 //   static Double_t x02[3] = { 18.8,    // C (cm)
295 //                           10.397,   // Concrete (cm)
296 //                           0.35 };    // W (cm) 
297 //   // z positions of the materials inside the absober outer part theta > 3 degres
298 //   static Double_t z1[10] = { 90, 315, 467, 472, 477, 482, 487, 492, 497, 502 };
299 //   // inner part theta < 3 degres
300 //   static Double_t z2[4] = { 90, 315, 467, 503 };
301 //   static Bool_t first = kTRUE;
302 //   static Double_t zBP1, zBP2, rLimit;
303 //   // Calculates z positions of the Branson's planes: zBP1 for outer part and zBP2 for inner part (only at the first call)
304 //   if (first) {
305 //     first = kFALSE;
306 //     Double_t aNBP = 0.0;
307 //     Double_t aDBP = 0.0;
308 //     Int_t iBound;
309     
310 //     for (iBound = 0; iBound < 9; iBound++) {
311 //       aNBP = aNBP +
312 //      (z1[iBound+1] * z1[iBound+1] * z1[iBound+1] -
313 //       z1[iBound]   * z1[iBound]   * z1[iBound]    ) / x01[iBound];
314 //       aDBP = aDBP +
315 //      (z1[iBound+1] * z1[iBound+1] - z1[iBound]   * z1[iBound]    ) / x01[iBound];
316 //     }
317 //     zBP1 = (2.0 * aNBP) / (3.0 * aDBP);
318 //     aNBP = 0.0;
319 //     aDBP = 0.0;
320 //     for (iBound = 0; iBound < 3; iBound++) {
321 //       aNBP = aNBP +
322 //      (z2[iBound+1] * z2[iBound+1] * z2[iBound+1] -
323 //       z2[iBound]   * z2[iBound ]  * z2[iBound]    ) / x02[iBound];
324 //       aDBP = aDBP +
325 //      (z2[iBound+1] * z2[iBound+1] - z2[iBound] * z2[iBound]) / x02[iBound];
326 //     }
327 //     zBP2 = (2.0 * aNBP) / (3.0 * aDBP);
328 //     rLimit = z2[3] * TMath::Tan(3.0 * (TMath::Pi()) / 180.);
329 //   }
330
331 //   pYZ = TMath::Abs(1.0 / fInverseBendingMomentum);
332 //   sign = 1;      
333 //   if (fInverseBendingMomentum < 0) sign = -1;     
334 //   pZ = pYZ / (TMath::Sqrt(1.0 + fBendingSlope * fBendingSlope)); 
335 //   pX = pZ * fNonBendingSlope; 
336 //   pY = pZ * fBendingSlope; 
337 //   pTotal = TMath::Sqrt(pYZ *pYZ + pX * pX);
338 //   xEndAbsorber = fNonBendingCoor; 
339 //   yEndAbsorber = fBendingCoor; 
340 //   radiusEndAbsorber2 = xEndAbsorber * xEndAbsorber + yEndAbsorber * yEndAbsorber;
341
342 //   if (radiusEndAbsorber2 > rLimit*rLimit) {
343 //     zEndAbsorber = z1[9];
344 //     zBP = zBP1;
345 //   } else {
346 //     zEndAbsorber = z2[3];
347 //     zBP = zBP2;
348 //   }
349
350 //   xBP = xEndAbsorber - (pX / pZ) * (zEndAbsorber - zBP);
351 //   yBP = yEndAbsorber - (pY / pZ) * (zEndAbsorber - zBP);
352
353 //   // new parameters after Branson and energy loss corrections
354 //   pZ = pTotal * zBP / TMath::Sqrt(xBP * xBP + yBP * yBP + zBP * zBP);
355 //   pX = pZ * xBP / zBP;
356 //   pY = pZ * yBP / zBP;
357 //   fBendingSlope = pY / pZ;
358 //   fNonBendingSlope = pX / pZ;
359   
360 //   pT = TMath::Sqrt(pX * pX + pY * pY);      
361 //   theta = TMath::ATan2(pT, pZ); 
362 //   pTotal =
363 //     TotalMomentumEnergyLoss(rLimit, pTotal, theta, xEndAbsorber, yEndAbsorber);
364
365 //   fInverseBendingMomentum = (sign / pTotal) *
366 //     TMath::Sqrt(1.0 +
367 //              fBendingSlope * fBendingSlope +
368 //              fNonBendingSlope * fNonBendingSlope) /
369 //     TMath::Sqrt(1.0 + fBendingSlope * fBendingSlope);
370
371 //   // vertex position at (0,0,0)
372 //   // should be taken from vertex measurement ???
373 //   fBendingCoor = 0.0;
374 //   fNonBendingCoor = 0;
375 //   fZ= 0;
376 // }
377
378 void AliMUONTrackParam::BransonCorrection()
379 {
380   // Branson correction of track parameters
381   // the entry parameters have to be calculated at the end of the absorber
382   // simplified version: the z positions of Branson's planes are no longer calculated
383   // but are given as inputs. One can use the macros MUONTestAbso.C and DrawTestAbso.C
384   // to test this correction. 
385   // Would it be possible to calculate all that from Geant configuration ????
386   // and to get the Branson parameters from a function in ABSO module ????
387   // with an eventual contribution from other detectors like START ????
388   Double_t  zBP, xBP, yBP;
389   Double_t  pYZ, pX, pY, pZ, pTotal, xEndAbsorber, yEndAbsorber, radiusEndAbsorber2, pT, theta;
390   Int_t sign;
391   static Bool_t first = kTRUE;
392   static Double_t zBP1, zBP2, rLimit, thetaLimit, zEndAbsorber;
393   // zBP1 for outer part and zBP2 for inner part (only at the first call)
394   if (first) {
395     first = kFALSE;
396   
397     zEndAbsorber = -503;  // spectro (z<0)
398     thetaLimit = 3.0 * (TMath::Pi()) / 180.;
399     rLimit = TMath::Abs(zEndAbsorber) * TMath::Tan(thetaLimit);
400     zBP1 = -450; // values close to those calculated with EvalAbso.C
401     zBP2 = -480;
402   }
403
404   pYZ = TMath::Abs(1.0 / fInverseBendingMomentum);
405   sign = 1;      
406   if (fInverseBendingMomentum < 0) sign = -1;     
407   pZ = -pYZ / (TMath::Sqrt(1.0 + fBendingSlope * fBendingSlope)); // spectro (z<0)
408   pX = pZ * fNonBendingSlope; 
409   pY = pZ * fBendingSlope; 
410   pTotal = TMath::Sqrt(pYZ *pYZ + pX * pX);
411   xEndAbsorber = fNonBendingCoor; 
412   yEndAbsorber = fBendingCoor; 
413   radiusEndAbsorber2 = xEndAbsorber * xEndAbsorber + yEndAbsorber * yEndAbsorber;
414
415   if (radiusEndAbsorber2 > rLimit*rLimit) {
416     zBP = zBP1;
417   } else {
418     zBP = zBP2;
419   }
420
421   xBP = xEndAbsorber - (pX / pZ) * (zEndAbsorber - zBP);
422   yBP = yEndAbsorber - (pY / pZ) * (zEndAbsorber - zBP);
423
424   // new parameters after Branson and energy loss corrections
425 //   Float_t zSmear = zBP - gRandom->Gaus(0.,2.);  // !!! possible smearing of Z vertex position
426   Float_t zSmear = zBP;
427   
428   pZ = pTotal * zSmear / TMath::Sqrt(xBP * xBP + yBP * yBP + zSmear * zSmear);
429   pX = pZ * xBP / zSmear;
430   pY = pZ * yBP / zSmear;
431   fBendingSlope = pY / pZ;
432   fNonBendingSlope = pX / pZ;
433
434   
435   pT = TMath::Sqrt(pX * pX + pY * pY);      
436   theta = TMath::ATan2(pT, TMath::Abs(pZ)); 
437   pTotal = TotalMomentumEnergyLoss(thetaLimit, pTotal, theta);
438
439   fInverseBendingMomentum = (sign / pTotal) *
440     TMath::Sqrt(1.0 +
441                 fBendingSlope * fBendingSlope +
442                 fNonBendingSlope * fNonBendingSlope) /
443     TMath::Sqrt(1.0 + fBendingSlope * fBendingSlope);
444
445   // vertex position at (0,0,0)
446   // should be taken from vertex measurement ???
447   fBendingCoor = 0.0;
448   fNonBendingCoor = 0;
449   fZ= 0;
450 }
451
452   //__________________________________________________________________________
453 Double_t AliMUONTrackParam::TotalMomentumEnergyLoss(Double_t thetaLimit, Double_t pTotal, Double_t theta)
454 {
455   // Returns the total momentum corrected from energy loss in the front absorber
456   // One can use the macros MUONTestAbso.C and DrawTestAbso.C
457   // to test this correction. 
458   // Momentum energy loss behaviour evaluated with the simulation of single muons (april 2002)
459   Double_t deltaP, pTotalCorrected;
460
461    // Parametrization to be redone according to change of absorber material ????
462   // See remark in function BransonCorrection !!!!
463   // The name is not so good, and there are many arguments !!!!
464   if (theta  < thetaLimit ) {
465     if (pTotal < 20) {
466       deltaP = 2.5938 + 0.0570 * pTotal - 0.001151 * pTotal * pTotal;
467     } else {
468       deltaP = 3.0714 + 0.011767 *pTotal;
469     }
470   } else {
471     if (pTotal < 20) {
472       deltaP  = 2.1207 + 0.05478 * pTotal - 0.00145079 * pTotal * pTotal;
473     } else { 
474       deltaP = 2.6069 + 0.0051705 * pTotal;
475     }
476   }
477   pTotalCorrected = pTotal + deltaP / TMath::Cos(theta);
478   return pTotalCorrected;
479 }
480
481   //__________________________________________________________________________
482 void AliMUONTrackParam::FieldCorrection(Double_t Z)
483 {
484   // 
485   // Correction of the effect of the magnetic field in the absorber
486   // Assume a constant field along Z axis.
487
488   Float_t b[3],x[3]; 
489   Double_t bZ;
490   Double_t pYZ,pX,pY,pZ,pT;
491   Double_t pXNew,pYNew;
492   Double_t c;
493
494   pYZ = TMath::Abs(1.0 / fInverseBendingMomentum);
495   c = TMath::Sign(1.0,fInverseBendingMomentum); // particle charge 
496  
497   pZ = -pYZ / (TMath::Sqrt(1.0 + fBendingSlope * fBendingSlope));  // spectro. (z<0)
498   pX = pZ * fNonBendingSlope; 
499   pY = pZ * fBendingSlope;
500   pT = TMath::Sqrt(pX*pX+pY*pY);
501
502   if (TMath::Abs(pZ) <= 0) return;
503   x[2] = Z/2;
504   x[0] = x[2]*fNonBendingSlope;  
505   x[1] = x[2]*fBendingSlope;
506
507   // Take magn. field value at position x.
508   gAlice->Field()->Field(x, b);
509   bZ =  b[2];
510  
511   // Transverse momentum rotation
512   // Parameterized with the study of DeltaPhi = phiReco - phiGen as a function of pZ.
513   Double_t phiShift = c*0.436*0.0003*bZ*Z/pZ; 
514  // Rotate momentum around Z axis.
515   pXNew = pX*TMath::Cos(phiShift) - pY*TMath::Sin(phiShift);
516   pYNew = pX*TMath::Sin(phiShift) + pY*TMath::Cos(phiShift);
517  
518   fBendingSlope = pYNew / pZ;
519   fNonBendingSlope = pXNew / pZ;
520   
521   fInverseBendingMomentum = c / TMath::Sqrt(pYNew*pYNew+pZ*pZ);
522  
523 }