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