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