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