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