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