]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONTrackExtrap.cxx
A new function for calculating the mean crossed material (A.Dainese, A.Gheata)
[u/mrichter/AliRoot.git] / MUON / AliMUONTrackExtrap.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //-----------------------------------------------------------------------------
19 // Class AliMUONTrackExtrap
20 // ------------------------
21 // Tools for track extrapolation in ALICE dimuon spectrometer
22 // Author: Philippe Pillot
23 //-----------------------------------------------------------------------------
24
25 #include "AliMUONTrackExtrap.h" 
26 #include "AliMUONTrackParam.h"
27 #include "AliMUONConstants.h"
28
29 #include "AliMagF.h" 
30
31 #include <TMath.h>
32 #include <TMatrixD.h>
33 #include <TGeoManager.h>
34
35 #include <Riostream.h>
36
37 /// \cond CLASSIMP
38 ClassImp(AliMUONTrackExtrap) // Class implementation in ROOT context
39 /// \endcond
40
41 const AliMagF* AliMUONTrackExtrap::fgkField = 0x0;
42 const Bool_t   AliMUONTrackExtrap::fgkUseHelix = kFALSE;
43 const Int_t    AliMUONTrackExtrap::fgkMaxStepNumber = 5000;
44 const Double_t AliMUONTrackExtrap::fgkHelixStepLength = 6.;
45 const Double_t AliMUONTrackExtrap::fgkRungeKuttaMaxResidue = 0.002;
46
47   //__________________________________________________________________________
48 Double_t AliMUONTrackExtrap::GetImpactParamFromBendingMomentum(Double_t bendingMomentum)
49 {
50   /// Returns impact parameter at vertex in bending plane (cm),
51   /// from the signed bending momentum "BendingMomentum" in bending plane (GeV/c),
52   /// using simple values for dipole magnetic field.
53   /// The sign of "BendingMomentum" is the sign of the charge.
54   
55   if (bendingMomentum == 0.) return 1.e10;
56   
57   Double_t simpleBPosition = 0.5 * (AliMUONConstants::CoilZ() + AliMUONConstants::YokeZ());
58   Double_t simpleBLength = 0.5 * (AliMUONConstants::CoilL() + AliMUONConstants::YokeL());
59   Float_t b[3], x[3] = {0.,0.,(Float_t) simpleBPosition};
60   if (fgkField) fgkField->Field(x,b);
61   else {
62     cout<<"F-AliMUONTrackExtrap::GetField: fgkField = 0x0"<<endl;
63     exit(-1);
64   }
65   Double_t simpleBValue = (Double_t) b[0];
66   
67   return (-0.0003 * simpleBValue * simpleBLength * simpleBPosition / bendingMomentum);
68 }
69
70   //__________________________________________________________________________
71 Double_t AliMUONTrackExtrap::GetBendingMomentumFromImpactParam(Double_t impactParam)
72 {
73   /// Returns signed bending momentum in bending plane (GeV/c),
74   /// the sign being the sign of the charge for particles moving forward in Z,
75   /// from the impact parameter "ImpactParam" at vertex in bending plane (cm),
76   /// using simple values for dipole magnetic field.
77   
78   if (impactParam == 0.) return 1.e10;
79   
80   Double_t simpleBPosition = 0.5 * (AliMUONConstants::CoilZ() + AliMUONConstants::YokeZ());
81   Double_t simpleBLength = 0.5 * (AliMUONConstants::CoilL() + AliMUONConstants::YokeL());
82   Float_t b[3], x[3] = {0.,0.,(Float_t) simpleBPosition};
83   if (fgkField) fgkField->Field(x,b);
84   else {
85     cout<<"F-AliMUONTrackExtrap::GetField: fgkField = 0x0"<<endl;
86     exit(-1);
87   }
88   Double_t simpleBValue = (Double_t) b[0];
89   
90   return (-0.0003 * simpleBValue * simpleBLength * simpleBPosition / impactParam);
91 }
92
93   //__________________________________________________________________________
94 void AliMUONTrackExtrap::LinearExtrapToZ(AliMUONTrackParam* trackParam, Double_t zEnd)
95 {
96   /// Track parameters (and their covariances if any) linearly extrapolated to the plane at "zEnd".
97   /// On return, results from the extrapolation are updated in trackParam.
98   
99   if (trackParam->GetZ() == zEnd) return; // nothing to be done if same z
100   
101   // Compute track parameters
102   Double_t dZ = zEnd - trackParam->GetZ();
103   trackParam->SetNonBendingCoor(trackParam->GetNonBendingCoor() + trackParam->GetNonBendingSlope() * dZ);
104   trackParam->SetBendingCoor(trackParam->GetBendingCoor() + trackParam->GetBendingSlope() * dZ);
105   trackParam->SetZ(zEnd);
106   
107   // Update track parameters covariances if any
108   if (trackParam->CovariancesExist()) {
109     TMatrixD paramCov(trackParam->GetCovariances());
110     paramCov(0,0) += dZ * dZ * paramCov(1,1) + 2. * dZ * paramCov(0,1);
111     paramCov(0,1) += dZ * paramCov(1,1);
112     paramCov(1,0) = paramCov(0,1);
113     paramCov(2,2) += dZ * dZ * paramCov(3,3) + 2. * dZ * paramCov(2,3);
114     paramCov(2,3) += dZ * paramCov(3,3);
115     paramCov(3,2) = paramCov(2,3);
116     trackParam->SetCovariances(paramCov);
117   }
118   
119 }
120
121   //__________________________________________________________________________
122 void AliMUONTrackExtrap::ExtrapToZ(AliMUONTrackParam* trackParam, Double_t zEnd)
123 {
124   /// Interface to track parameter extrapolation to the plane at "Z" using Helix or Rungekutta algorithm.
125   /// On return, the track parameters resulting from the extrapolation are updated in trackParam.
126   if (fgkUseHelix) AliMUONTrackExtrap::ExtrapToZHelix(trackParam,zEnd);
127   else AliMUONTrackExtrap::ExtrapToZRungekutta(trackParam,zEnd);
128 }
129
130   //__________________________________________________________________________
131 void AliMUONTrackExtrap::ExtrapToZHelix(AliMUONTrackParam* trackParam, Double_t zEnd)
132 {
133   /// Track parameter extrapolation to the plane at "Z" using Helix algorithm.
134   /// On return, the track parameters resulting from the extrapolation are updated in trackParam.
135   if (trackParam->GetZ() == zEnd) return; // nothing to be done if same Z
136   Double_t forwardBackward; // +1 if forward, -1 if backward
137   if (zEnd < trackParam->GetZ()) forwardBackward = 1.0; // spectro. z<0 
138   else forwardBackward = -1.0;
139   Double_t v3[7], v3New[7]; // 7 in parameter ????
140   Int_t i3, stepNumber;
141   // For safety: return kTRUE or kFALSE ????
142   // Parameter vector for calling EXTRAP_ONESTEP
143   ConvertTrackParamForExtrap(trackParam, forwardBackward, v3);
144   // sign of charge (sign of fInverseBendingMomentum if forward motion)
145   // must be changed if backward extrapolation
146   Double_t chargeExtrap = forwardBackward * TMath::Sign(Double_t(1.0), trackParam->GetInverseBendingMomentum());
147   // Extrapolation loop
148   stepNumber = 0;
149   while (((-forwardBackward * (v3[2] - zEnd)) <= 0.0) && (stepNumber < fgkMaxStepNumber)) { // spectro. z<0
150     stepNumber++;
151     ExtrapOneStepHelix(chargeExtrap, fgkHelixStepLength, v3, v3New);
152     if ((-forwardBackward * (v3New[2] - zEnd)) > 0.0) break; // one is beyond Z spectro. z<0
153     // better use TArray ????
154     for (i3 = 0; i3 < 7; i3++) {v3[i3] = v3New[i3];}
155   }
156   // check fgkMaxStepNumber ????
157   // Interpolation back to exact Z (2nd order)
158   // should be in function ???? using TArray ????
159   Double_t dZ12 = v3New[2] - v3[2]; // 1->2
160   if (TMath::Abs(dZ12) > 0) {
161     Double_t dZ1i = zEnd - v3[2]; // 1-i
162     Double_t dZi2 = v3New[2] - zEnd; // i->2
163     Double_t xPrime = (v3New[0] - v3[0]) / dZ12;
164     Double_t xSecond = ((v3New[3] / v3New[5]) - (v3[3] / v3[5])) / dZ12;
165     Double_t yPrime = (v3New[1] - v3[1]) / dZ12;
166     Double_t ySecond = ((v3New[4] / v3New[5]) - (v3[4] / v3[5])) / dZ12;
167     v3[0] = v3[0] + xPrime * dZ1i - 0.5 * xSecond * dZ1i * dZi2; // X
168     v3[1] = v3[1] + yPrime * dZ1i - 0.5 * ySecond * dZ1i * dZi2; // Y
169     v3[2] = zEnd; // Z
170     Double_t xPrimeI = xPrime - 0.5 * xSecond * (dZi2 - dZ1i);
171     Double_t yPrimeI = yPrime - 0.5 * ySecond * (dZi2 - dZ1i);
172     // (PX, PY, PZ)/PTOT assuming forward motion
173     v3[5] = 1.0 / TMath::Sqrt(1.0 + xPrimeI * xPrimeI + yPrimeI * yPrimeI); // PZ/PTOT
174     v3[3] = xPrimeI * v3[5]; // PX/PTOT
175     v3[4] = yPrimeI * v3[5]; // PY/PTOT
176   } else {
177     cout<<"W-AliMUONTrackExtrap::ExtrapToZHelix: Extrap. to Z not reached, Z = "<<zEnd<<endl;
178   }
179   // Recover track parameters (charge back for forward motion)
180   RecoverTrackParam(v3, chargeExtrap * forwardBackward, trackParam);
181 }
182
183   //__________________________________________________________________________
184 void AliMUONTrackExtrap::ExtrapToZRungekutta(AliMUONTrackParam* trackParam, Double_t zEnd)
185 {
186   /// Track parameter extrapolation to the plane at "Z" using Rungekutta algorithm.
187   /// On return, the track parameters resulting from the extrapolation are updated in trackParam.
188   if (trackParam->GetZ() == zEnd) return; // nothing to be done if same Z
189   Double_t forwardBackward; // +1 if forward, -1 if backward
190   if (zEnd < trackParam->GetZ()) forwardBackward = 1.0; // spectro. z<0 
191   else forwardBackward = -1.0;
192   // sign of charge (sign of fInverseBendingMomentum if forward motion)
193   // must be changed if backward extrapolation
194   Double_t chargeExtrap = forwardBackward * TMath::Sign(Double_t(1.0), trackParam->GetInverseBendingMomentum());
195   Double_t v3[7], v3New[7];
196   Double_t dZ, step;
197   Int_t stepNumber = 0;
198   
199   // Extrapolation loop (until within tolerance)
200   Double_t residue = zEnd - trackParam->GetZ();
201   while (TMath::Abs(residue) > fgkRungeKuttaMaxResidue && stepNumber <= fgkMaxStepNumber) {
202     dZ = zEnd - trackParam->GetZ();
203     // step lenght assuming linear trajectory
204     step = dZ * TMath::Sqrt(1.0 + trackParam->GetBendingSlope()*trackParam->GetBendingSlope() +
205                             trackParam->GetNonBendingSlope()*trackParam->GetNonBendingSlope());
206     ConvertTrackParamForExtrap(trackParam, forwardBackward, v3);
207     do { // reduce step lenght while zEnd oversteped
208       if (stepNumber > fgkMaxStepNumber) {
209         cout<<"W-AliMUONTrackExtrap::ExtrapToZRungekutta: Too many trials: "<<stepNumber<<endl;
210         break;
211       }
212       stepNumber ++;
213       step = TMath::Abs(step);
214       AliMUONTrackExtrap::ExtrapOneStepRungekutta(chargeExtrap,step,v3,v3New);
215       residue = zEnd - v3New[2];
216       step *= dZ/(v3New[2]-trackParam->GetZ());
217     } while (residue*dZ < 0 && TMath::Abs(residue) > fgkRungeKuttaMaxResidue);
218     RecoverTrackParam(v3New, chargeExtrap * forwardBackward, trackParam);
219   }
220   
221   // terminate the extropolation with a straight line up to the exact "zEnd" value
222   trackParam->SetNonBendingCoor(trackParam->GetNonBendingCoor() + residue * trackParam->GetNonBendingSlope());
223   trackParam->SetBendingCoor(trackParam->GetBendingCoor() + residue * trackParam->GetBendingSlope());
224   trackParam->SetZ(zEnd);
225 }
226
227   //__________________________________________________________________________
228 void AliMUONTrackExtrap::ConvertTrackParamForExtrap(AliMUONTrackParam* trackParam, Double_t forwardBackward, Double_t *v3)
229 {
230   /// Set vector of Geant3 parameters pointed to by "v3" from track parameters in trackParam.
231   /// Since AliMUONTrackParam is only geometry, one uses "forwardBackward"
232   /// to know whether the particle is going forward (+1) or backward (-1).
233   v3[0] = trackParam->GetNonBendingCoor(); // X
234   v3[1] = trackParam->GetBendingCoor(); // Y
235   v3[2] = trackParam->GetZ(); // Z
236   Double_t pYZ = TMath::Abs(1.0 / trackParam->GetInverseBendingMomentum());
237   Double_t pZ = pYZ / TMath::Sqrt(1.0 + trackParam->GetBendingSlope() * trackParam->GetBendingSlope());
238   v3[6] = TMath::Sqrt(pYZ * pYZ + pZ * pZ * trackParam->GetNonBendingSlope() * trackParam->GetNonBendingSlope()); // PTOT
239   v3[5] = -forwardBackward * pZ / v3[6]; // PZ/PTOT spectro. z<0
240   v3[3] = trackParam->GetNonBendingSlope() * v3[5]; // PX/PTOT
241   v3[4] = trackParam->GetBendingSlope() * v3[5]; // PY/PTOT
242 }
243
244   //__________________________________________________________________________
245 void AliMUONTrackExtrap::RecoverTrackParam(Double_t *v3, Double_t charge, AliMUONTrackParam* trackParam)
246 {
247   /// Set track parameters in trackParam from Geant3 parameters pointed to by "v3",
248   /// assumed to be calculated for forward motion in Z.
249   /// "InverseBendingMomentum" is signed with "charge".
250   trackParam->SetNonBendingCoor(v3[0]); // X
251   trackParam->SetBendingCoor(v3[1]); // Y
252   trackParam->SetZ(v3[2]); // Z
253   Double_t pYZ = v3[6] * TMath::Sqrt(1.0 - v3[3] * v3[3]);
254   trackParam->SetInverseBendingMomentum(charge/pYZ);
255   trackParam->SetBendingSlope(v3[4]/v3[5]);
256   trackParam->SetNonBendingSlope(v3[3]/v3[5]);
257 }
258
259   //__________________________________________________________________________
260 void AliMUONTrackExtrap::ExtrapToZCov(AliMUONTrackParam* trackParam, Double_t zEnd, Bool_t updatePropagator)
261 {
262   /// Track parameters and their covariances extrapolated to the plane at "zEnd".
263   /// On return, results from the extrapolation are updated in trackParam.
264   
265   if (trackParam->GetZ() == zEnd) return; // nothing to be done if same z
266   
267   // No need to propagate the covariance matrix if it does not exist
268   if (!trackParam->CovariancesExist()) {
269     cout<<"W-AliMUONTrackExtrap::ExtrapToZCov: Covariance matrix does not exist"<<endl;
270     // Extrapolate track parameters to "zEnd"
271     ExtrapToZ(trackParam,zEnd);
272     return;
273   }
274   
275   // Save the actual track parameters
276   AliMUONTrackParam trackParamSave(*trackParam);
277   TMatrixD paramSave(trackParamSave.GetParameters());
278   Double_t zBegin = trackParamSave.GetZ();
279   
280   // Get reference to the parameter covariance matrix
281   const TMatrixD& kParamCov = trackParam->GetCovariances();
282   
283   // Extrapolate track parameters to "zEnd"
284   ExtrapToZ(trackParam,zEnd);
285   
286   // Get reference to the extrapolated parameters
287   const TMatrixD& extrapParam = trackParam->GetParameters();
288   
289   // Calculate the jacobian related to the track parameters extrapolation to "zEnd"
290   TMatrixD jacob(5,5);
291   jacob.Zero();
292   TMatrixD dParam(5,1);
293   for (Int_t i=0; i<5; i++) {
294     // Skip jacobian calculation for parameters with no associated error
295     if (kParamCov(i,i) == 0.) continue;
296     
297     // Small variation of parameter i only
298     for (Int_t j=0; j<5; j++) {
299       if (j==i) {
300         dParam(j,0) = TMath::Sqrt(kParamCov(i,i));
301         if (j == 4) dParam(j,0) *= TMath::Sign(1.,-paramSave(4,0)); // variation always in the same direction
302       } else dParam(j,0) = 0.;
303     }
304     
305     // Set new parameters
306     trackParamSave.SetParameters(paramSave);
307     trackParamSave.AddParameters(dParam);
308     trackParamSave.SetZ(zBegin);
309     
310     // Extrapolate new track parameters to "zEnd"
311     ExtrapToZ(&trackParamSave,zEnd);
312     
313     // Calculate the jacobian
314     TMatrixD jacobji(trackParamSave.GetParameters(),TMatrixD::kMinus,extrapParam);
315     jacobji *= 1. / dParam(i,0);
316     jacob.SetSub(0,i,jacobji);
317   }
318   
319   // Extrapolate track parameter covariances to "zEnd"
320   TMatrixD tmp(kParamCov,TMatrixD::kMultTranspose,jacob);
321   TMatrixD tmp2(jacob,TMatrixD::kMult,tmp);
322   trackParam->SetCovariances(tmp2);
323   
324   // Update the propagator if required
325   if (updatePropagator) trackParam->UpdatePropagator(jacob);
326   
327 }
328
329   //__________________________________________________________________________
330 void AliMUONTrackExtrap::ExtrapToStation(AliMUONTrackParam* trackParamIn, Int_t station, AliMUONTrackParam *trackParamOut)
331 {
332   /// Track parameters extrapolated from "trackParamIn" to both chambers of the station(0..) "station"
333   /// are returned in the array (dimension 2) of track parameters pointed to by "TrackParamOut"
334   /// (index 0 and 1 for first and second chambers).
335   Double_t extZ[2], z1, z2;
336   Int_t i1 = -1, i2 = -1; // = -1 to avoid compilation warnings
337   // range of station to be checked ????
338   z1 = AliMUONConstants::DefaultChamberZ(2 * station);
339   z2 = AliMUONConstants::DefaultChamberZ(2 * station + 1);
340   // First and second Z to extrapolate at
341   if ((z1 > trackParamIn->GetZ()) && (z2 > trackParamIn->GetZ())) {i1 = 0; i2 = 1;}
342   else if ((z1 < trackParamIn->GetZ()) && (z2 < trackParamIn->GetZ())) {i1 = 1; i2 = 0;}
343   else {
344     cout<<"E-AliMUONTrackExtrap::ExtrapToStation: Starting Z ("<<trackParamIn->GetZ()
345         <<") in between z1 ("<<z1<<") and z2 ("<<z2<<") of station(0..)"<<station<<endl;
346     exit(-1);
347   }
348   extZ[i1] = z1;
349   extZ[i2] = z2;
350   // copy of track parameters
351   trackParamOut[i1] = *trackParamIn;
352   // first extrapolation
353   ExtrapToZ(&(trackParamOut[i1]),extZ[0]);
354   trackParamOut[i2] = trackParamOut[i1];
355   // second extrapolation
356   ExtrapToZ(&(trackParamOut[i2]),extZ[1]);
357   return;
358 }
359
360   //__________________________________________________________________________
361 void AliMUONTrackExtrap::ExtrapToVertexUncorrected(AliMUONTrackParam* trackParam, Double_t zVtx)
362 {
363   /// Extrapolation to the vertex (at the z position "zVtx") without Branson and energy loss corrections.
364   /// Returns the track parameters resulting from the extrapolation in the current TrackParam.
365   /// Include multiple Coulomb scattering effects in trackParam covariances.
366   
367   if (trackParam->GetZ() == zVtx) return; // nothing to be done if already at vertex
368   
369   if (trackParam->GetZ() > zVtx) { // spectro. (z<0)
370     cout<<"W-AliMUONTrackExtrap::ExtrapToVertexUncorrected: Starting Z ("<<trackParam->GetZ()
371         <<") upstream the vertex (zVtx = "<<zVtx<<")"<<endl;
372     exit(-1);
373   }
374   
375   // Check the vertex position relatively to the absorber
376   if (zVtx < AliMUONConstants::AbsZBeg() && zVtx > AliMUONConstants::AbsZEnd()) { // spectro. (z<0)
377     cout<<"W-AliMUONTrackExtrap::ExtrapToVertex: Ending Z ("<<zVtx
378         <<") inside the front absorber ("<<AliMUONConstants::AbsZBeg()<<","<<AliMUONConstants::AbsZEnd()<<")"<<endl;
379   } else if (zVtx < AliMUONConstants::AbsZEnd() ) { // spectro. (z<0)
380     cout<<"W-AliMUONTrackExtrap::ExtrapToVertex: Ending Z ("<<zVtx
381         <<") downstream the front absorber (zAbsorberEnd = "<<AliMUONConstants::AbsZEnd()<<")"<<endl;
382     if (trackParam->CovariancesExist()) ExtrapToZCov(trackParam,zVtx);
383     else ExtrapToZ(trackParam,zVtx);
384     return;
385   }
386   
387   // Check the track position relatively to the absorber and extrapolate track parameters to the end of the absorber if needed
388   if (trackParam->GetZ() > AliMUONConstants::AbsZBeg()) { // spectro. (z<0)
389     cout<<"W-AliMUONTrackExtrap::ExtrapToVertex: Starting Z ("<<trackParam->GetZ()
390         <<") upstream the front absorber (zAbsorberBegin = "<<AliMUONConstants::AbsZBeg()<<")"<<endl;
391     if (trackParam->CovariancesExist()) ExtrapToZCov(trackParam,zVtx);
392     else ExtrapToZ(trackParam,zVtx);
393     return;
394   } else if (trackParam->GetZ() > AliMUONConstants::AbsZEnd()) { // spectro. (z<0)
395     cout<<"W-AliMUONTrackExtrap::ExtrapToVertex: Starting Z ("<<trackParam->GetZ()
396         <<") inside the front absorber ("<<AliMUONConstants::AbsZBeg()<<","<<AliMUONConstants::AbsZEnd()<<")"<<endl;
397   } else {
398     if (trackParam->CovariancesExist()) ExtrapToZCov(trackParam,AliMUONConstants::AbsZEnd());
399     else ExtrapToZ(trackParam,AliMUONConstants::AbsZEnd());
400   }
401   
402   // Then add MCS effect in absorber to the parameters covariances
403   AliMUONTrackParam trackParamIn(*trackParam);
404   ExtrapToZ(&trackParamIn, TMath::Min(zVtx, AliMUONConstants::AbsZBeg()));
405   Double_t trackXYZIn[3];
406   trackXYZIn[0] = trackParamIn.GetNonBendingCoor();
407   trackXYZIn[1] = trackParamIn.GetBendingCoor();
408   trackXYZIn[2] = trackParamIn.GetZ();
409   Double_t trackXYZOut[3];
410   trackXYZOut[0] = trackParam->GetNonBendingCoor();
411   trackXYZOut[1] = trackParam->GetBendingCoor();
412   trackXYZOut[2] = trackParam->GetZ();
413   Double_t pathLength = 0.;
414   Double_t f0 = 0.;
415   Double_t f1 = 0.;
416   Double_t f2 = 0.;
417   Double_t meanRho = 0.;
418   GetAbsorberCorrectionParam(trackXYZIn,trackXYZOut,pathLength,f0,f1,f2,meanRho);
419   AddMCSEffectInAbsorber(trackParam,pathLength,f0,f1,f2);
420   
421   // finally go to the vertex
422   ExtrapToZCov(trackParam,zVtx);
423   
424 }
425
426   //__________________________________________________________________________
427 void AliMUONTrackExtrap::AddMCSEffectInAbsorber(AliMUONTrackParam* param, Double_t pathLength, Double_t f0, Double_t f1, Double_t f2)
428 {
429   /// Add to the track parameter covariances the effects of multiple Coulomb scattering
430   /// at the end of the front absorber using the absorber correction parameters
431   
432   // absorber related covariance parameters
433   Double_t bendingSlope = param->GetBendingSlope();
434   Double_t nonBendingSlope = param->GetNonBendingSlope();
435   Double_t inverseBendingMomentum = param->GetInverseBendingMomentum();
436   Double_t alpha2 = 0.0136 * 0.0136 * inverseBendingMomentum * inverseBendingMomentum * (1.0 + bendingSlope * bendingSlope) /
437                                         (1.0 + bendingSlope *bendingSlope + nonBendingSlope * nonBendingSlope); // velocity = 1
438   Double_t varCoor = alpha2 * (pathLength * pathLength * f0 - 2. * pathLength * f1 + f2);
439   Double_t covCorrSlope = alpha2 * (pathLength * f0 - f1);
440   Double_t varSlop = alpha2 * f0;
441   
442   TMatrixD newParamCov(param->GetCovariances());
443   // Non bending plane
444   newParamCov(0,0) += varCoor;       newParamCov(0,1) += covCorrSlope;
445   newParamCov(1,0) += covCorrSlope;  newParamCov(1,1) += varSlop;
446   // Bending plane
447   newParamCov(2,2) += varCoor;       newParamCov(2,3) += covCorrSlope;
448   newParamCov(3,2) += covCorrSlope;  newParamCov(3,3) += varSlop;
449   
450   // Set new covariances
451   param->SetCovariances(newParamCov);
452   
453 }
454
455   //__________________________________________________________________________
456 void AliMUONTrackExtrap::GetAbsorberCorrectionParam(Double_t trackXYZIn[3], Double_t trackXYZOut[3], Double_t &pathLength,
457                                                     Double_t &f0, Double_t &f1, Double_t &f2, Double_t &meanRho)
458 {
459   /// Parameters used to correct for Multiple Coulomb Scattering and energy loss in absorber
460   /// Calculated assuming a linear propagation between track positions trackXYZIn and trackXYZOut
461   // pathLength: path length between trackXYZIn and trackXYZOut (cm)
462   // f0:         0th moment of z calculated with the inverse radiation-length distribution
463   // f1:         1st moment of z calculated with the inverse radiation-length distribution
464   // f2:         2nd moment of z calculated with the inverse radiation-length distribution
465   // meanRho:    average density of crossed material (g/cm3)
466   
467   // Reset absorber's parameters
468   pathLength = 0.;
469   f0 = 0.;
470   f1 = 0.;
471   f2 = 0.;
472   meanRho = 0.;
473   
474   // Check whether the geometry is available
475   if (!gGeoManager) {
476     cout<<"E-AliMUONTrackExtrap::GetAbsorberCorrectionParam: no TGeo"<<endl;
477     return;
478   }
479   
480   // Initialize starting point and direction
481   pathLength = TMath::Sqrt((trackXYZOut[0] - trackXYZIn[0])*(trackXYZOut[0] - trackXYZIn[0])+
482                            (trackXYZOut[1] - trackXYZIn[1])*(trackXYZOut[1] - trackXYZIn[1])+
483                            (trackXYZOut[2] - trackXYZIn[2])*(trackXYZOut[2] - trackXYZIn[2]));
484   if (pathLength < TGeoShape::Tolerance()) return;
485   Double_t b[3];
486   b[0] = (trackXYZOut[0] - trackXYZIn[0]) / pathLength;
487   b[1] = (trackXYZOut[1] - trackXYZIn[1]) / pathLength;
488   b[2] = (trackXYZOut[2] - trackXYZIn[2]) / pathLength;
489   TGeoNode *currentnode = gGeoManager->InitTrack(trackXYZIn, b);
490   if (!currentnode) {
491     cout<<"E-AliMUONTrackExtrap::GetAbsorberCorrectionParam: start point out of geometry"<<endl;
492     return;
493   }
494   
495   // loop over absorber slices and calculate absorber's parameters
496   Double_t rho = 0.; // material density (g/cm3)
497   Double_t x0 = 0.;  // radiation-length (cm-1)
498   Double_t localPathLength = 0;
499   Double_t remainingPathLength = pathLength;
500   Double_t zB = trackXYZIn[2];
501   Double_t zE, dzB, dzE;
502   do {
503     // Get material properties
504     TGeoMaterial *material = currentnode->GetVolume()->GetMedium()->GetMaterial();
505     rho = material->GetDensity();
506     x0 = material->GetRadLen();
507     if (!material->IsMixture()) x0 /= rho; // different normalization in the modeler for mixture
508     
509     // Get path length within this material
510     gGeoManager->FindNextBoundary(remainingPathLength);
511     localPathLength = gGeoManager->GetStep() + 1.e-6;
512     // Check if boundary within remaining path length. If so, make sure to cross the boundary to prepare the next step
513     if (localPathLength >= remainingPathLength) localPathLength = remainingPathLength;
514     else {
515       currentnode = gGeoManager->Step();
516       if (!currentnode) {
517         cout<<"E-AliMUONTrackExtrap::GetAbsorberCorrectionParam: navigation failed"<<endl;
518         f0 = f1 = f2 = meanRho = 0.;
519         return;
520       }
521       if (!gGeoManager->IsEntering()) {
522         // make another small step to try to enter in new absorber slice
523         gGeoManager->SetStep(0.001);
524         currentnode = gGeoManager->Step();
525         if (!gGeoManager->IsEntering() || !currentnode) {
526           cout<<"E-AliMUONTrackExtrap::GetAbsorberCorrectionParam: navigation failed"<<endl;
527           f0 = f1 = f2 = meanRho = 0.;
528           return;
529         }
530         localPathLength += 0.001;
531       }
532     }
533     
534     // calculate absorber's parameters
535     zE = b[2] * localPathLength + zB;
536     dzB = zB - trackXYZIn[2];
537     dzE = zE - trackXYZIn[2];
538     f0 += localPathLength / x0;
539     f1 += (dzE*dzE - dzB*dzB) / b[2] / b[2] / x0 / 2.;
540     f2 += (dzE*dzE*dzE - dzB*dzB*dzB) / b[2] / b[2] / b[2] / x0 / 3.;
541     meanRho += localPathLength * rho;
542     
543     // prepare next step
544     zB = zE;
545     remainingPathLength -= localPathLength;
546   } while (remainingPathLength > TGeoShape::Tolerance());
547   
548   meanRho /= pathLength;
549 }
550
551   //__________________________________________________________________________
552 Double_t AliMUONTrackExtrap::GetMCSAngle2(const AliMUONTrackParam& param, Double_t dZ, Double_t x0)
553 {
554   /// Return the angular dispersion square due to multiple Coulomb scattering
555   /// through a material of thickness "dZ" and of radiation length "x0"
556   /// assuming linear propagation and using the small angle approximation.
557   
558   Double_t bendingSlope = param.GetBendingSlope();
559   Double_t nonBendingSlope = param.GetNonBendingSlope();
560   Double_t inverseTotalMomentum2 = param.GetInverseBendingMomentum() * param.GetInverseBendingMomentum() *
561                                    (1.0 + bendingSlope * bendingSlope) /
562                                    (1.0 + bendingSlope *bendingSlope + nonBendingSlope * nonBendingSlope); 
563   // Path length in the material
564   Double_t pathLength = TMath::Abs(dZ) * TMath::Sqrt(1.0 + bendingSlope*bendingSlope + nonBendingSlope*nonBendingSlope);
565   // relativistic velocity
566   Double_t velo = 1.;
567   // Angular dispersion square of the track (variance) in a plane perpendicular to the trajectory
568   Double_t theta02 = 0.0136 / velo * (1 + 0.038 * TMath::Log(pathLength/x0));
569   
570   return theta02 * theta02 * inverseTotalMomentum2 * pathLength / x0;
571 }
572
573   //__________________________________________________________________________
574 void AliMUONTrackExtrap::AddMCSEffect(AliMUONTrackParam *param, Double_t dZ, Double_t x0)
575 {
576   /// Add to the track parameter covariances the effects of multiple Coulomb scattering
577   /// through a material of thickness "dZ" and of radiation length "x0"
578   /// assuming linear propagation and using the small angle approximation.
579   
580   Double_t bendingSlope = param->GetBendingSlope();
581   Double_t nonBendingSlope = param->GetNonBendingSlope();
582   Double_t inverseTotalMomentum2 = param->GetInverseBendingMomentum() * param->GetInverseBendingMomentum() *
583                                    (1.0 + bendingSlope * bendingSlope) /
584                                    (1.0 + bendingSlope *bendingSlope + nonBendingSlope * nonBendingSlope); 
585   // Path length in the material
586   Double_t pathLength = TMath::Abs(dZ) * TMath::Sqrt(1.0 + bendingSlope*bendingSlope + nonBendingSlope*nonBendingSlope);
587   Double_t pathLength2 = pathLength * pathLength;
588   // relativistic velocity
589   Double_t velo = 1.;
590   // Angular dispersion square of the track (variance) in a plane perpendicular to the trajectory
591   Double_t theta02 = 0.0136 / velo * (1 + 0.038 * TMath::Log(pathLength/x0));
592   theta02 *= theta02 * inverseTotalMomentum2 * pathLength / x0;
593   
594   Double_t varCoor      = pathLength2 * theta02 / 3.;
595   Double_t varSlop      = theta02;
596   Double_t covCorrSlope = pathLength * theta02 / 2.;
597   
598   // Add effects of multiple Coulomb scattering in track parameter covariances
599   TMatrixD newParamCov(param->GetCovariances());
600   // Non bending plane
601   newParamCov(0,0) += varCoor;       newParamCov(0,1) += covCorrSlope;
602   newParamCov(1,0) += covCorrSlope;  newParamCov(1,1) += varSlop;
603   // Bending plane
604   newParamCov(2,2) += varCoor;       newParamCov(2,3) += covCorrSlope;
605   newParamCov(3,2) += covCorrSlope;  newParamCov(3,3) += varSlop;
606   
607   // Set new covariances
608   param->SetCovariances(newParamCov);
609 }
610
611   //__________________________________________________________________________
612 void AliMUONTrackExtrap::ExtrapToVertex(AliMUONTrackParam* trackParam, Double_t xVtx, Double_t yVtx, Double_t zVtx,
613                                         Bool_t CorrectForMCS, Bool_t CorrectForEnergyLoss)
614 {
615   /// Extrapolation to the vertex.
616   /// Returns the track parameters resulting from the extrapolation of the current TrackParam.
617   /// Changes parameters according to Branson correction through the absorber and energy loss
618   
619   if (trackParam->GetZ() == zVtx) return; // nothing to be done if already at vertex
620   
621   if (trackParam->GetZ() > zVtx) { // spectro. (z<0)
622     cout<<"F-AliMUONTrackExtrap::ExtrapToVertex: Starting Z ("<<trackParam->GetZ()
623         <<") upstream the vertex (zVtx = "<<zVtx<<")"<<endl;
624     exit(-1);
625   }
626   
627   // Check if correction required
628   if (!CorrectForMCS && !CorrectForEnergyLoss) {
629     ExtrapToZ(trackParam,zVtx);
630     return;
631   }
632   
633   // Check the vertex position relatively to the absorber
634   if (zVtx < AliMUONConstants::AbsZBeg() && zVtx > AliMUONConstants::AbsZEnd()) { // spectro. (z<0)
635     cout<<"W-AliMUONTrackExtrap::ExtrapToVertex: Ending Z ("<<zVtx
636         <<") inside the front absorber ("<<AliMUONConstants::AbsZBeg()<<","<<AliMUONConstants::AbsZEnd()<<")"<<endl;
637   } else if (zVtx < AliMUONConstants::AbsZEnd() ) { // spectro. (z<0)
638     cout<<"W-AliMUONTrackExtrap::ExtrapToVertex: Ending Z ("<<zVtx
639         <<") downstream the front absorber (zAbsorberEnd = "<<AliMUONConstants::AbsZEnd()<<")"<<endl;
640     ExtrapToZ(trackParam,zVtx);
641     return;
642   }
643   
644   // Check the track position relatively to the absorber and extrapolate track parameters to the end of the absorber if needed
645   if (trackParam->GetZ() > AliMUONConstants::AbsZBeg()) { // spectro. (z<0)
646     cout<<"W-AliMUONTrackExtrap::ExtrapToVertex: Starting Z ("<<trackParam->GetZ()
647         <<") upstream the front absorber (zAbsorberBegin = "<<AliMUONConstants::AbsZBeg()<<")"<<endl;
648     ExtrapToZ(trackParam,zVtx);
649     return;
650   } else if (trackParam->GetZ() > AliMUONConstants::AbsZEnd()) { // spectro. (z<0)
651     cout<<"W-AliMUONTrackExtrap::ExtrapToVertex: Starting Z ("<<trackParam->GetZ()
652         <<") inside the front absorber ("<<AliMUONConstants::AbsZBeg()<<","<<AliMUONConstants::AbsZEnd()<<")"<<endl;
653   } else {
654     ExtrapToZ(trackParam,AliMUONConstants::AbsZEnd());
655   }
656   
657   // Get absorber correction parameters assuming linear propagation from vertex to the track position
658   Double_t trackXYZOut[3];
659   trackXYZOut[0] = trackParam->GetNonBendingCoor();
660   trackXYZOut[1] = trackParam->GetBendingCoor();
661   trackXYZOut[2] = trackParam->GetZ();
662   Double_t trackXYZIn[3];
663   trackXYZIn[2] = TMath::Min(zVtx, AliMUONConstants::AbsZBeg()); // spectro. (z<0)
664   trackXYZIn[0] = trackXYZOut[0] + (xVtx - trackXYZOut[0]) / (zVtx - trackXYZOut[2]) * (trackXYZIn[2] - trackXYZOut[2]);
665   trackXYZIn[1] = trackXYZOut[1] + (yVtx - trackXYZOut[1]) / (zVtx - trackXYZOut[2]) * (trackXYZIn[2] - trackXYZOut[2]);
666   Double_t pathLength = 0.;
667   Double_t f0 = 0.;
668   Double_t f1 = 0.;
669   Double_t f2 = 0.;
670   Double_t meanRho = 0.;
671   GetAbsorberCorrectionParam(trackXYZIn,trackXYZOut,pathLength,f0,f1,f2,meanRho);
672   
673   // Calculate energy loss
674   Double_t pTot = trackParam->P();
675   Double_t charge = TMath::Sign(Double_t(1.0), trackParam->GetInverseBendingMomentum());      
676   Double_t deltaP = TotalMomentumEnergyLoss(pTot,pathLength,meanRho);
677   
678   // Correct for half of energy loss
679   Double_t nonBendingSlope, bendingSlope;
680   if (CorrectForEnergyLoss) {
681     pTot += 0.5 * deltaP;
682     nonBendingSlope = trackParam->GetNonBendingSlope();
683     bendingSlope = trackParam->GetBendingSlope();
684     trackParam->SetInverseBendingMomentum(charge / pTot *
685           TMath::Sqrt(1.0 + nonBendingSlope*nonBendingSlope + bendingSlope*bendingSlope) /
686           TMath::Sqrt(1.0 + bendingSlope*bendingSlope));
687   }
688   
689   if (CorrectForMCS) {
690     // Position of the Branson plane (spectro. (z<0))
691     Double_t zB = (f1>0.) ? trackXYZIn[2] - f2/f1 : 0.;
692     
693     // Get track position in the Branson plane corrected for magnetic field effect
694     ExtrapToZ(trackParam,zVtx);
695     Double_t xB = trackParam->GetNonBendingCoor() + (zB - zVtx) * trackParam->GetNonBendingSlope();
696     Double_t yB = trackParam->GetBendingCoor()    + (zB - zVtx) * trackParam->GetBendingSlope();
697     
698     // Get track slopes corrected for multiple scattering (spectro. (z<0))
699     nonBendingSlope = (zB<0.) ? (xB - xVtx) / (zB - zVtx) : trackParam->GetNonBendingSlope();
700     bendingSlope    = (zB<0.) ? (yB - yVtx) / (zB - zVtx) : trackParam->GetBendingSlope();
701     
702     // Set track parameters at vertex
703     trackParam->SetNonBendingCoor(xVtx);
704     trackParam->SetBendingCoor(yVtx);
705     trackParam->SetZ(zVtx);
706     trackParam->SetNonBendingSlope(nonBendingSlope);
707     trackParam->SetBendingSlope(bendingSlope);
708   } else {
709     ExtrapToZ(trackParam,zVtx);
710     nonBendingSlope = trackParam->GetNonBendingSlope();
711     bendingSlope = trackParam->GetBendingSlope();
712   }
713   
714   // Correct for second half of energy loss
715   if (CorrectForEnergyLoss) pTot += 0.5 * deltaP;
716   
717   // Set track parameters at vertex
718   trackParam->SetInverseBendingMomentum(charge / pTot *
719         TMath::Sqrt(1.0 + nonBendingSlope*nonBendingSlope + bendingSlope*bendingSlope) /
720         TMath::Sqrt(1.0 + bendingSlope*bendingSlope));
721   
722 }
723
724   //__________________________________________________________________________
725 Double_t AliMUONTrackExtrap::TotalMomentumEnergyLoss(AliMUONTrackParam* trackParam, Double_t xVtx, Double_t yVtx, Double_t zVtx)
726 {
727   /// Calculate the total momentum energy loss in-between the track position and the vertex assuming a linear propagation
728   
729   if (trackParam->GetZ() == zVtx) return 0.; // nothing to be done if already at vertex
730   
731   // Check whether the geometry is available
732   if (!gGeoManager) {
733     cout<<"E-AliMUONTrackExtrap::TotalMomentumEnergyLoss: no TGeo"<<endl;
734     return 0.;
735   }
736   
737   // Get encountered material correction parameters assuming linear propagation from vertex to the track position
738   Double_t trackXYZOut[3];
739   trackXYZOut[0] = trackParam->GetNonBendingCoor();
740   trackXYZOut[1] = trackParam->GetBendingCoor();
741   trackXYZOut[2] = trackParam->GetZ();
742   Double_t trackXYZIn[3];
743   trackXYZIn[0] = xVtx;
744   trackXYZIn[1] = yVtx;
745   trackXYZIn[2] = zVtx;
746   Double_t pathLength = 0.;
747   Double_t f0 = 0.;
748   Double_t f1 = 0.;
749   Double_t f2 = 0.;
750   Double_t meanRho = 0.;
751   GetAbsorberCorrectionParam(trackXYZIn,trackXYZOut,pathLength,f0,f1,f2,meanRho);
752   
753   // Calculate energy loss
754   Double_t pTot = trackParam->P();
755   return TotalMomentumEnergyLoss(pTot,pathLength,meanRho);
756 }
757
758   //__________________________________________________________________________
759 Double_t AliMUONTrackExtrap::TotalMomentumEnergyLoss(Double_t pTotal, Double_t pathLength, Double_t rho)
760 {
761   /// Returns the total momentum energy loss in the front absorber
762   Double_t muMass = 0.105658369;
763   Double_t p2=pTotal*pTotal;
764   Double_t beta2=p2/(p2 + muMass*muMass);
765   Double_t dE=ApproximateBetheBloch(beta2)*pathLength*rho;
766   
767   return dE;
768 }
769
770   //__________________________________________________________________________
771 Double_t AliMUONTrackExtrap::ApproximateBetheBloch(Double_t beta2) 
772 {
773 /// This is an approximation of the Bethe-Bloch formula with 
774 /// the density effect taken into account at beta*gamma > 3.5
775 /// (the approximation is reasonable only for solid materials) 
776
777   if (beta2/(1-beta2)>3.5*3.5)
778      return 0.153e-3/beta2*(log(3.5*5940)+0.5*log(beta2/(1-beta2)) - beta2);
779
780   return 0.153e-3/beta2*(log(5940*beta2/(1-beta2)) - beta2);
781 }
782
783  //__________________________________________________________________________
784 void AliMUONTrackExtrap::ExtrapOneStepHelix(Double_t charge, Double_t step, Double_t *vect, Double_t *vout)
785 {
786 /// <pre>
787 ///    ******************************************************************
788 ///    *                                                                *
789 ///    *  Performs the tracking of one step in a magnetic field         *
790 ///    *  The trajectory is assumed to be a helix in a constant field   *
791 ///    *  taken at the mid point of the step.                           *
792 ///    *  Parameters:                                                   *
793 ///    *   input                                                        *
794 ///    *     STEP =arc length of the step asked                         *
795 ///    *     VECT =input vector (position,direction cos and momentum)   *
796 ///    *     CHARGE=  electric charge of the particle                   *
797 ///    *   output                                                       *
798 ///    *     VOUT = same as VECT after completion of the step           *
799 ///    *                                                                *
800 ///    *    ==>Called by : <USER>, GUSWIM                               *
801 ///    *       Author    m.hansroul  *********                          *
802 ///    *       modified  s.egli, s.v.levonian                           *
803 ///    *       modified  v.perevoztchikov
804 ///    *                                                                *
805 ///    ******************************************************************
806 /// </pre>
807
808 // modif: everything in double precision
809
810     Double_t xyz[3], h[4], hxp[3];
811     Double_t h2xy, hp, rho, tet;
812     Double_t sint, sintt, tsint, cos1t;
813     Double_t f1, f2, f3, f4, f5, f6;
814
815     const Int_t kix  = 0;
816     const Int_t kiy  = 1;
817     const Int_t kiz  = 2;
818     const Int_t kipx = 3;
819     const Int_t kipy = 4;
820     const Int_t kipz = 5;
821     const Int_t kipp = 6;
822
823     const Double_t kec = 2.9979251e-4;
824     //
825     //    ------------------------------------------------------------------
826     //
827     //       units are kgauss,centimeters,gev/c
828     //
829     vout[kipp] = vect[kipp];
830     if (TMath::Abs(charge) < 0.00001) {
831       for (Int_t i = 0; i < 3; i++) {
832         vout[i] = vect[i] + step * vect[i+3];
833         vout[i+3] = vect[i+3];
834       }
835       return;
836     }
837     xyz[0]    = vect[kix] + 0.5 * step * vect[kipx];
838     xyz[1]    = vect[kiy] + 0.5 * step * vect[kipy];
839     xyz[2]    = vect[kiz] + 0.5 * step * vect[kipz];
840
841     //cmodif: call gufld (xyz, h) changed into:
842     GetField (xyz, h);
843  
844     h2xy = h[0]*h[0] + h[1]*h[1];
845     h[3] = h[2]*h[2]+ h2xy;
846     if (h[3] < 1.e-12) {
847       for (Int_t i = 0; i < 3; i++) {
848         vout[i] = vect[i] + step * vect[i+3];
849         vout[i+3] = vect[i+3];
850       }
851       return;
852     }
853     if (h2xy < 1.e-12*h[3]) {
854       ExtrapOneStepHelix3(charge*h[2], step, vect, vout);
855       return;
856     }
857     h[3] = TMath::Sqrt(h[3]);
858     h[0] /= h[3];
859     h[1] /= h[3];
860     h[2] /= h[3];
861     h[3] *= kec;
862
863     hxp[0] = h[1]*vect[kipz] - h[2]*vect[kipy];
864     hxp[1] = h[2]*vect[kipx] - h[0]*vect[kipz];
865     hxp[2] = h[0]*vect[kipy] - h[1]*vect[kipx];
866  
867     hp = h[0]*vect[kipx] + h[1]*vect[kipy] + h[2]*vect[kipz];
868
869     rho = -charge*h[3]/vect[kipp];
870     tet = rho * step;
871
872     if (TMath::Abs(tet) > 0.15) {
873       sint = TMath::Sin(tet);
874       sintt = (sint/tet);
875       tsint = (tet-sint)/tet;
876       cos1t = 2.*(TMath::Sin(0.5*tet))*(TMath::Sin(0.5*tet))/tet;
877     } else {
878       tsint = tet*tet/36.;
879       sintt = (1. - tsint);
880       sint = tet*sintt;
881       cos1t = 0.5*tet;
882     }
883
884     f1 = step * sintt;
885     f2 = step * cos1t;
886     f3 = step * tsint * hp;
887     f4 = -tet*cos1t;
888     f5 = sint;
889     f6 = tet * cos1t * hp;
890  
891     vout[kix] = vect[kix] + f1*vect[kipx] + f2*hxp[0] + f3*h[0];
892     vout[kiy] = vect[kiy] + f1*vect[kipy] + f2*hxp[1] + f3*h[1];
893     vout[kiz] = vect[kiz] + f1*vect[kipz] + f2*hxp[2] + f3*h[2];
894  
895     vout[kipx] = vect[kipx] + f4*vect[kipx] + f5*hxp[0] + f6*h[0];
896     vout[kipy] = vect[kipy] + f4*vect[kipy] + f5*hxp[1] + f6*h[1];
897     vout[kipz] = vect[kipz] + f4*vect[kipz] + f5*hxp[2] + f6*h[2];
898  
899     return;
900 }
901
902  //__________________________________________________________________________
903 void AliMUONTrackExtrap::ExtrapOneStepHelix3(Double_t field, Double_t step, Double_t *vect, Double_t *vout)
904 {
905 /// <pre>
906 ///     ******************************************************************
907 ///     *                                                                *
908 ///     *       Tracking routine in a constant field oriented            *
909 ///     *       along axis 3                                             *
910 ///     *       Tracking is performed with a conventional                *
911 ///     *       helix step method                                        *
912 ///     *                                                                *
913 ///     *    ==>Called by : <USER>, GUSWIM                               *
914 ///     *       Authors    R.Brun, M.Hansroul  *********                 *
915 ///     *       Rewritten  V.Perevoztchikov
916 ///     *                                                                *
917 ///     ******************************************************************
918 /// </pre>
919
920     Double_t hxp[3];
921     Double_t h4, hp, rho, tet;
922     Double_t sint, sintt, tsint, cos1t;
923     Double_t f1, f2, f3, f4, f5, f6;
924
925     const Int_t kix  = 0;
926     const Int_t kiy  = 1;
927     const Int_t kiz  = 2;
928     const Int_t kipx = 3;
929     const Int_t kipy = 4;
930     const Int_t kipz = 5;
931     const Int_t kipp = 6;
932
933     const Double_t kec = 2.9979251e-4;
934
935 // 
936 //     ------------------------------------------------------------------
937 // 
938 //       units are kgauss,centimeters,gev/c
939 // 
940     vout[kipp] = vect[kipp];
941     h4 = field * kec;
942
943     hxp[0] = - vect[kipy];
944     hxp[1] = + vect[kipx];
945  
946     hp = vect[kipz];
947
948     rho = -h4/vect[kipp];
949     tet = rho * step;
950     if (TMath::Abs(tet) > 0.15) {
951       sint = TMath::Sin(tet);
952       sintt = (sint/tet);
953       tsint = (tet-sint)/tet;
954       cos1t = 2.* TMath::Sin(0.5*tet) * TMath::Sin(0.5*tet)/tet;
955     } else {
956       tsint = tet*tet/36.;
957       sintt = (1. - tsint);
958       sint = tet*sintt;
959       cos1t = 0.5*tet;
960     }
961
962     f1 = step * sintt;
963     f2 = step * cos1t;
964     f3 = step * tsint * hp;
965     f4 = -tet*cos1t;
966     f5 = sint;
967     f6 = tet * cos1t * hp;
968  
969     vout[kix] = vect[kix] + f1*vect[kipx] + f2*hxp[0];
970     vout[kiy] = vect[kiy] + f1*vect[kipy] + f2*hxp[1];
971     vout[kiz] = vect[kiz] + f1*vect[kipz] + f3;
972  
973     vout[kipx] = vect[kipx] + f4*vect[kipx] + f5*hxp[0];
974     vout[kipy] = vect[kipy] + f4*vect[kipy] + f5*hxp[1];
975     vout[kipz] = vect[kipz] + f4*vect[kipz] + f6;
976
977     return;
978 }
979
980  //__________________________________________________________________________
981 void AliMUONTrackExtrap::ExtrapOneStepRungekutta(Double_t charge, Double_t step, Double_t* vect, Double_t* vout)
982 {
983 /// <pre>
984 ///     ******************************************************************
985 ///     *                                                                *
986 ///     *  Runge-Kutta method for tracking a particle through a magnetic *
987 ///     *  field. Uses Nystroem algorithm (See Handbook Nat. Bur. of     *
988 ///     *  Standards, procedure 25.5.20)                                 *
989 ///     *                                                                *
990 ///     *  Input parameters                                              *
991 ///     *       CHARGE    Particle charge                                *
992 ///     *       STEP      Step size                                      *
993 ///     *       VECT      Initial co-ords,direction cosines,momentum     *
994 ///     *  Output parameters                                             *
995 ///     *       VOUT      Output co-ords,direction cosines,momentum      *
996 ///     *  User routine called                                           *
997 ///     *       CALL GUFLD(X,F)                                          *
998 ///     *                                                                *
999 ///     *    ==>Called by : <USER>, GUSWIM                               *
1000 ///     *       Authors    R.Brun, M.Hansroul  *********                 *
1001 ///     *                  V.Perevoztchikov (CUT STEP implementation)    *
1002 ///     *                                                                *
1003 ///     *                                                                *
1004 ///     ******************************************************************
1005 /// </pre>
1006
1007     Double_t h2, h4, f[4];
1008     Double_t xyzt[3], a, b, c, ph,ph2;
1009     Double_t secxs[4],secys[4],seczs[4],hxp[3];
1010     Double_t g1, g2, g3, g4, g5, g6, ang2, dxt, dyt, dzt;
1011     Double_t est, at, bt, ct, cba;
1012     Double_t f1, f2, f3, f4, rho, tet, hnorm, hp, rho1, sint, cost;
1013     
1014     Double_t x;
1015     Double_t y;
1016     Double_t z;
1017     
1018     Double_t xt;
1019     Double_t yt;
1020     Double_t zt;
1021
1022     Double_t maxit = 1992;
1023     Double_t maxcut = 11;
1024
1025     const Double_t kdlt   = 1e-4;
1026     const Double_t kdlt32 = kdlt/32.;
1027     const Double_t kthird = 1./3.;
1028     const Double_t khalf  = 0.5;
1029     const Double_t kec = 2.9979251e-4;
1030
1031     const Double_t kpisqua = 9.86960440109;
1032     const Int_t kix  = 0;
1033     const Int_t kiy  = 1;
1034     const Int_t kiz  = 2;
1035     const Int_t kipx = 3;
1036     const Int_t kipy = 4;
1037     const Int_t kipz = 5;
1038   
1039     // *.
1040     // *.    ------------------------------------------------------------------
1041     // *.
1042     // *             this constant is for units cm,gev/c and kgauss
1043     // *
1044     Int_t iter = 0;
1045     Int_t ncut = 0;
1046     for(Int_t j = 0; j < 7; j++)
1047       vout[j] = vect[j];
1048
1049     Double_t  pinv   = kec * charge / vect[6];
1050     Double_t tl = 0.;
1051     Double_t h = step;
1052     Double_t rest;
1053
1054  
1055     do {
1056       rest  = step - tl;
1057       if (TMath::Abs(h) > TMath::Abs(rest)) h = rest;
1058       //cmodif: call gufld(vout,f) changed into:
1059
1060       GetField(vout,f);
1061
1062       // *
1063       // *             start of integration
1064       // *
1065       x      = vout[0];
1066       y      = vout[1];
1067       z      = vout[2];
1068       a      = vout[3];
1069       b      = vout[4];
1070       c      = vout[5];
1071
1072       h2     = khalf * h;
1073       h4     = khalf * h2;
1074       ph     = pinv * h;
1075       ph2    = khalf * ph;
1076       secxs[0] = (b * f[2] - c * f[1]) * ph2;
1077       secys[0] = (c * f[0] - a * f[2]) * ph2;
1078       seczs[0] = (a * f[1] - b * f[0]) * ph2;
1079       ang2 = (secxs[0]*secxs[0] + secys[0]*secys[0] + seczs[0]*seczs[0]);
1080       if (ang2 > kpisqua) break;
1081
1082       dxt    = h2 * a + h4 * secxs[0];
1083       dyt    = h2 * b + h4 * secys[0];
1084       dzt    = h2 * c + h4 * seczs[0];
1085       xt     = x + dxt;
1086       yt     = y + dyt;
1087       zt     = z + dzt;
1088       // *
1089       // *              second intermediate point
1090       // *
1091
1092       est = TMath::Abs(dxt) + TMath::Abs(dyt) + TMath::Abs(dzt);
1093       if (est > h) {
1094         if (ncut++ > maxcut) break;
1095         h *= khalf;
1096         continue;
1097       }
1098  
1099       xyzt[0] = xt;
1100       xyzt[1] = yt;
1101       xyzt[2] = zt;
1102
1103       //cmodif: call gufld(xyzt,f) changed into:
1104       GetField(xyzt,f);
1105
1106       at     = a + secxs[0];
1107       bt     = b + secys[0];
1108       ct     = c + seczs[0];
1109
1110       secxs[1] = (bt * f[2] - ct * f[1]) * ph2;
1111       secys[1] = (ct * f[0] - at * f[2]) * ph2;
1112       seczs[1] = (at * f[1] - bt * f[0]) * ph2;
1113       at     = a + secxs[1];
1114       bt     = b + secys[1];
1115       ct     = c + seczs[1];
1116       secxs[2] = (bt * f[2] - ct * f[1]) * ph2;
1117       secys[2] = (ct * f[0] - at * f[2]) * ph2;
1118       seczs[2] = (at * f[1] - bt * f[0]) * ph2;
1119       dxt    = h * (a + secxs[2]);
1120       dyt    = h * (b + secys[2]);
1121       dzt    = h * (c + seczs[2]);
1122       xt     = x + dxt;
1123       yt     = y + dyt;
1124       zt     = z + dzt;
1125       at     = a + 2.*secxs[2];
1126       bt     = b + 2.*secys[2];
1127       ct     = c + 2.*seczs[2];
1128
1129       est = TMath::Abs(dxt)+TMath::Abs(dyt)+TMath::Abs(dzt);
1130       if (est > 2.*TMath::Abs(h)) {
1131         if (ncut++ > maxcut) break;
1132         h *= khalf;
1133         continue;
1134       }
1135  
1136       xyzt[0] = xt;
1137       xyzt[1] = yt;
1138       xyzt[2] = zt;
1139
1140       //cmodif: call gufld(xyzt,f) changed into:
1141       GetField(xyzt,f);
1142
1143       z      = z + (c + (seczs[0] + seczs[1] + seczs[2]) * kthird) * h;
1144       y      = y + (b + (secys[0] + secys[1] + secys[2]) * kthird) * h;
1145       x      = x + (a + (secxs[0] + secxs[1] + secxs[2]) * kthird) * h;
1146
1147       secxs[3] = (bt*f[2] - ct*f[1])* ph2;
1148       secys[3] = (ct*f[0] - at*f[2])* ph2;
1149       seczs[3] = (at*f[1] - bt*f[0])* ph2;
1150       a      = a+(secxs[0]+secxs[3]+2. * (secxs[1]+secxs[2])) * kthird;
1151       b      = b+(secys[0]+secys[3]+2. * (secys[1]+secys[2])) * kthird;
1152       c      = c+(seczs[0]+seczs[3]+2. * (seczs[1]+seczs[2])) * kthird;
1153
1154       est    = TMath::Abs(secxs[0]+secxs[3] - (secxs[1]+secxs[2]))
1155         + TMath::Abs(secys[0]+secys[3] - (secys[1]+secys[2]))
1156         + TMath::Abs(seczs[0]+seczs[3] - (seczs[1]+seczs[2]));
1157
1158       if (est > kdlt && TMath::Abs(h) > 1.e-4) {
1159         if (ncut++ > maxcut) break;
1160         h *= khalf;
1161         continue;
1162       }
1163
1164       ncut = 0;
1165       // *               if too many iterations, go to helix
1166       if (iter++ > maxit) break;
1167
1168       tl += h;
1169       if (est < kdlt32) 
1170         h *= 2.;
1171       cba    = 1./ TMath::Sqrt(a*a + b*b + c*c);
1172       vout[0] = x;
1173       vout[1] = y;
1174       vout[2] = z;
1175       vout[3] = cba*a;
1176       vout[4] = cba*b;
1177       vout[5] = cba*c;
1178       rest = step - tl;
1179       if (step < 0.) rest = -rest;
1180       if (rest < 1.e-5*TMath::Abs(step)) return;
1181
1182     } while(1);
1183
1184     // angle too big, use helix
1185
1186     f1  = f[0];
1187     f2  = f[1];
1188     f3  = f[2];
1189     f4  = TMath::Sqrt(f1*f1+f2*f2+f3*f3);
1190     rho = -f4*pinv;
1191     tet = rho * step;
1192  
1193     hnorm = 1./f4;
1194     f1 = f1*hnorm;
1195     f2 = f2*hnorm;
1196     f3 = f3*hnorm;
1197
1198     hxp[0] = f2*vect[kipz] - f3*vect[kipy];
1199     hxp[1] = f3*vect[kipx] - f1*vect[kipz];
1200     hxp[2] = f1*vect[kipy] - f2*vect[kipx];
1201  
1202     hp = f1*vect[kipx] + f2*vect[kipy] + f3*vect[kipz];
1203
1204     rho1 = 1./rho;
1205     sint = TMath::Sin(tet);
1206     cost = 2.*TMath::Sin(khalf*tet)*TMath::Sin(khalf*tet);
1207
1208     g1 = sint*rho1;
1209     g2 = cost*rho1;
1210     g3 = (tet-sint) * hp*rho1;
1211     g4 = -cost;
1212     g5 = sint;
1213     g6 = cost * hp;
1214  
1215     vout[kix] = vect[kix] + g1*vect[kipx] + g2*hxp[0] + g3*f1;
1216     vout[kiy] = vect[kiy] + g1*vect[kipy] + g2*hxp[1] + g3*f2;
1217     vout[kiz] = vect[kiz] + g1*vect[kipz] + g2*hxp[2] + g3*f3;
1218  
1219     vout[kipx] = vect[kipx] + g4*vect[kipx] + g5*hxp[0] + g6*f1;
1220     vout[kipy] = vect[kipy] + g4*vect[kipy] + g5*hxp[1] + g6*f2;
1221     vout[kipz] = vect[kipz] + g4*vect[kipz] + g5*hxp[2] + g6*f3;
1222
1223     return;
1224 }
1225
1226 //___________________________________________________________
1227  void  AliMUONTrackExtrap::GetField(Double_t *Position, Double_t *Field)
1228 {
1229   /// interface for arguments in double precision (Why ? ChF)
1230   Float_t x[3], b[3];
1231
1232   x[0] = Position[0]; x[1] = Position[1]; x[2] = Position[2];
1233
1234   if (fgkField) fgkField->Field(x,b);
1235   else {
1236     cout<<"F-AliMUONTrackExtrap::GetField: fgkField = 0x0"<<endl;
1237     exit(-1);
1238   }
1239   
1240   Field[0] = b[0]; Field[1] = b[1]; Field[2] = b[2];
1241
1242   return;
1243 }