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