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