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