]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONTrackParam.cxx
Using AliTracker (external mag field) instead of gAlice->Field (Nicolas & Christian)
[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 #include "AliMUON.h"
30 #include "AliMUONTrackParam.h" 
31 #include "AliMUONConstants.h"
32 #include "AliESDMuonTrack.h"
33 #include "AliMagF.h" 
34 #include "AliLog.h" 
35 #include "AliTracker.h"
36
37 ClassImp(AliMUONTrackParam) // Class implementation in ROOT context
38
39   //_________________________________________________________________________
40 AliMUONTrackParam::AliMUONTrackParam()
41   : TObject(),
42     fInverseBendingMomentum(0.),
43     fBendingSlope(0.),
44     fNonBendingSlope(0.),
45     fZ(0.),
46     fBendingCoor(0.),
47     fNonBendingCoor(0.)
48 {
49 // Constructor
50   // get field from outside
51   fkField = AliTracker::GetFieldMap();
52   if (!fkField) AliFatal("No field available");
53 }
54
55   //_________________________________________________________________________
56 AliMUONTrackParam& 
57 AliMUONTrackParam::operator=(const AliMUONTrackParam& theMUONTrackParam)
58 {
59   // Asignment operator
60   if (this == &theMUONTrackParam)
61     return *this;
62
63   // base class assignement
64   TObject::operator=(theMUONTrackParam);
65
66   fInverseBendingMomentum =  theMUONTrackParam.fInverseBendingMomentum; 
67   fBendingSlope           =  theMUONTrackParam.fBendingSlope; 
68   fNonBendingSlope        =  theMUONTrackParam.fNonBendingSlope; 
69   fZ                      =  theMUONTrackParam.fZ; 
70   fBendingCoor            =  theMUONTrackParam.fBendingCoor; 
71   fNonBendingCoor         =  theMUONTrackParam.fNonBendingCoor;
72
73   return *this;
74 }
75   //_________________________________________________________________________
76 AliMUONTrackParam::AliMUONTrackParam(const AliMUONTrackParam& theMUONTrackParam)
77   : TObject(theMUONTrackParam),
78     fInverseBendingMomentum(theMUONTrackParam.fInverseBendingMomentum), 
79     fBendingSlope(theMUONTrackParam.fBendingSlope),
80     fNonBendingSlope(theMUONTrackParam.fNonBendingSlope),
81     fZ(theMUONTrackParam.fZ),
82     fBendingCoor(theMUONTrackParam.fBendingCoor),
83     fNonBendingCoor(theMUONTrackParam.fNonBendingCoor)
84 {
85   // Copy constructor
86
87 }
88
89   //_________________________________________________________________________
90 void AliMUONTrackParam::GetParamFrom(const AliESDMuonTrack& esdMuonTrack)
91 {
92   // assigned value form ESD track.
93   fInverseBendingMomentum =  esdMuonTrack.GetInverseBendingMomentum();
94   fBendingSlope           =  TMath::Tan(esdMuonTrack.GetThetaY());
95   fNonBendingSlope        =  TMath::Tan(esdMuonTrack.GetThetaX());
96   fZ                      =  esdMuonTrack.GetZ(); 
97   fBendingCoor            =  esdMuonTrack.GetBendingCoor(); 
98   fNonBendingCoor         =  esdMuonTrack.GetNonBendingCoor();
99 }
100
101   //_________________________________________________________________________
102 void AliMUONTrackParam::SetParamFor(AliESDMuonTrack& esdMuonTrack)
103 {
104   // assigned value form ESD track.
105   esdMuonTrack.SetInverseBendingMomentum(fInverseBendingMomentum);
106   esdMuonTrack.SetThetaX(TMath::ATan(fNonBendingSlope));
107   esdMuonTrack.SetThetaY(TMath::ATan(fBendingSlope));
108   esdMuonTrack.SetZ(fZ); 
109   esdMuonTrack.SetBendingCoor(fBendingCoor); 
110   esdMuonTrack.SetNonBendingCoor(fNonBendingCoor);
111 }
112
113   //__________________________________________________________________________
114 void AliMUONTrackParam::ExtrapToZ(Double_t Z)
115 {
116   // Track parameter extrapolation to the plane at "Z".
117   // On return, the track parameters resulting from the extrapolation
118   // replace the current track parameters.
119   if (this->fZ == Z) return; // nothing to be done if same Z
120   Double_t forwardBackward; // +1 if forward, -1 if backward
121   if (Z < this->fZ) forwardBackward = 1.0; // spectro. z<0 
122   else forwardBackward = -1.0;
123   Double_t vGeant3[7], vGeant3New[7]; // 7 in parameter ????
124   Int_t iGeant3, stepNumber;
125   Int_t maxStepNumber = 5000; // in parameter ????
126   // For safety: return kTRUE or kFALSE ????
127   // Parameter vector for calling EXTRAP_ONESTEP
128   SetGeant3Parameters(vGeant3, forwardBackward);
129   // sign of charge (sign of fInverseBendingMomentum if forward motion)
130   // must be changed if backward extrapolation
131   Double_t chargeExtrap = forwardBackward *
132     TMath::Sign(Double_t(1.0), this->fInverseBendingMomentum);
133   Double_t stepLength = 6.0; // in parameter ????
134   // Extrapolation loop
135   stepNumber = 0;
136   while (((-forwardBackward * (vGeant3[2] - Z)) <= 0.0) &&  // spectro. z<0
137          (stepNumber < maxStepNumber)) {
138     stepNumber++;
139     // Option for switching between helix and Runge-Kutta ???? 
140     //ExtrapOneStepRungekutta(chargeExtrap, stepLength, vGeant3, vGeant3New);
141     ExtrapOneStepHelix(chargeExtrap, stepLength, vGeant3, vGeant3New);
142     if ((-forwardBackward * (vGeant3New[2] - Z)) > 0.0) break; // one is beyond Z spectro. z<0
143     // better use TArray ????
144     for (iGeant3 = 0; iGeant3 < 7; iGeant3++)
145       {vGeant3[iGeant3] = vGeant3New[iGeant3];}
146   }
147   // check maxStepNumber ????
148   // Interpolation back to exact Z (2nd order)
149   // should be in function ???? using TArray ????
150   Double_t dZ12 = vGeant3New[2] - vGeant3[2]; // 1->2
151   if (TMath::Abs(dZ12) > 0) {
152     Double_t dZ1i = Z - vGeant3[2]; // 1-i
153     Double_t dZi2 = vGeant3New[2] - Z; // i->2
154     Double_t xPrime = (vGeant3New[0] - vGeant3[0]) / dZ12;
155     Double_t xSecond =
156       ((vGeant3New[3] / vGeant3New[5]) - (vGeant3[3] / vGeant3[5])) / dZ12;
157     Double_t yPrime = (vGeant3New[1] - vGeant3[1]) / dZ12;
158     Double_t ySecond =
159       ((vGeant3New[4] / vGeant3New[5]) - (vGeant3[4] / vGeant3[5])) / dZ12;
160     vGeant3[0] = vGeant3[0] + xPrime * dZ1i - 0.5 * xSecond * dZ1i * dZi2; // X
161     vGeant3[1] = vGeant3[1] + yPrime * dZ1i - 0.5 * ySecond * dZ1i * dZi2; // Y
162     vGeant3[2] = Z; // Z
163     Double_t xPrimeI = xPrime - 0.5 * xSecond * (dZi2 - dZ1i);
164     Double_t yPrimeI = yPrime - 0.5 * ySecond * (dZi2 - dZ1i);
165     // (PX, PY, PZ)/PTOT assuming forward motion
166     vGeant3[5] =
167       1.0 / TMath::Sqrt(1.0 + xPrimeI * xPrimeI + yPrimeI * yPrimeI); // PZ/PTOT
168     vGeant3[3] = xPrimeI * vGeant3[5]; // PX/PTOT
169     vGeant3[4] = yPrimeI * vGeant3[5]; // PY/PTOT
170   } else {
171     AliWarning(Form("Extrap. to Z not reached, Z = %f",Z));    
172   }
173   // Track parameters from Geant3 parameters,
174   // with charge back for forward motion
175   GetFromGeant3Parameters(vGeant3, chargeExtrap * forwardBackward);
176 }
177
178   //__________________________________________________________________________
179 void AliMUONTrackParam::SetGeant3Parameters(Double_t *VGeant3, Double_t ForwardBackward)
180 {
181   // Set vector of Geant3 parameters pointed to by "VGeant3"
182   // from track parameters in current AliMUONTrackParam.
183   // Since AliMUONTrackParam is only geometry, one uses "ForwardBackward"
184   // to know whether the particle is going forward (+1) or backward (-1).
185   VGeant3[0] = this->fNonBendingCoor; // X
186   VGeant3[1] = this->fBendingCoor; // Y
187   VGeant3[2] = this->fZ; // Z
188   Double_t pYZ = TMath::Abs(1.0 / this->fInverseBendingMomentum);
189   Double_t pZ =
190     pYZ / TMath::Sqrt(1.0 + this->fBendingSlope * this->fBendingSlope);
191   VGeant3[6] =
192     TMath::Sqrt(pYZ * pYZ +
193                 pZ * pZ * this->fNonBendingSlope * this->fNonBendingSlope); // PTOT
194   VGeant3[5] = -ForwardBackward * pZ / VGeant3[6]; // PZ/PTOT spectro. z<0
195   VGeant3[3] = this->fNonBendingSlope * VGeant3[5]; // PX/PTOT
196   VGeant3[4] = this->fBendingSlope * VGeant3[5]; // PY/PTOT
197 }
198
199   //__________________________________________________________________________
200 void AliMUONTrackParam::GetFromGeant3Parameters(Double_t *VGeant3, Double_t Charge)
201 {
202   // Get track parameters in current AliMUONTrackParam
203   // from Geant3 parameters pointed to by "VGeant3",
204   // assumed to be calculated for forward motion in Z.
205   // "InverseBendingMomentum" is signed with "Charge".
206   this->fNonBendingCoor = VGeant3[0]; // X
207   this->fBendingCoor = VGeant3[1]; // Y
208   this->fZ = VGeant3[2]; // Z
209   Double_t pYZ = VGeant3[6] * TMath::Sqrt(1.0 - VGeant3[3] * VGeant3[3]);
210   this->fInverseBendingMomentum = Charge / pYZ;
211   this->fBendingSlope = VGeant3[4] / VGeant3[5];
212   this->fNonBendingSlope = VGeant3[3] / VGeant3[5];
213 }
214
215   //__________________________________________________________________________
216 void AliMUONTrackParam::ExtrapToStation(Int_t Station, AliMUONTrackParam *TrackParam)
217 {
218   // Track parameters extrapolated from current track parameters ("this")
219   // to both chambers of the station(0..) "Station"
220   // are returned in the array (dimension 2) of track parameters
221   // pointed to by "TrackParam" (index 0 and 1 for first and second chambers).
222   Double_t extZ[2], z1, z2;
223   Int_t i1 = -1, i2 = -1; // = -1 to avoid compilation warnings
224   // range of Station to be checked ????
225   z1 = AliMUONConstants::DefaultChamberZ(2 * Station);
226   z2 = AliMUONConstants::DefaultChamberZ(2 * Station + 1);
227   // First and second Z to extrapolate at
228   if ((z1 > this->fZ) && (z2 > this->fZ)) {i1 = 0; i2 = 1;}
229   else if ((z1 < this->fZ) && (z2 < this->fZ)) {i1 = 1; i2 = 0;}
230   else {
231         AliError(Form("Starting Z (%f) in between z1 (%f) and z2 (%f) of station(0..)%d",this->fZ,z1,z2,Station));
232 //     cout << "ERROR in AliMUONTrackParam::CreateExtrapSegmentInStation" << endl;
233 //     cout << "Starting Z (" << this->fZ << ") in between z1 (" << z1 <<
234 //       ") and z2 (" << z2 << ") of station(0..) " << Station << endl;
235   }
236   extZ[i1] = z1;
237   extZ[i2] = z2;
238   // copy of track parameters
239   TrackParam[i1] = *this;
240   // first extrapolation
241   (&(TrackParam[i1]))->ExtrapToZ(extZ[0]);
242   TrackParam[i2] = TrackParam[i1];
243   // second extrapolation
244   (&(TrackParam[i2]))->ExtrapToZ(extZ[1]);
245   return;
246 }
247
248   //__________________________________________________________________________
249 void AliMUONTrackParam::ExtrapToVertex(Double_t xVtx, Double_t yVtx, Double_t zVtx)
250 {
251   // Extrapolation to the vertex.
252   // Returns the track parameters resulting from the extrapolation,
253   // in the current TrackParam.
254   // Changes parameters according to Branson correction through the absorber 
255   
256   Double_t zAbsorber = -503.0; // to be coherent with the Geant absorber geometry !!!!
257                                // spectro. (z<0) 
258   // Extrapolates track parameters upstream to the "Z" end of the front absorber
259   ExtrapToZ(zAbsorber); // !!!
260   // Makes Branson correction (multiple scattering + energy loss)
261   BransonCorrection(xVtx,yVtx,zVtx);
262   // Makes a simple magnetic field correction through the absorber
263   FieldCorrection(zAbsorber);
264 }
265
266
267 //  Keep this version for future developments
268   //__________________________________________________________________________
269 // void AliMUONTrackParam::BransonCorrection()
270 // {
271 //   // Branson correction of track parameters
272 //   // the entry parameters have to be calculated at the end of the absorber
273 //   Double_t zEndAbsorber, zBP, xBP, yBP;
274 //   Double_t  pYZ, pX, pY, pZ, pTotal, xEndAbsorber, yEndAbsorber, radiusEndAbsorber2, pT, theta;
275 //   Int_t sign;
276 //   // Would it be possible to calculate all that from Geant configuration ????
277 //   // and to get the Branson parameters from a function in ABSO module ????
278 //   // with an eventual contribution from other detectors like START ????
279 //   // Radiation lengths outer part theta > 3 degres
280 //   static Double_t x01[9] = { 18.8,    // C (cm)
281 //                           10.397,   // Concrete (cm)
282 //                           0.56,    // Plomb (cm)
283 //                           47.26,   // Polyethylene (cm)
284 //                           0.56,   // Plomb (cm)
285 //                           47.26,   // Polyethylene (cm)
286 //                           0.56,   // Plomb (cm)
287 //                           47.26,   // Polyethylene (cm)
288 //                           0.56 };   // Plomb (cm)
289 //   // inner part theta < 3 degres
290 //   static Double_t x02[3] = { 18.8,    // C (cm)
291 //                           10.397,   // Concrete (cm)
292 //                           0.35 };    // W (cm) 
293 //   // z positions of the materials inside the absober outer part theta > 3 degres
294 //   static Double_t z1[10] = { 90, 315, 467, 472, 477, 482, 487, 492, 497, 502 };
295 //   // inner part theta < 3 degres
296 //   static Double_t z2[4] = { 90, 315, 467, 503 };
297 //   static Bool_t first = kTRUE;
298 //   static Double_t zBP1, zBP2, rLimit;
299 //   // Calculates z positions of the Branson's planes: zBP1 for outer part and zBP2 for inner part (only at the first call)
300 //   if (first) {
301 //     first = kFALSE;
302 //     Double_t aNBP = 0.0;
303 //     Double_t aDBP = 0.0;
304 //     Int_t iBound;
305     
306 //     for (iBound = 0; iBound < 9; iBound++) {
307 //       aNBP = aNBP +
308 //      (z1[iBound+1] * z1[iBound+1] * z1[iBound+1] -
309 //       z1[iBound]   * z1[iBound]   * z1[iBound]    ) / x01[iBound];
310 //       aDBP = aDBP +
311 //      (z1[iBound+1] * z1[iBound+1] - z1[iBound]   * z1[iBound]    ) / x01[iBound];
312 //     }
313 //     zBP1 = (2.0 * aNBP) / (3.0 * aDBP);
314 //     aNBP = 0.0;
315 //     aDBP = 0.0;
316 //     for (iBound = 0; iBound < 3; iBound++) {
317 //       aNBP = aNBP +
318 //      (z2[iBound+1] * z2[iBound+1] * z2[iBound+1] -
319 //       z2[iBound]   * z2[iBound ]  * z2[iBound]    ) / x02[iBound];
320 //       aDBP = aDBP +
321 //      (z2[iBound+1] * z2[iBound+1] - z2[iBound] * z2[iBound]) / x02[iBound];
322 //     }
323 //     zBP2 = (2.0 * aNBP) / (3.0 * aDBP);
324 //     rLimit = z2[3] * TMath::Tan(3.0 * (TMath::Pi()) / 180.);
325 //   }
326
327 //   pYZ = TMath::Abs(1.0 / fInverseBendingMomentum);
328 //   sign = 1;      
329 //   if (fInverseBendingMomentum < 0) sign = -1;     
330 //   pZ = pYZ / (TMath::Sqrt(1.0 + fBendingSlope * fBendingSlope)); 
331 //   pX = pZ * fNonBendingSlope; 
332 //   pY = pZ * fBendingSlope; 
333 //   pTotal = TMath::Sqrt(pYZ *pYZ + pX * pX);
334 //   xEndAbsorber = fNonBendingCoor; 
335 //   yEndAbsorber = fBendingCoor; 
336 //   radiusEndAbsorber2 = xEndAbsorber * xEndAbsorber + yEndAbsorber * yEndAbsorber;
337
338 //   if (radiusEndAbsorber2 > rLimit*rLimit) {
339 //     zEndAbsorber = z1[9];
340 //     zBP = zBP1;
341 //   } else {
342 //     zEndAbsorber = z2[3];
343 //     zBP = zBP2;
344 //   }
345
346 //   xBP = xEndAbsorber - (pX / pZ) * (zEndAbsorber - zBP);
347 //   yBP = yEndAbsorber - (pY / pZ) * (zEndAbsorber - zBP);
348
349 //   // new parameters after Branson and energy loss corrections
350 //   pZ = pTotal * zBP / TMath::Sqrt(xBP * xBP + yBP * yBP + zBP * zBP);
351 //   pX = pZ * xBP / zBP;
352 //   pY = pZ * yBP / zBP;
353 //   fBendingSlope = pY / pZ;
354 //   fNonBendingSlope = pX / pZ;
355   
356 //   pT = TMath::Sqrt(pX * pX + pY * pY);      
357 //   theta = TMath::ATan2(pT, pZ); 
358 //   pTotal =
359 //     TotalMomentumEnergyLoss(rLimit, pTotal, theta, xEndAbsorber, yEndAbsorber);
360
361 //   fInverseBendingMomentum = (sign / pTotal) *
362 //     TMath::Sqrt(1.0 +
363 //              fBendingSlope * fBendingSlope +
364 //              fNonBendingSlope * fNonBendingSlope) /
365 //     TMath::Sqrt(1.0 + fBendingSlope * fBendingSlope);
366
367 //   // vertex position at (0,0,0)
368 //   // should be taken from vertex measurement ???
369 //   fBendingCoor = 0.0;
370 //   fNonBendingCoor = 0;
371 //   fZ= 0;
372 // }
373
374 void AliMUONTrackParam::BransonCorrection(Double_t xVtx,Double_t yVtx,Double_t zVtx)
375 {
376   // Branson correction of track parameters
377   // the entry parameters have to be calculated at the end of the absorber
378   // simplified version: the z positions of Branson's planes are no longer calculated
379   // but are given as inputs. One can use the macros MUONTestAbso.C and DrawTestAbso.C
380   // to test this correction. 
381   // Would it be possible to calculate all that from Geant configuration ????
382   // and to get the Branson parameters from a function in ABSO module ????
383   // with an eventual contribution from other detectors like START ????
384   //change to take into account the vertex postition (real, reconstruct,....)
385
386   Double_t  zBP, xBP, yBP;
387   Double_t  pYZ, pX, pY, pZ, pTotal, xEndAbsorber, yEndAbsorber, radiusEndAbsorber2, pT, theta;
388   Int_t sign;
389   static Bool_t first = kTRUE;
390   static Double_t zBP1, zBP2, rLimit, thetaLimit, zEndAbsorber;
391   // zBP1 for outer part and zBP2 for inner part (only at the first call)
392   if (first) {
393     first = kFALSE;
394   
395     zEndAbsorber = -503;  // spectro (z<0)
396     thetaLimit = 3.0 * (TMath::Pi()) / 180.;
397     rLimit = TMath::Abs(zEndAbsorber) * TMath::Tan(thetaLimit);
398     zBP1 = -450; // values close to those calculated with EvalAbso.C
399     zBP2 = -480;
400   }
401
402   pYZ = TMath::Abs(1.0 / fInverseBendingMomentum);
403   sign = 1;      
404   if (fInverseBendingMomentum < 0) sign = -1;  
405   pZ = Pz();
406   pX = Px(); 
407   pY = Py(); 
408   pTotal = TMath::Sqrt(pYZ *pYZ + pX * pX);
409   xEndAbsorber = fNonBendingCoor; 
410   yEndAbsorber = fBendingCoor; 
411   radiusEndAbsorber2 = xEndAbsorber * xEndAbsorber + yEndAbsorber * yEndAbsorber;
412
413   if (radiusEndAbsorber2 > rLimit*rLimit) {
414     zBP = zBP1;
415   } else {
416     zBP = zBP2;
417   }
418
419   xBP = xEndAbsorber - (pX / pZ) * (zEndAbsorber - zBP);
420   yBP = yEndAbsorber - (pY / pZ) * (zEndAbsorber - zBP);
421
422   // new parameters after Branson and energy loss corrections
423 //   Float_t zSmear = zBP - gRandom->Gaus(0.,2.);  // !!! possible smearing of Z vertex position
424
425   Float_t zSmear = zBP ;
426   
427    pZ = pTotal * (zSmear-zVtx) / TMath::Sqrt((xBP-xVtx) * (xBP-xVtx) + (yBP-yVtx) * (yBP-yVtx) +( zSmear-zVtx) * (zSmear-zVtx) );
428    pX = pZ * (xBP - xVtx)/ (zSmear-zVtx);
429    pY = pZ * (yBP - yVtx) / (zSmear-zVtx);
430   fBendingSlope = pY / pZ;
431   fNonBendingSlope = pX / pZ;
432
433   
434   pT = TMath::Sqrt(pX * pX + pY * pY);      
435   theta = TMath::ATan2(pT, TMath::Abs(pZ)); 
436   pTotal = TotalMomentumEnergyLoss(thetaLimit, pTotal, theta);
437
438   fInverseBendingMomentum = (sign / pTotal) *
439     TMath::Sqrt(1.0 +
440                 fBendingSlope * fBendingSlope +
441                 fNonBendingSlope * fNonBendingSlope) /
442     TMath::Sqrt(1.0 + fBendingSlope * fBendingSlope);
443
444   // vertex position at (0,0,0)
445   // should be taken from vertex measurement ???
446
447   fBendingCoor = xVtx;
448   fNonBendingCoor = yVtx;
449   fZ= zVtx;
450
451 }
452
453   //__________________________________________________________________________
454 Double_t AliMUONTrackParam::TotalMomentumEnergyLoss(Double_t thetaLimit, Double_t pTotal, Double_t theta)
455 {
456   // Returns the total momentum corrected from energy loss in the front absorber
457   // One can use the macros MUONTestAbso.C and DrawTestAbso.C
458   // to test this correction. 
459   // Momentum energy loss behaviour evaluated with the simulation of single muons (april 2002)
460   Double_t deltaP, pTotalCorrected;
461
462    // Parametrization to be redone according to change of absorber material ????
463   // See remark in function BransonCorrection !!!!
464   // The name is not so good, and there are many arguments !!!!
465   if (theta  < thetaLimit ) {
466     if (pTotal < 20) {
467       deltaP = 2.5938 + 0.0570 * pTotal - 0.001151 * pTotal * pTotal;
468     } else {
469       deltaP = 3.0714 + 0.011767 *pTotal;
470     }
471     deltaP *= 0.75; // AZ
472   } else {
473     if (pTotal < 20) {
474       deltaP  = 2.1207 + 0.05478 * pTotal - 0.00145079 * pTotal * pTotal;
475     } else { 
476       deltaP = 2.6069 + 0.0051705 * pTotal;
477     }
478     deltaP *= 0.9; // AZ
479   }
480   pTotalCorrected = pTotal + deltaP / TMath::Cos(theta);
481   return pTotalCorrected;
482 }
483
484   //__________________________________________________________________________
485 void AliMUONTrackParam::FieldCorrection(Double_t Z)
486 {
487   // 
488   // Correction of the effect of the magnetic field in the absorber
489   // Assume a constant field along Z axis.
490
491   Float_t b[3],x[3]; 
492   Double_t bZ;
493   Double_t pYZ,pX,pY,pZ,pT;
494   Double_t pXNew,pYNew;
495   Double_t c;
496
497   pYZ = TMath::Abs(1.0 / fInverseBendingMomentum);
498   c = TMath::Sign(1.0,fInverseBendingMomentum); // particle charge 
499  
500   pZ = Pz();
501   pX = Px(); 
502   pY = Py();
503   pT = TMath::Sqrt(pX*pX+pY*pY);
504
505   if (TMath::Abs(pZ) <= 0) return;
506   x[2] = Z/2;
507   x[0] = x[2]*fNonBendingSlope;  
508   x[1] = x[2]*fBendingSlope;
509
510   // Take magn. field value at position x.
511   fkField->Field(x, b);
512   bZ =  b[2];
513  
514   // Transverse momentum rotation
515   // Parameterized with the study of DeltaPhi = phiReco - phiGen as a function of pZ.
516   Double_t phiShift = c*0.436*0.0003*bZ*Z/pZ; 
517  // Rotate momentum around Z axis.
518   pXNew = pX*TMath::Cos(phiShift) - pY*TMath::Sin(phiShift);
519   pYNew = pX*TMath::Sin(phiShift) + pY*TMath::Cos(phiShift);
520  
521   fBendingSlope = pYNew / pZ;
522   fNonBendingSlope = pXNew / pZ;
523   
524   fInverseBendingMomentum = c / TMath::Sqrt(pYNew*pYNew+pZ*pZ);
525  
526 }
527   //__________________________________________________________________________
528 Double_t AliMUONTrackParam::Px() const
529 {
530   // return px from track paramaters
531   Double_t pYZ, pZ, pX;
532   pYZ = 0;
533   if (  TMath::Abs(fInverseBendingMomentum) > 0 )
534     pYZ = TMath::Abs(1.0 / fInverseBendingMomentum);
535   pZ = -pYZ / (TMath::Sqrt(1.0 + fBendingSlope * fBendingSlope));  // spectro. (z<0)
536   pX = pZ * fNonBendingSlope; 
537   return pX;
538 }
539   //__________________________________________________________________________
540 Double_t AliMUONTrackParam::Py() const
541 {
542   // return px from track paramaters
543   Double_t pYZ, pZ, pY;
544   pYZ = 0;
545   if (  TMath::Abs(fInverseBendingMomentum) > 0 )
546     pYZ = TMath::Abs(1.0 / fInverseBendingMomentum);
547   pZ = -pYZ / (TMath::Sqrt(1.0 + fBendingSlope * fBendingSlope));  // spectro. (z<0)
548   pY = pZ * fBendingSlope; 
549   return pY;
550 }
551   //__________________________________________________________________________
552 Double_t AliMUONTrackParam::Pz() const
553 {
554   // return px from track paramaters
555   Double_t pYZ, pZ;
556   pYZ = 0;
557   if (  TMath::Abs(fInverseBendingMomentum) > 0 )
558     pYZ = TMath::Abs(1.0 / fInverseBendingMomentum);
559   pZ = -pYZ / (TMath::Sqrt(1.0 + fBendingSlope * fBendingSlope));  // spectro. (z<0)
560   return pZ;
561 }
562   //__________________________________________________________________________
563 Double_t AliMUONTrackParam::P() const
564 {
565   // return p from track paramaters
566   Double_t  pYZ, pZ, p;
567   pYZ = 0;
568   if (  TMath::Abs(fInverseBendingMomentum) > 0 )
569     pYZ = TMath::Abs(1.0 / fInverseBendingMomentum);
570   pZ = -pYZ / (TMath::Sqrt(1.0 + fBendingSlope * fBendingSlope));  // spectro. (z<0)
571   p = TMath::Abs(pZ) * 
572     TMath::Sqrt(1.0 + fBendingSlope * fBendingSlope + fNonBendingSlope * fNonBendingSlope);
573   return p;
574   
575 }
576  //__________________________________________________________________________
577 void AliMUONTrackParam::ExtrapOneStepHelix(Double_t charge, Double_t step, 
578                                          Double_t *vect, Double_t *vout) const
579 {
580 //    ******************************************************************
581 //    *                                                                *
582 //    *  Performs the tracking of one step in a magnetic field         *
583 //    *  The trajectory is assumed to be a helix in a constant field   *
584 //    *  taken at the mid point of the step.                           *
585 //    *  Parameters:                                                   *
586 //    *   input                                                        *
587 //    *     STEP =arc length of the step asked                         *
588 //    *     VECT =input vector (position,direction cos and momentum)   *
589 //    *     CHARGE=  electric charge of the particle                   *
590 //    *   output                                                       *
591 //    *     VOUT = same as VECT after completion of the step           *
592 //    *                                                                *
593 //    *    ==>Called by : <USER>, GUSWIM                               *
594 //    *       Author    m.hansroul  *********                          *
595 //    *       modified  s.egli, s.v.levonian                           *
596 //    *       modified  v.perevoztchikov
597 //    *                                                                *
598 //    ******************************************************************
599 //
600
601 // modif: everything in double precision
602
603     Double_t xyz[3], h[4], hxp[3];
604     Double_t h2xy, hp, rho, tet;
605     Double_t sint, sintt, tsint, cos1t;
606     Double_t f1, f2, f3, f4, f5, f6;
607
608     const Int_t kix  = 0;
609     const Int_t kiy  = 1;
610     const Int_t kiz  = 2;
611     const Int_t kipx = 3;
612     const Int_t kipy = 4;
613     const Int_t kipz = 5;
614     const Int_t kipp = 6;
615
616     const Double_t kec = 2.9979251e-4;
617     //
618     //    ------------------------------------------------------------------
619     //
620     //       units are kgauss,centimeters,gev/c
621     //
622     vout[kipp] = vect[kipp];
623     if (TMath::Abs(charge) < 0.00001) {
624       for (Int_t i = 0; i < 3; i++) {
625         vout[i] = vect[i] + step * vect[i+3];
626         vout[i+3] = vect[i+3];
627       }
628       return;
629     }
630     xyz[0]    = vect[kix] + 0.5 * step * vect[kipx];
631     xyz[1]    = vect[kiy] + 0.5 * step * vect[kipy];
632     xyz[2]    = vect[kiz] + 0.5 * step * vect[kipz];
633
634     //cmodif: call gufld (xyz, h) changed into:
635     GetField (xyz, h);
636  
637     h2xy = h[0]*h[0] + h[1]*h[1];
638     h[3] = h[2]*h[2]+ h2xy;
639     if (h[3] < 1.e-12) {
640       for (Int_t i = 0; i < 3; i++) {
641         vout[i] = vect[i] + step * vect[i+3];
642         vout[i+3] = vect[i+3];
643       }
644       return;
645     }
646     if (h2xy < 1.e-12*h[3]) {
647       ExtrapOneStepHelix3(charge*h[2], step, vect, vout);
648       return;
649     }
650     h[3] = TMath::Sqrt(h[3]);
651     h[0] /= h[3];
652     h[1] /= h[3];
653     h[2] /= h[3];
654     h[3] *= kec;
655
656     hxp[0] = h[1]*vect[kipz] - h[2]*vect[kipy];
657     hxp[1] = h[2]*vect[kipx] - h[0]*vect[kipz];
658     hxp[2] = h[0]*vect[kipy] - h[1]*vect[kipx];
659  
660     hp = h[0]*vect[kipx] + h[1]*vect[kipy] + h[2]*vect[kipz];
661
662     rho = -charge*h[3]/vect[kipp];
663     tet = rho * step;
664
665     if (TMath::Abs(tet) > 0.15) {
666       sint = TMath::Sin(tet);
667       sintt = (sint/tet);
668       tsint = (tet-sint)/tet;
669       cos1t = 2.*(TMath::Sin(0.5*tet))*(TMath::Sin(0.5*tet))/tet;
670     } else {
671       tsint = tet*tet/36.;
672       sintt = (1. - tsint);
673       sint = tet*sintt;
674       cos1t = 0.5*tet;
675     }
676
677     f1 = step * sintt;
678     f2 = step * cos1t;
679     f3 = step * tsint * hp;
680     f4 = -tet*cos1t;
681     f5 = sint;
682     f6 = tet * cos1t * hp;
683  
684     vout[kix] = vect[kix] + f1*vect[kipx] + f2*hxp[0] + f3*h[0];
685     vout[kiy] = vect[kiy] + f1*vect[kipy] + f2*hxp[1] + f3*h[1];
686     vout[kiz] = vect[kiz] + f1*vect[kipz] + f2*hxp[2] + f3*h[2];
687  
688     vout[kipx] = vect[kipx] + f4*vect[kipx] + f5*hxp[0] + f6*h[0];
689     vout[kipy] = vect[kipy] + f4*vect[kipy] + f5*hxp[1] + f6*h[1];
690     vout[kipz] = vect[kipz] + f4*vect[kipz] + f5*hxp[2] + f6*h[2];
691  
692     return;
693 }
694
695  //__________________________________________________________________________
696 void AliMUONTrackParam::ExtrapOneStepHelix3(Double_t field, Double_t step, 
697                                                Double_t *vect, Double_t *vout) const
698 {
699 // 
700 //     ******************************************************************
701 //     *                                                                *
702 //     *       Tracking routine in a constant field oriented            *
703 //     *       along axis 3                                             *
704 //     *       Tracking is performed with a conventional                *
705 //     *       helix step method                                        *
706 //     *                                                                *
707 //     *    ==>Called by : <USER>, GUSWIM                               *
708 //     *       Authors    R.Brun, M.Hansroul  *********                 *
709 //     *       Rewritten  V.Perevoztchikov
710 //     *                                                                *
711 //     ******************************************************************
712 // 
713
714     Double_t hxp[3];
715     Double_t h4, hp, rho, tet;
716     Double_t sint, sintt, tsint, cos1t;
717     Double_t f1, f2, f3, f4, f5, f6;
718
719     const Int_t kix  = 0;
720     const Int_t kiy  = 1;
721     const Int_t kiz  = 2;
722     const Int_t kipx = 3;
723     const Int_t kipy = 4;
724     const Int_t kipz = 5;
725     const Int_t kipp = 6;
726
727     const Double_t kec = 2.9979251e-4;
728
729 // 
730 //     ------------------------------------------------------------------
731 // 
732 //       units are kgauss,centimeters,gev/c
733 // 
734     vout[kipp] = vect[kipp];
735     h4 = field * kec;
736
737     hxp[0] = - vect[kipy];
738     hxp[1] = + vect[kipx];
739  
740     hp = vect[kipz];
741
742     rho = -h4/vect[kipp];
743     tet = rho * step;
744     if (TMath::Abs(tet) > 0.15) {
745       sint = TMath::Sin(tet);
746       sintt = (sint/tet);
747       tsint = (tet-sint)/tet;
748       cos1t = 2.* TMath::Sin(0.5*tet) * TMath::Sin(0.5*tet)/tet;
749     } else {
750       tsint = tet*tet/36.;
751       sintt = (1. - tsint);
752       sint = tet*sintt;
753       cos1t = 0.5*tet;
754     }
755
756     f1 = step * sintt;
757     f2 = step * cos1t;
758     f3 = step * tsint * hp;
759     f4 = -tet*cos1t;
760     f5 = sint;
761     f6 = tet * cos1t * hp;
762  
763     vout[kix] = vect[kix] + f1*vect[kipx] + f2*hxp[0];
764     vout[kiy] = vect[kiy] + f1*vect[kipy] + f2*hxp[1];
765     vout[kiz] = vect[kiz] + f1*vect[kipz] + f3;
766  
767     vout[kipx] = vect[kipx] + f4*vect[kipx] + f5*hxp[0];
768     vout[kipy] = vect[kipy] + f4*vect[kipy] + f5*hxp[1];
769     vout[kipz] = vect[kipz] + f4*vect[kipz] + f6;
770
771     return;
772 }
773  //__________________________________________________________________________
774 void AliMUONTrackParam::ExtrapOneStepRungekutta(Double_t charge, Double_t step, 
775                                                      Double_t* vect, Double_t* vout) const
776 {
777 // 
778 //     ******************************************************************
779 //     *                                                                *
780 //     *  Runge-Kutta method for tracking a particle through a magnetic *
781 //     *  field. Uses Nystroem algorithm (See Handbook Nat. Bur. of     *
782 //     *  Standards, procedure 25.5.20)                                 *
783 //     *                                                                *
784 //     *  Input parameters                                              *
785 //     *       CHARGE    Particle charge                                *
786 //     *       STEP      Step size                                      *
787 //     *       VECT      Initial co-ords,direction cosines,momentum     *
788 //     *  Output parameters                                             *
789 //     *       VOUT      Output co-ords,direction cosines,momentum      *
790 //     *  User routine called                                           *
791 //     *       CALL GUFLD(X,F)                                          *
792 //     *                                                                *
793 //     *    ==>Called by : <USER>, GUSWIM                               *
794 //     *       Authors    R.Brun, M.Hansroul  *********                 *
795 //     *                  V.Perevoztchikov (CUT STEP implementation)    *
796 //     *                                                                *
797 //     *                                                                *
798 //     ******************************************************************
799 // 
800
801     Double_t h2, h4, f[4];
802     Double_t xyzt[3], a, b, c, ph,ph2;
803     Double_t secxs[4],secys[4],seczs[4],hxp[3];
804     Double_t g1, g2, g3, g4, g5, g6, ang2, dxt, dyt, dzt;
805     Double_t est, at, bt, ct, cba;
806     Double_t f1, f2, f3, f4, rho, tet, hnorm, hp, rho1, sint, cost;
807     
808     Double_t x;
809     Double_t y;
810     Double_t z;
811     
812     Double_t xt;
813     Double_t yt;
814     Double_t zt;
815
816     Double_t maxit = 1992;
817     Double_t maxcut = 11;
818
819     const Double_t kdlt   = 1e-4;
820     const Double_t kdlt32 = kdlt/32.;
821     const Double_t kthird = 1./3.;
822     const Double_t khalf  = 0.5;
823     const Double_t kec = 2.9979251e-4;
824
825     const Double_t kpisqua = 9.86960440109;
826     const Int_t kix  = 0;
827     const Int_t kiy  = 1;
828     const Int_t kiz  = 2;
829     const Int_t kipx = 3;
830     const Int_t kipy = 4;
831     const Int_t kipz = 5;
832   
833     // *.
834     // *.    ------------------------------------------------------------------
835     // *.
836     // *             this constant is for units cm,gev/c and kgauss
837     // *
838     Int_t iter = 0;
839     Int_t ncut = 0;
840     for(Int_t j = 0; j < 7; j++)
841       vout[j] = vect[j];
842
843     Double_t  pinv   = kec * charge / vect[6];
844     Double_t tl = 0.;
845     Double_t h = step;
846     Double_t rest;
847
848  
849     do {
850       rest  = step - tl;
851       if (TMath::Abs(h) > TMath::Abs(rest)) h = rest;
852       //cmodif: call gufld(vout,f) changed into:
853
854       GetField(vout,f);
855
856       // *
857       // *             start of integration
858       // *
859       x      = vout[0];
860       y      = vout[1];
861       z      = vout[2];
862       a      = vout[3];
863       b      = vout[4];
864       c      = vout[5];
865
866       h2     = khalf * h;
867       h4     = khalf * h2;
868       ph     = pinv * h;
869       ph2    = khalf * ph;
870       secxs[0] = (b * f[2] - c * f[1]) * ph2;
871       secys[0] = (c * f[0] - a * f[2]) * ph2;
872       seczs[0] = (a * f[1] - b * f[0]) * ph2;
873       ang2 = (secxs[0]*secxs[0] + secys[0]*secys[0] + seczs[0]*seczs[0]);
874       if (ang2 > kpisqua) break;
875
876       dxt    = h2 * a + h4 * secxs[0];
877       dyt    = h2 * b + h4 * secys[0];
878       dzt    = h2 * c + h4 * seczs[0];
879       xt     = x + dxt;
880       yt     = y + dyt;
881       zt     = z + dzt;
882       // *
883       // *              second intermediate point
884       // *
885
886       est = TMath::Abs(dxt) + TMath::Abs(dyt) + TMath::Abs(dzt);
887       if (est > h) {
888         if (ncut++ > maxcut) break;
889         h *= khalf;
890         continue;
891       }
892  
893       xyzt[0] = xt;
894       xyzt[1] = yt;
895       xyzt[2] = zt;
896
897       //cmodif: call gufld(xyzt,f) changed into:
898       GetField(xyzt,f);
899
900       at     = a + secxs[0];
901       bt     = b + secys[0];
902       ct     = c + seczs[0];
903
904       secxs[1] = (bt * f[2] - ct * f[1]) * ph2;
905       secys[1] = (ct * f[0] - at * f[2]) * ph2;
906       seczs[1] = (at * f[1] - bt * f[0]) * ph2;
907       at     = a + secxs[1];
908       bt     = b + secys[1];
909       ct     = c + seczs[1];
910       secxs[2] = (bt * f[2] - ct * f[1]) * ph2;
911       secys[2] = (ct * f[0] - at * f[2]) * ph2;
912       seczs[2] = (at * f[1] - bt * f[0]) * ph2;
913       dxt    = h * (a + secxs[2]);
914       dyt    = h * (b + secys[2]);
915       dzt    = h * (c + seczs[2]);
916       xt     = x + dxt;
917       yt     = y + dyt;
918       zt     = z + dzt;
919       at     = a + 2.*secxs[2];
920       bt     = b + 2.*secys[2];
921       ct     = c + 2.*seczs[2];
922
923       est = TMath::Abs(dxt)+TMath::Abs(dyt)+TMath::Abs(dzt);
924       if (est > 2.*TMath::Abs(h)) {
925         if (ncut++ > maxcut) break;
926         h *= khalf;
927         continue;
928       }
929  
930       xyzt[0] = xt;
931       xyzt[1] = yt;
932       xyzt[2] = zt;
933
934       //cmodif: call gufld(xyzt,f) changed into:
935       GetField(xyzt,f);
936
937       z      = z + (c + (seczs[0] + seczs[1] + seczs[2]) * kthird) * h;
938       y      = y + (b + (secys[0] + secys[1] + secys[2]) * kthird) * h;
939       x      = x + (a + (secxs[0] + secxs[1] + secxs[2]) * kthird) * h;
940
941       secxs[3] = (bt*f[2] - ct*f[1])* ph2;
942       secys[3] = (ct*f[0] - at*f[2])* ph2;
943       seczs[3] = (at*f[1] - bt*f[0])* ph2;
944       a      = a+(secxs[0]+secxs[3]+2. * (secxs[1]+secxs[2])) * kthird;
945       b      = b+(secys[0]+secys[3]+2. * (secys[1]+secys[2])) * kthird;
946       c      = c+(seczs[0]+seczs[3]+2. * (seczs[1]+seczs[2])) * kthird;
947
948       est    = TMath::Abs(secxs[0]+secxs[3] - (secxs[1]+secxs[2]))
949         + TMath::Abs(secys[0]+secys[3] - (secys[1]+secys[2]))
950         + TMath::Abs(seczs[0]+seczs[3] - (seczs[1]+seczs[2]));
951
952       if (est > kdlt && TMath::Abs(h) > 1.e-4) {
953         if (ncut++ > maxcut) break;
954         h *= khalf;
955         continue;
956       }
957
958       ncut = 0;
959       // *               if too many iterations, go to helix
960       if (iter++ > maxit) break;
961
962       tl += h;
963       if (est < kdlt32) 
964         h *= 2.;
965       cba    = 1./ TMath::Sqrt(a*a + b*b + c*c);
966       vout[0] = x;
967       vout[1] = y;
968       vout[2] = z;
969       vout[3] = cba*a;
970       vout[4] = cba*b;
971       vout[5] = cba*c;
972       rest = step - tl;
973       if (step < 0.) rest = -rest;
974       if (rest < 1.e-5*TMath::Abs(step)) return;
975
976     } while(1);
977
978     // angle too big, use helix
979
980     f1  = f[0];
981     f2  = f[1];
982     f3  = f[2];
983     f4  = TMath::Sqrt(f1*f1+f2*f2+f3*f3);
984     rho = -f4*pinv;
985     tet = rho * step;
986  
987     hnorm = 1./f4;
988     f1 = f1*hnorm;
989     f2 = f2*hnorm;
990     f3 = f3*hnorm;
991
992     hxp[0] = f2*vect[kipz] - f3*vect[kipy];
993     hxp[1] = f3*vect[kipx] - f1*vect[kipz];
994     hxp[2] = f1*vect[kipy] - f2*vect[kipx];
995  
996     hp = f1*vect[kipx] + f2*vect[kipy] + f3*vect[kipz];
997
998     rho1 = 1./rho;
999     sint = TMath::Sin(tet);
1000     cost = 2.*TMath::Sin(khalf*tet)*TMath::Sin(khalf*tet);
1001
1002     g1 = sint*rho1;
1003     g2 = cost*rho1;
1004     g3 = (tet-sint) * hp*rho1;
1005     g4 = -cost;
1006     g5 = sint;
1007     g6 = cost * hp;
1008  
1009     vout[kix] = vect[kix] + g1*vect[kipx] + g2*hxp[0] + g3*f1;
1010     vout[kiy] = vect[kiy] + g1*vect[kipy] + g2*hxp[1] + g3*f2;
1011     vout[kiz] = vect[kiz] + g1*vect[kipz] + g2*hxp[2] + g3*f3;
1012  
1013     vout[kipx] = vect[kipx] + g4*vect[kipx] + g5*hxp[0] + g6*f1;
1014     vout[kipy] = vect[kipy] + g4*vect[kipy] + g5*hxp[1] + g6*f2;
1015     vout[kipz] = vect[kipz] + g4*vect[kipz] + g5*hxp[2] + g6*f3;
1016
1017     return;
1018 }
1019 //___________________________________________________________
1020  void  AliMUONTrackParam::GetField(Double_t *Position, Double_t *Field) const
1021 {
1022     // interface for arguments in double precision (Why ? ChF)
1023
1024     Float_t x[3], b[3];
1025
1026     x[0] = Position[0]; x[1] = Position[1]; x[2] = Position[2];
1027
1028     fkField->Field(x, b);
1029     Field[0] = b[0]; Field[1] = b[1]; Field[2] = b[2];
1030
1031     return;
1032   }
1033 //_____________________________________________-
1034 void AliMUONTrackParam::Print(Option_t* opt) const
1035 {
1036 //
1037   // Printing TrackParam information 
1038   // "full" option for printing all the information about the TrackParam
1039   //
1040   TString sopt(opt);
1041   sopt.ToUpper();
1042  
1043   if ( sopt.Contains("FULL") ) { 
1044     cout << "<AliMUONTrackParam> Bending P=" << setw(5) << setprecision(3)  << 1./GetInverseBendingMomentum() << 
1045       ", NonBendSlope=" << setw(5) << setprecision(3)  << GetNonBendingSlope()*180./TMath::Pi() <<
1046       ", BendSlope=" << setw(5) << setprecision(3)     << GetBendingSlope()*180./TMath::Pi()  << 
1047       ", (x,y,z)_IP=(" <<  setw(5) << setprecision(3) << GetNonBendingCoor() <<
1048       "," <<  setw(5) << setprecision(3) << GetBendingCoor() <<
1049       "," <<  setw(5) << setprecision(3) << GetZ() <<
1050       ") cm, (px,py,pz)=(" << setw(5) << setprecision(3) << Px() <<
1051       "," << setw(5) << setprecision(3) << Py() <<
1052       "," << setw(5) << setprecision(3) << Pz() << ") GeV/c" << endl;
1053   }
1054   else {
1055     cout << "<AliMUONTrackParam>"  << endl;
1056   }
1057     
1058 }