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