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