ATLAS Offline Software
EnergyLossUpdator.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 // EnergyLossUpdator.cxx, (c) ATLAS Detector software
8 
9 // Trk include
13 // Gaudi
14 #include "GaudiKernel/MsgStream.h"
15 #include "GaudiKernel/SystemOfUnits.h"
16 
17 
18 #include <cmath>
19 
20 
21 // constructor
23  const std::string& n,
24  const IInterface* p)
25  : AthAlgTool(t, n, p)
26  , m_detailedEloss(true)
27  , m_optimalRadiation(true)
28 {
29  declareInterface<Trk::IEnergyLossUpdator>(this);
30  // scale for the most probable value
31  declareProperty("DetailedEloss", m_detailedEloss);
32  declareProperty("OptimalRadiation", m_optimalRadiation);
33 }
34 
35 // public interface method
36 double
38  double p,
40 {
42  return 0.;
43  }
44 
45  // preparation of kinetic constants
47  double E = std::sqrt(p * p + m * m);
48  double beta = p / E;
49  double gamma = E / m;
50 
51  // add ionization and radiation
52  double dEdX =
55 
56  // add e+e- pair production and photonuclear effect for muons at energies
57  // above 8 GeV
58  if ((particle == Trk::muon) && (E > 8000.)) {
59  if (E < 1.e6) {
60  dEdX += -0.5345 / mat.x0() + 6.803e-5 * E / mat.x0() +
61  2.278e-11 * E * E / mat.x0() -
62  9.899e-18 * E * E * E / mat.x0(); // E below 1 TeV
63  } else {
64  dEdX += -2.986 / mat.x0() + 9.253e-5 * E / mat.x0(); // E above 1 TeV
65  }
66  }
67  return dEdX;
68 }
69 
70 // public interface method
73  double p,
74  double pathcorrection,
77  bool useMPV) const
78 {
79  if (particle == Trk::undefined) {
81  "undefined ParticleHypothesis, energy loss calculation cancelled");
82  return {};
83  }
84 
85  if (useMPV) {
86  return ionizationEnergyLoss(mat, p, pathcorrection, dir, particle);
87  }
88 
89  double deltaE = 0.;
90  // preparation
91  double sign = (dir == Trk::oppositeMomentum) ? -1. : 1.;
92 
93  double pathLength = pathcorrection * mat.thicknessInX0() * mat.x0();
94 
95  double sigIoni = 0.;
96  double sigRad = 0.;
97  double kazL = 0.;
99  p, (mat.material()), particle, sigIoni, kazL);
101  p, (mat.material()), particle, sigRad);
102 
103  meanIoni = sign * pathLength * meanIoni;
104  meanRad = sign * pathLength * meanRad;
105  sigIoni = pathLength * sigIoni;
106  sigRad = pathLength * sigRad;
107  kazL = pathLength * kazL;
108 
109  //
110  // include pathlength dependence of Landau ionization
111  //
112  sigIoni = sigIoni - kazL * std::log(pathLength);
113 
114  deltaE = meanIoni + meanRad;
115 
116  double sigmaDeltaE = std::sqrt(sigIoni * sigIoni + sigRad * sigRad);
117  ATH_MSG_DEBUG(" Energy loss updator deltaE "
118  << deltaE << " meanIoni " << meanIoni << " meanRad " << meanRad
119  << " sigIoni " << sigIoni << " sigRad " << sigRad << " sign "
120  << sign << " pathLength " << pathLength);
121  return (!m_detailedEloss ? Trk::EnergyLoss(deltaE, sigmaDeltaE)
122  : Trk::EnergyLoss(deltaE, sigmaDeltaE, sigmaDeltaE,
123  sigmaDeltaE, meanIoni, sigIoni,
124  meanRad, sigRad, pathLength));
125 }
126 
127 // public interface method
128 
131  double caloEnergy,
132  double caloEnergyError,
133  double pCaloEntry,
134  double momentumError,
135  int& elossFlag) const
136 {
137  //
138  // Input: the detailed EnergyLoss object in the Calorimeter that contains the
139  // different Eloss terms and their uncertainties; caloEnergy and error; and
140  // the muon momentumError (all in MeV)
141  //
142  // For use in the MuonSystem
143  // Input: caloEnergy = 0. caloEnergyError = 0. and pCaloEntry = pMuonEntry
144  // momentum at MuonEntry
145  //
146  // Output: an updated Energy loss values deltaE()
147  // that can be used in the track fit and corresponds to the Most
148  // Probable EnergyLoss value taking into account the ionization,
149  // radiation and smearing due to the errors including the
150  // momentumError (in MeV)
151  //
152  // elossFlag = false if Calorimeter Energy is NOT stored (and later
153  // fitted) on the Eloss object
154  // = true Calorimeter Energy is stored and will be fitted
155  //
156  // deltaE is used in the final fit
157  //
158 
159  elossFlag = 0;
160 
161  int isign = 1;
162  if (eLoss.deltaE() < 0) {
163  isign = -1;
164  }
165 
166  double deltaE = eLoss.deltaE();
167  double sigmaDeltaE = eLoss.sigmaDeltaE();
168  // Detailed Eloss
169  double deltaE_ioni = eLoss.meanIoni();
170  double sigmaDeltaE_ioni = 0.45 * eLoss.sigmaIoni(); // sigma Landau
171  double deltaE_rad = eLoss.meanRad();
172  double sigmaDeltaE_rad =
173  eLoss.sigmaRad(); // rms and mean of steep exponential
174  double depth = eLoss.length();
175 
176  // Eloss radiative protection
177 
178  if (eLoss.meanRad() > 100000.) {
179  deltaE_rad = 100000.;
180  sigmaDeltaE_rad = eLoss.sigmaRad() * 100000. / eLoss.meanRad();
181  }
182 
183  double sigmaPlusDeltaE = eLoss.sigmaPlusDeltaE();
184  double sigmaMinusDeltaE = eLoss.sigmaMinusDeltaE();
185 
186  double MOP = deltaE_ioni - isign * 3.59524 * sigmaDeltaE_ioni;
187 
188  //
189  // MOP shift due to ionization and radiation
190  //
191  double MOPshift =
192  isign * 50 * 10000. / pCaloEntry +
193  isign * 0.75 * std::sqrt(sigmaDeltaE_ioni * sigmaDeltaE_rad);
194  double MOPshiftNoRad = isign * 50 * 10000. / pCaloEntry;
195  //
196  // define sigmas for Landau convoluted with exponential
197  //
198  double fracErad = sigmaDeltaE_rad +
199  std::abs(deltaE_rad) * pCaloEntry / (800000. + pCaloEntry);
200  double sigmaL = sigmaDeltaE_ioni + 0.8 * fracErad / 3.59524;
201  // double sigmaLNoRad = sigmaDeltaE_ioni;
202  double sigmaMinus = 1.02 * sigmaDeltaE_ioni + 0.08 * sigmaDeltaE_rad;
203  double sigmaPlus = 4.65 * sigmaDeltaE_ioni + 1.16 * sigmaDeltaE_rad;
204  // double sigmaMinusNoRad = 1.02 * sigmaDeltaE_ioni;
205  // double sigmaPlusNoRad = 4.65 * sigmaDeltaE_ioni;
206  double xc = momentumError / (sigmaL > 0. ? sigmaL : 1.);
207  double correction =
208  (0.3849 * xc * xc + 7.76672e-03 * xc * xc * xc) /
209  (1 + 2.8631 * xc + 0.3849 * xc * xc + 7.76672e-03 * xc * xc * xc);
210 
211  //
212  // Case where the measured Calorimeter energy is not available (e.g. low pT or
213  // not isolated)
214  //
215 
216  if (caloEnergyError <= 0) {
217  //
218  // Shift of MOP due to momentum resolution
219  //
220  double MOPreso = isign * 3.59524 * sigmaL * correction;
221 
222  deltaE = MOP + MOPshift + MOPreso;
223  sigmaMinusDeltaE = sigmaMinus;
224  sigmaPlusDeltaE = sigmaPlus;
225  sigmaDeltaE = std::sqrt(0.5 * sigmaMinusDeltaE * sigmaMinusDeltaE +
226  0.5 * sigmaPlusDeltaE * sigmaPlusDeltaE);
227  //
228  if (m_optimalRadiation && std::abs(deltaE) < caloEnergy &&
229  pCaloEntry > 100000) {
230  //
231  // Calorimeter measurement can be used as veto to say there was no
232  // significant radiation
233  //
234  // In that case the Eloss is taken as the ionization Eloss
235  // Use MOP after correction for landau tail (MOPshiftNoRad) and momentum
236  // resolution smearing (MOPreso)
237  //
238  sigmaL = sigmaDeltaE_ioni + 0.3 * fracErad / 3.59524;
239  xc = momentumError / (sigmaL > 0. ? sigmaL : 1.);
240  correction =
241  (0.3849 * xc * xc + 7.76672e-03 * xc * xc * xc) /
242  (1 + 2.8631 * xc + 0.3849 * xc * xc + 7.76672e-03 * xc * xc * xc);
243 
244  MOPreso = isign * 3.59524 * sigmaL * correction;
245  deltaE = MOP + MOPshift + MOPreso;
246  sigmaMinusDeltaE = sigmaMinus;
247  sigmaPlusDeltaE = sigmaPlus;
248  sigmaDeltaE = std::sqrt(0.5 * sigmaMinusDeltaE * sigmaMinusDeltaE +
249  0.5 * sigmaPlusDeltaE * sigmaPlusDeltaE);
250  }
251  } else {
252  double sigmaPlusTot =
253  std::sqrt(sigmaPlus * sigmaPlus + caloEnergyError * caloEnergyError);
254  if (m_optimalRadiation) {
255  sigmaPlusTot =
256  std::sqrt(4.65 * sigmaDeltaE_ioni * 4.65 * sigmaDeltaE_ioni +
257  caloEnergyError * caloEnergyError);
258  }
259  double MOPtot = std::abs(MOP + MOPshift);
260  if (m_optimalRadiation) {
261  MOPtot = std::abs(MOP + MOPshiftNoRad);
262  }
263 
264  if (caloEnergy > MOPtot + 2 * sigmaPlusTot) {
265  //
266  // Use measured Calorimeter energy
267  //
268  //
269  // take into account the tail in the Measured Eloss
270  //
271  double MOPreso = isign * 3.59524 * sigmaL * correction;
272  deltaE = isign * caloEnergy + MOPreso;
273  sigmaMinusDeltaE = caloEnergyError + 0.08 * sigmaDeltaE_rad;
274  sigmaPlusDeltaE = caloEnergyError + 1.16 * sigmaDeltaE_rad;
275  sigmaDeltaE = std::sqrt(0.5 * sigmaMinusDeltaE * sigmaMinusDeltaE +
276  0.5 * sigmaPlusDeltaE * sigmaPlusDeltaE);
277  elossFlag = 1;
278  } else {
279  // Use MOP after corrections
280 
281  //
282  // Shift of MOP due to momentum resolution smearing
283  //
284  sigmaL = sigmaDeltaE_ioni + 0.3 * fracErad / 3.59524;
285  xc = momentumError / (sigmaL > 0. ? sigmaL : 1);
286  correction =
287  (0.3849 * xc * xc + 7.76672e-03 * xc * xc * xc) /
288  (1 + 2.8631 * xc + 0.3849 * xc * xc + 7.76672e-03 * xc * xc * xc);
289  double MOPreso = isign * 3.59524 * sigmaL * correction;
290  //
291  // Use MOP after correction for landau tail (MOPshiftNoRad) and radiation
292  // (MOPshift) and momentum resolution smearing (MOPreso)
293  //
294  deltaE = MOP + MOPshift + MOPreso;
295  sigmaMinusDeltaE = sigmaMinus;
296  sigmaPlusDeltaE = sigmaPlus;
297  sigmaDeltaE = std::sqrt(0.5 * sigmaMinusDeltaE * sigmaMinusDeltaE +
298  0.5 * sigmaPlusDeltaE * sigmaPlusDeltaE);
299  }
300  }
301 
302  return {deltaE,
303  sigmaDeltaE,
304  sigmaMinusDeltaE,
305  sigmaPlusDeltaE,
306  deltaE_ioni,
307  sigmaDeltaE_ioni,
308  deltaE_rad,
309  sigmaDeltaE_rad,
310  depth};
311 }
312 
313 // public interface method
314 void
316  double eta,
317  double phi,
318  double& X0Scale,
319  double& ElossScale) const
320 {
321  //
322  // for Calorimeter icalo = 1
323  // Muon System icalo = 0
324  // convention eta, phi is at Calorimeter Exit (or Muon Entry)
325  // eta and phi are from the position (not direction)
326  //
327  // input X0 and ElossScale = 1
328  // output updated X0Scale and ElossScale
329  //
330 
331  double X0CaloGirder[4] = {
332  -1.02877e-01, -2.74322e-02, 8.12989e-02, 9.73551e-01
333  };
334 
335  // R2012 to R2015 determined on ATLAS-R2-2015-03-01-00 to be used in rel 21
336  constexpr double X0CaloScale[60] = {
337  1.01685, 1.02092, 1.01875, 1.01812, 1.01791, 1.01345, 1.01354,
338  1.02145, 1.01645, 1.01585, 1.0172, 1.02262, 1.01464, 0.990931,
339  0.971953, 0.99845, 1.01433, 0.982143, 0.974015, 0.978742, 0.960029,
340  0.966766, 0.980199, 0.989586, 0.997144, 1.00169, 0.994166, 0.966332,
341  0.93671, 0.935656, 0.921994, 0.901489, 0.897799, 0.89638, 0.905629,
342  0.903374, 0.925922, 0.941203, 0.956273, 0.968618, 0.976883, 0.988349,
343  0.99855, 1.00212, 1.01456, 1.01541, 1.02532, 1.03238, 1.03688,
344  1.03783, 1.02078, 1.01529, 1.0156, 1.02212, 1.02226, 1.02406,
345  1.02188, 1.00661, 1.00661, 1.00661
346  };
347 
348  // R2012 to R2015 determined on ATLAS-R2-2015-03-01-00 to be used in rel 21
349  constexpr double ElossCaloScale[30] = {
350  1.06921, 1.06828, 1.06734, 1.06092, 1.06638, 1.06335, 1.07421, 1.05885,
351  1.07351, 1.07435, 1.06902, 1.07704, 1.08782, 1.09844, 1.115, 1.07609,
352  1.08233, 1.08764, 1.08209, 1.08255, 1.08008, 1.07573, 1.077, 1.07271,
353  1.07343, 1.07769, 1.07794, 1.08377, 1.08377, 1.08377
354  };
355 
356  //
357  constexpr double X0MuonScale[60] = {
358  -0.0320612, -0.0320612, -0.0320612, -0.0320612, -0.0693796, -0.0389677,
359  -0.0860891, -0.124606, -0.0882329, -0.100014, -0.0790912, -0.0745538,
360  -0.099088, -0.0933711, -0.0618782, -0.0619762, -0.0658361, -0.109704,
361  -0.129547, -0.143364, -0.0774768, -0.0739859, -0.0417835, -0.022119,
362  0.00308797, 0.0197657, -0.0137871, -0.036848, -0.0643794, -0.0514949,
363  -0.0317105, 0.016539, 0.0308435, -0.00056883, -0.00756813, -0.00760612,
364  -0.0234571, -0.0980915, -0.101175, -0.102354, -0.0920337, -0.100337,
365  -0.0887628, 0.0660931, 0.228999, 0.260675, 0.266301, 0.267907,
366  0.281668, 0.194433, 0.132954, 0.20707, 0.220466, 0.20936,
367  0.191441, 0.191441, 0.191441, 0.191441, 0.191441, 0.191441
368  };
369 
370  int i60 = std::abs(eta) * 20.;
371 
372  if (i60 < 0) {
373  i60 = 0;
374  }
375  if (i60 > 59) {
376  i60 = 59;
377  }
378 
379  if (icalo == 1) {
380  //
381  // Girder parametrization
382  //
383  double x =
384  phi + 3.1416 - 3.1416 / 32. * int((3.1416 + phi) / (3.1416 / 32.));
385  double scale = 0.;
386  if (x > M_PI / 64.) {
387  x = M_PI / 32. - x;
388  }
389 
390  if (x < 0.005) {
391  scale =
392  X0CaloGirder[0] * (1 - x / 0.005) + X0CaloGirder[1] + X0CaloGirder[3];
393  } else if (x < 0.017) {
394  scale = X0CaloGirder[1] + X0CaloGirder[3];
395  } else if (x < 0.028) {
396  scale = X0CaloGirder[2] + X0CaloGirder[3];
397  } else {
398  scale = X0CaloGirder[3];
399  }
400 
401  if (std::abs(eta) > 1.3) {
402  scale = 1.;
403  }
404  //
405  // eta dependence of X0
406  //
407  scale *= X0CaloScale[i60];
408  X0Scale = scale;
409  //
410  // eta dependence of Eloss
411  //
412  int i30 = std::abs(eta) * 10.;
413  if (i30 < 0) {
414  i30 = 0;
415  }
416  if (i30 > 29) {
417  i30 = 29;
418  }
419 
420  double nfactor = 0.987363 / 1.05471;
421 
422  ElossScale = nfactor * ElossCaloScale[i30] * scale;
423  } else {
424  //
425  // Muon system
426  //
427  // eta dependence of X0
428  //
429  double scale = 1. + X0MuonScale[i60];
430  //
431  // Muon scale is now 1 with MuonTrackingGeometry and TrkDetDescrGeoModelCnv
432  // fixes
433  //
434  scale = 1.0;
435  X0Scale = scale;
436  ElossScale = 0.93 * scale;
437  }
438 }
439 
442  double p,
443  double pathcorrection,
446 {
447  // preparation
448  double sign = (dir == Trk::oppositeMomentum) ? -1. : 1.;
449  double pathLength = pathcorrection * mat.thicknessInX0() * mat.x0();
450 
451  double sigIoni = 0.;
452  double kazL = 0.;
453 
454  double meanIoni =
456  p, (mat.material()), particle, sigIoni, kazL, pathLength);
457 
458  return (!m_detailedEloss
459  ? Trk::EnergyLoss(meanIoni, sigIoni)
460  : Trk::EnergyLoss(meanIoni,
461  sigIoni,
462  sigIoni,
463  sigIoni,
464  meanIoni,
465  sigIoni,
466  0.,
467  0.,
468  pathLength));
469 }
Trk::EnergyLoss::sigmaMinusDeltaE
double sigmaMinusDeltaE() const
returns the negative side
EnergyLoss.h
Trk::MaterialInteraction::dEdl_ionization
static double dEdl_ionization(double p, const Material &mat, ParticleHypothesis particle, double &sigma, double &kazL)
dE/dl ionization energy loss per path unit
Definition: MaterialInteraction.cxx:117
Trk::ParticleSwitcher::particle
constexpr ParticleHypothesis particle[PARTICLEHYPOTHESES]
the array of masses
Definition: ParticleHypothesis.h:76
egammaParameters::depth
@ depth
pointing depth of the shower as calculated in egammaqgcld
Definition: egammaParamDefs.h:276
python.SystemOfUnits.m
int m
Definition: SystemOfUnits.py:91
python.PerfMonSerializer.p
def p
Definition: PerfMonSerializer.py:743
Trk::EnergyLossUpdator::m_optimalRadiation
bool m_optimalRadiation
use calorimeter more optimal for radiation detection
Definition: EnergyLossUpdator.h:188
Trk::EnergyLossUpdator::dEdX
virtual double dEdX(const MaterialProperties &mat, double p, ParticleHypothesis particle=pion) const override final
dEdX calculation when providing MaterialProperties, a momentum, a pathlength, and a ParicleHypothesis...
Definition: EnergyLossUpdator.cxx:37
CaloCellPos2Ntuple.int
int
Definition: CaloCellPos2Ntuple.py:24
MaterialProperties.h
eta
Scalar eta() const
pseudorapidity method
Definition: AmgMatrixBasePlugin.h:79
Trk::oppositeMomentum
@ oppositeMomentum
Definition: PropDirection.h:21
AthCommonDataStore< AthCommonMsg< AlgTool > >::declareProperty
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T > &t)
Definition: AthCommonDataStore.h:145
Trk::MaterialInteraction::dE_MPV_ionization
static double dE_MPV_ionization(double p, const Trk::Material &mat, Trk::ParticleHypothesis particle, double &sigma, double &kazL, double path)
Most Propable dE ionization energly loss.
Definition: MaterialInteraction.cxx:205
mat
GeoMaterial * mat
Definition: LArDetectorConstructionTBEC.cxx:53
Trk::EnergyLoss::sigmaDeltaE
double sigmaDeltaE() const
returns the symmatric error
M_PI
#define M_PI
Definition: ActiveFraction.h:11
Trk::EnergyLossUpdator::energyLoss
virtual EnergyLoss energyLoss(const MaterialProperties &mat, double p, double pathcorrection, PropDirection dir=alongMomentum, ParticleHypothesis particle=pion, bool useMPV=false) const override final
deltaE calculation using dEdX and integrating along pathlength, assuming constant dEdX during for the...
Definition: EnergyLossUpdator.cxx:72
Trk::EnergyLoss::sigmaRad
double sigmaRad() const
Trk::undefined
@ undefined
Definition: ParticleHypothesis.h:38
read_hist_ntuple.t
t
Definition: read_hist_ntuple.py:5
Trk::EnergyLoss::meanIoni
double meanIoni() const
yodamerge_tmp.scale
scale
Definition: yodamerge_tmp.py:138
Trk::ParticleHypothesis
ParticleHypothesis
Definition: ParticleHypothesis.h:25
Trk::EnergyLossUpdator::m_detailedEloss
bool m_detailedEloss
provide extended EnergyLoss info
Definition: EnergyLossUpdator.h:187
tools.zlumi_mc_cf.correction
def correction(mu, runmode, campaign, run=None)
Definition: zlumi_mc_cf.py:4
Trk::EnergyLoss::length
double length() const
Trk::PropDirection
PropDirection
Definition: PropDirection.h:19
Trk::EnergyLossUpdator::EnergyLossUpdator
EnergyLossUpdator(const std::string &, const std::string &, const IInterface *)
AlgTool like constructor.
Definition: EnergyLossUpdator.cxx:22
TrigVtx::gamma
@ gamma
Definition: TrigParticleTable.h:26
beamspotman.n
n
Definition: beamspotman.py:731
ATH_MSG_DEBUG
#define ATH_MSG_DEBUG(x)
Definition: AthMsgStreamMacros.h:29
pdg_comparison.kazL
kazL
Definition: pdg_comparison.py:325
Trk::EnergyLoss::deltaE
double deltaE() const
returns the
sign
int sign(int a)
Definition: TRT_StrawNeighbourSvc.h:127
Trk::muon
@ muon
Definition: ParticleHypothesis.h:28
Trk::EnergyLoss::meanRad
double meanRad() const
beamspotman.dir
string dir
Definition: beamspotman.py:623
Trk::MaterialInteraction::dEdl_radiation
static double dEdl_radiation(double p, const Material &mat, ParticleHypothesis particle, double &sigma)
dE/dl radiation energy loss per path unit
Definition: MaterialInteraction.cxx:229
Trk::ParticleMasses::mass
constexpr double mass[PARTICLEHYPOTHESES]
the array of masses
Definition: ParticleHypothesis.h:53
Trk::EnergyLossUpdator::updateEnergyLoss
virtual EnergyLoss updateEnergyLoss(EnergyLoss &eLoss, double caloEnergy, double caloEnergyError, double pCaloEntry, double momentumError, int &elossFlag) const override final
Method to recalculate Eloss values for the fit setting an elossFlag using as an input the detailed El...
Definition: EnergyLossUpdator.cxx:130
VP1PartSpect::E
@ E
Definition: VP1PartSpectFlags.h:21
Trk::nonInteracting
@ nonInteracting
Definition: ParticleHypothesis.h:25
Trk::EnergyLoss
This class describes energy loss material effects in the ATLAS tracking EDM.
Definition: EnergyLoss.h:34
Trk::MaterialProperties
Definition: MaterialProperties.h:40
DiTauMassTools::MaxHistStrategyV2::e
e
Definition: PhysicsAnalysis/TauID/DiTauMassTools/DiTauMassTools/HelperFunctions.h:26
EnergyLossUpdator.h
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
python.CaloCondTools.log
log
Definition: CaloCondTools.py:20
Trk::MaterialInteraction::dEdXBetheBloch
static double dEdXBetheBloch(const Trk::MaterialProperties &mat, double beta, double gamma, Trk::ParticleHypothesis particle)
dE/dl ionization energy loss per path unit
Definition: MaterialInteraction.cxx:165
Trk::phi
@ phi
Definition: ParamDefs.h:81
Trk::MaterialInteraction::dEdXBetheHeitler
static double dEdXBetheHeitler(const Trk::MaterialProperties &mat, double initialE, Trk::ParticleHypothesis particle)
Definition: MaterialInteraction.cxx:266
AthAlgTool
Definition: AthAlgTool.h:26
Trk::x
@ x
Definition: ParamDefs.h:61
MuonParameters::beta
@ beta
Definition: MuonParamDefs.h:144
Trk::EnergyLossUpdator::ionizationEnergyLoss
Trk::EnergyLoss ionizationEnergyLoss(const MaterialProperties &mat, double p, double pathcorrection, PropDirection dir=alongMomentum, ParticleHypothesis particle=pion) const
Definition: EnergyLossUpdator.cxx:441
Trk::EnergyLoss::sigmaPlusDeltaE
double sigmaPlusDeltaE() const
returns the positive side
Trk::EnergyLossUpdator::getX0ElossScales
virtual void getX0ElossScales(int icalo, double eta, double phi, double &X0Scale, double &ElossScale) const override final
Routine to calculate X0 and Eloss scale factors for the Calorimeter and Muon System.
Definition: EnergyLossUpdator.cxx:315
Trk::EnergyLoss::sigmaIoni
double sigmaIoni() const