EM Physics¶
The physics models in Celeritas are primarily derived from references cited by Geant4, including the Geant4 physics reference manual. Undocumented adjustments to those models in Geant4 may also be implemented.
Processes and models¶
The following table summarizes the EM processes and models in Celeritas.
Particle |
Processes |
Models |
Celeritas Implementation |
Applicability |
---|---|---|---|---|
\(e^-\) |
Ionization |
Møller |
0–100 TeV |
|
Bremsstrahlung |
Seltzer–Berger |
0–1 GeV |
||
Relativistic |
1 GeV – 100 TeV |
|||
Coulomb scattering |
Urban |
|
100 eV – 100 TeV |
|
Coulomb |
0–100 TeV |
|||
\(e^+\) |
Ionization |
Bhabha |
0–100 TeV |
|
Bremsstrahlung |
Seltzer-Berger |
0–1 GeV |
||
Relativistic |
1 GeV – 100 TeV |
|||
Coulomb scattering |
Urban |
|
100 eV – 100 TeV |
|
Coulomb |
0–100 TeV |
|||
Annihilation |
\(e^+,e^- \to 2\gamma\) |
0–100 TeV |
||
\(\gamma\) |
Photoelectric |
Livermore |
0–100 TeV |
|
Compton scattering |
Klein–Nishina |
0–100 TeV |
||
Pair production |
Bethe–Heitler |
0–100 TeV |
||
Rayleigh scattering |
Livermore |
0–100 TeV |
||
\(\mu^-\) |
Ionization |
ICRU73QO |
0–200 keV |
|
Bethe–Bloch |
200 keV–1 GeV |
|||
Mu Bethe–Bloch |
200 keV–100 TeV |
|||
Bremsstrahlung |
Mu bremsstrahlung |
0–100 TeV |
||
\(\mu^+\) |
Ionization |
Bragg |
0–200 keV |
|
Bethe–Bloch |
200 keV–1 GeV |
|||
Mu Bethe–Bloch |
200 keV–100 TeV |
|||
Bremsstrahlung |
Mu bremsstrahlung |
0–100 TeV |
The implemented physics models are meant to match the defaults constructed in
G4EmStandardPhysics
. Known differences are:
Particles other than electrons, positrons, and gammas are not currently supported.
As with the AdePT project, Celeritas currently extends the range of Urban MSC to higher energies rather than implementing the Wentzel-VI and discrete Coulomb scattering.
Celeritas imports tracking cutoffs and other parameters from
G4EmParameters
, but some custom model cutoffs are not accessible to Celeritas.
As extension to the various random distributions, Celeritas expresses many physics operations as distributions of updated track states based on original track states. For example, the Tsai-Urban distribution used for sampling exiting angles of bremsstrahlung and pair production has parameters of incident particle energy and mass, and it samples the exiting polar angle cosine.
All discrete interactions (in Geant4 parlance, “post-step do-it”s) use distributions to sample an Interaction based on incident particle properties. The sampled result contains the updated particle direction and energy, as well as properties of any secondary particles produced.
Ionization¶
-
class MollerBhabhaInteractor¶
Perform Moller (e-e-) and Bhabha (e+e-) scattering.
This interaction, part of the ionization process, is when an incident electron or positron ejects an electron from the surrounding matter.
Note
This performs the same sampling routine as in Geant4’s G4MollerBhabhaModel class, as documented in section 10.1.4 of the Geant4 Physics Reference (release 10.6).
-
template<class EnergySampler>
class MuHadIonizationInteractor¶ Perform the discrete part of the muon or hadron ionization process.
This simulates the production of delta rays by incident muons or hadrons. The same basic sampling routine is used by multiple models, but the energy of the secondary is sampled from a distribution unique to the model.
Note
This performs the same sampling routine as in Geant4’s
G4BetheBlochModel
,G4MuBetheBlochModel
,G4BraggModel
, andG4ICRU73QOModel
, as documented in the Geant4 Physics Reference Manual release 11.2 sections 11.1 and 12.1.5.
The exiting energy distribution from most of these ionization models are sampled using external helper distributions.
-
class BetheBlochEnergyDistribution¶
Sample the energy of the delta ray for muon or hadron ionization.
This samples the energy according to the Bethe-Bloch model, as described in the Geant4 Physics Reference Manual release 11.2 section 12.1.5. The Bethe-Bloch differential cross section can be written as
\[ \difd{\sigma}{T} = 2\pi r_e^2 mc^2 Z \frac{z_p^2}{\beta^2}\frac{1}{T^2} \left[1 - \beta^2 \frac{T}{T_{max}} + s \frac{T^2}{2E^2} \right] \]\[ \difd{\sigma}{T} = C f(T) g(T) \]
-
class BraggICRU73QOEnergyDistribution¶
Sample the energy of the delta ray for muon or hadron ionization.
This samples the energy according to the Bragg and ICRU73QO models, as described in the Geant4 Physics Reference Manual release 11.2 section 11.1.
-
class BhabhaEnergyDistribution¶
Helper class for
MollerBhabhaInteractor
.Sample the exiting energy fraction for Bhabha scattering.
-
class MollerEnergyDistribution¶
Helper class for
MollerBhabhaInteractor
.Sample the exiting energy fraction for Moller scattering.
-
class MuBBEnergyDistribution¶
Sample delta ray energy for the muon Bethe-Bloch ionization model.
8 This samples the energy according to the muon Bethe-Bloch model, as described in the Geant4 Physics Reference Manual release 11.2 section 11.1. At the higher energies for which this model is applied, leading radiative corrections are taken into account. The differential cross section can be written as
\[ \sigma(E, \epsilon) = \sigma_{BB}(E, \epsilon)\left[1 + \frac{\alpha}{2\pi} \log \left(1 + \frac{2\epsilon}{m_e} \log \left(\frac{4 m_e E(E - \epsilon}{m_{\mu}^2(2\epsilon + m_e)} \right) \right) \right]. \]As in the Bethe-Bloch model, the energy is sampled by factorizing the cross section as \( \sigma = C f(T) g(T) \), where \( f(T) = \frac{1}{T^2} \) and \( T \in [T_{cut}, T_{max}] \). The energy is sampled from \( f(T) \) and accepted with probability \( g(T) \).
Bremsstrahlung¶
-
class RelativisticBremInteractor¶
Perform a high-energy Bremsstrahlung interaction.
This is a relativistic Bremsstrahlung model for high-energy (> 1 GeV) electrons and positrons.
Note
This performs the same sampling routine as in Geant4’s G4eBremsstrahlungRelModel class, as documented in section 10.2.2 of the Geant4 Physics Reference (release 10.7).
-
class SeltzerBergerInteractor¶
Seltzer-Berger model for electron and positron bremsstrahlung processes.
Given an incoming electron or positron of sufficient energy (as per CutOffView), this class provides the energy loss of these particles due to radiation of photons in the field of a nucleus. This model improves accuracy using cross sections based on interpolation of published tables from Seltzer and Berger given in Nucl. Instr. and Meth. in Phys. Research B, 12(1):95–134 (1985) and Atomic Data and Nuclear Data Tables, 35():345 (1986). The cross sections are obtained from SBEnergyDistribution and are appropriately scaled in the case of positrons via SBPositronXsCorrector.
Note
This interactor performs an analogous sampling as in Geant4’s G4SeltzerBergerModel, documented in 10.2.1 of the Geant Physics Reference (release 10.6). The implementation is based on Geant4 10.4.3.
-
class MuBremsstrahlungInteractor¶
Perform muon bremsstrahlung interaction.
This is a model for the Bremsstrahlung process for muons. Given an incident muon, the class computes the change to the incident muon direction and energy, and it adds a single secondary gamma to the secondary stack.
Note
This performs the same sampling routine as in Geant4’s G4MuBremsstrahlungModel class, as documented in section 11.2 of the Geant4 Physics Reference (release 10.6).
The Seltzer–Berger interactions are sampled with the help of an energy distribution and cross section correction:
-
template<class XSCorrector>
class SBEnergyDistribution¶ Sample exiting photon energy from Bremsstrahlung.
The SB energy distribution uses tabulated Seltzer-Berger differential cross section data from G4EMLOW, which stores scaled cross sections as a function of incident particle energy and exiting gamma energy (see SeltzerBergerModel for details). The sampling procedure is roughly laid out in section [PHYS341] of the GEANT3 physics reference manual, although like Geant4 we use raw tabulated SB data rather than a parameter fit. Also like Geant4 we include the extra density correction factor.
Once an element and incident energy have been selected, the exiting energy distribution is the differential cross section, which is stored as a scaled tabulated value. The reconstructed cross section gives the pdf
\[ p(\kappa) \propto \difd{\sigma}{k} s \propto \frac{1}{\kappa} \chi_Z(E, \kappa) \]\[ p_1(\kappa) \frac{1}{\ln 1/\kappa_c} \frac{1}{\kappa} \]\[ \kappa = \exp( \xi \ln \kappa_c ) \,, \]\[ p_2(\kappa) = \frac{\chi_Z(E, \kappa)}{\max_\kappa \chi_Z(E, \kappa)} \]The above algorithm is idealized; in practice, the minimum and maximum values are adjusted for a constant factor \(d_\rho\), which depends on the incident particle mass + energy and the material’s electron density:
\[ \kappa_\mathrm{min} = \ln (E_c^2 + d_\rho E^2) \]\[ \kappa_\mathrm{max} = \ln (E^2 + d_\rho E^2) \, . \]\[ k = \sqrt{ \exp(\kappa_\mathrm{min} + \xi [\kappa_\mathrm{max} - \kappa_\mathrm{min}]) - d_\rho E^2} \]Most of the mechanics of the sampling are in the template-free
SBEnergyDistHelper
, which is passed as a construction argument to this sampler. The separate class exists here to minimize duplication of templated code, which is required to provide for an on-the-fly correction of the cross section sampling.
-
class SBPositronXsCorrector¶
Scale SB differential cross sections for positrons.
This correction factor is the ratio of the positron to electron cross sections. It appears in the bowels of
G4SeltzerBergerModel::SampleEnergyTransfer
in Geant4 and scales the SB cross section by a factor\[ \frac{\sigma^+}{\sigma^-} = \exp( 2 \pi \alpha Z [ \beta^{-1}(E - k_c) - \beta^{-1}(E - k) ] ) \]\[ \beta^{-1}(E) = \frac{c}{v} = \sqrt{1 - \left( \frac{m_e c^2}{E + m_e c^2} \right)^2} \,, \]\( \frac{\sigma^+}{\sigma^-} = 1 \) at the low end of the spectrum (where \( k / E = 0 \)) and is 0 at the tip of the spectrum where \( k / E \to 1 \).
The correction factor is described in:
Kim, Longhuan, R. H. Pratt, S. M. Seltzer, and M. J. Berger. “Ratio of Positron to Electron Bremsstrahlung Energy Loss: An Approximate Scaling Law.” Physical Review A 33, no. 5 (May 1, 1986): 3002–9. https://doi.org/10.1103/PhysRevA.33.3002.
A simple distribution is used to sample exiting polar angles from electron bremsstrahlung (and gamma conversion).
-
class TsaiUrbanDistribution¶
Polar angular distribution for pair-production and bremsstrahlung processes.
For pair production, the polar angle of the electron (or positron) is defined with respect to the direction of the parent photon. The energy- angle distribution given by Tsai is quite complicated to sample and can be approximated by a density function suggested by Urban.
The angular distribution of the emitted photons is obtained from a simplified formula based on the Tsai cross-section, which is expected to become isotropic in the low energy limit.
Note
This performs the same sampling routine as in Geant4’s ModifiedTsai class, based on derivation from Tsai (Rev Mod Phys 49,421(1977) and documented in section 6.5.2 (pair-production), and 10.2.1 and 10.2.4 (bremsstrahlung) of the Geant4 Physics Reference (release 10.6).
Relativistic bremsstrahlung and relativistic Bethe-Heitler sampling both use a helper class to calculate LPM factors.
-
class LPMCalculator¶
Calculate the Landau-Pomeranchuk-Migdal (LPM) suppression functions.
The LPM effect is the suppression of low-energy photon production due to electron multiple scattering. At high energies and in high density materials, the cross sections for pair production and bremsstrahlung are reduced. The differential cross sections accounting for the LPM effect are expressed in terms of the LPM suppression functions \( xi(s) \), \( G(s) \), and \( \phi(s) \).
Here \( \epsilon \) is the ratio of the electron (or positron) energy to the photon energy, \( \epsilon = E / k \).
For small energies, the suppression factors all approach unity.
See section 10.2.2 of the Geant4 Physics Reference Manual and ComputeLPMfunctions and GetLPMFunctions in G4eBremsstrahlungRelModel and G4PairProductionRelModel. Also see T. Stanev, Ch. Vankov, Development of ultrahigh-energy electromagnetic cascades in water and lead including the Landau-Pomeranchuk-Migdal effect, Phys. Rev. D, 25 (1982), p. 1291.
Muon bremsstrahlung calculates the differential cross section as part of rejection sampling.
-
class MuBremsDiffXsCalculator¶
Calculate the differential cross section for muon bremsstrahlung.
The differential cross section can be written as
\[ \difd{\sigma}{\epsilon} = \frac{16}{3} \alpha N_A (\frac{m}{\mu} r_e)^2 \frac{1}{\epsilon A} Z(Z \Phi_n + \Phi_e) (1 - v + \frac{3}{4} v^2), \]The contribution to the cross section from the nucleus is given by
\[ \Phi_n = \ln \frac{B Z^{-1/3} (\mu + \delta(D'_n \sqrt{e} - 2))}{D'_n (m + \delta \sqrt{e} B Z^{-1/3})} \ , \]The contribution to the cross section from electrons is given by
\[ \Phi_e = \ln \frac{B' Z^{-2/3} \mu}{\left(1 + \frac{\delta \mu}{m^2 \sqrt{e}}\right)(m + \delta \sqrt{e} B' Z^{-2/3})} \ . \]The constants \( B \) and \( B' \) were calculated using the Thomas-Fermi model. In the case of hydrogen, where the Thomas-Fermi model does not serve as a good approximation, the exact values of the constants were calculated analytically.
This performs the same calculation as in Geant4’s
G4MuBremsstrahlungModel::ComputeDMicroscopicCrossSection()
and as described in section 11.2.1 of the Physics Reference Manual. The formulae are taken mainly from SR Kelner, RP Kokoulin, and AA Petrukhin. About cross section for high-energy muon bremsstrahlung. Technical Report, MEphI, 1995. Preprint MEPhI 024-95, Moscow, 1995, CERN SCAN-9510048.
Muon bremsstrahlung and pair production use a simple distribution to sample the exiting polar angles.
Warning
doxygenclass: Cannot find class “celeritas::MuBremsPPAngularDistribution” in doxygen xml output for project “celeritas” from directory: /home/runner/work/celeritas/celeritas/build/doc/xml
Photon scattering¶
-
class KleinNishinaInteractor¶
Perform Compton scattering, neglecting atomic binding energy.
This is a model for the discrete Compton inelastic scattering process. Given an incident gamma, it adds a single secondary (electron) to the secondary stack and returns an interaction for the change to the incident gamma direction and energy. No cutoffs are performed for the incident energy or the exiting gamma energy. A secondary production cutoff is applied to the outgoing electron.
Note
This performs the same sampling routine as in Geant4’s G4KleinNishinaCompton, as documented in section 6.4.2 of the Geant4 Physics Reference (release 10.6).
-
class RayleighInteractor¶
Apply the Livermore model of Rayleigh scattering to photons.
Note
This performs the same sampling routine as in Geant4’s G4LivermoreRayleighModel class, as documented in section 6.2.2 of the Geant4 Physics Reference (release 10.6).
Conversion/annihilation/photoelectric¶
-
class BetheHeitlerInteractor¶
Relativistic model for electron-positron pair production.
The energies of the secondary electron and positron are sampled using the Bethe-Heitler cross sections with a Coulomb correction. The LPM effect is taken into account for incident gamma energies above 100 GeV. Exiting particle directions are sampled with the
TsaiUrbanDistribution
. Note that energy is not exactly conserved.Note
This performs the same sampling routine as in Geant4’s G4PairProductionRelModel, as documented in sections 6.5 (gamma conversion) and 10.2.2 (LPM effect) of the Geant4 Physics Reference Manual (release 10.7)
-
class EPlusGGInteractor¶
Annihilate a positron to create two gammas.
This is a model for the discrete positron-electron annihilation process which simulates the in-flight annihilation of a positron with an atomic electron and produces into two photons. It is assumed that the atomic electron is initially free and at rest.
Note
This performs the same sampling routine as in Geant4’s G4eeToTwoGammaModel class, as documented in section 10.3 of the Geant4 Physics Reference (release 10.6). The maximum energy for the model is specified in Geant4.
-
class LivermorePEInteractor¶
Livermore model for the photoelectric effect.
A parameterization of the photoelectric cross sections in two different energy intervals, formulated as
\[ \sigma(E) = a_1 / E + a_2 / E^2 + a_3 / E^3 + a_4 / E^4 + a_5 / E^5 + a_6 / E^6 \: , \]Note
This performs the same sampling routine as in Geant4’s G4LivermorePhotoElectricModel class, as documented in section 6.3.5 of the Geant4 Physics Reference (release 10.6).
-
class AtomicRelaxation¶
Simulate particle emission from atomic deexcitation.
The EADL radiative and non-radiative transition data is used to simulate the emission of fluorescence photons and (optionally) Auger electrons given an initial shell vacancy created by a primary process.
Positron annihilation and Livermore photoelectric cross sections are calculated on the fly (as opposed to pre-tabulated cross sections).
-
class EPlusGGMacroXsCalculator¶
Calculates the macroscopic cross section of positron annihilation.
The Heitler formula (section 10.3.2 of the Geant4 Physics Reference Manual, Release 10.6) is used to compute the macroscopic cross section for positron annihilation on the fly at the given energy.
-
class LivermorePEMicroXsCalculator¶
Calculate photoelectric effect cross sections using the Livermore data.
Coulomb scattering¶
Elastic scattering of charged particles off atoms can be simulated in three ways:
A detailed single scattering model in which each scattering interaction is sampled
A multiple scattering approach which calculates global effects from many collisions
A combination of the two
Though it is the most accurate, the single Coulomb scattering model is too computationally expensive to be used in most applications as the number of collisions can be extremely large. Instead, a “condensed” simulation algorithm is typically used to determine the net energy loss, displacement, and direction change from many collisions after a given path length. The Urban model is the default multiple scattering model in Celeritas for all energies and in Geant4 below 100 MeV. A third “mixed” simulation approach uses multiple scattering to simulate interactions with scattering angles below a given polar angle limit and single scattering for large angles. The Wentzel VI model, used together with the single Coulomb scattering model, is an implementation of the mixed simulation algorithm. It is the default model in Geant4 above 100 MeV and currently under development in Celeritas.
-
class CoulombScatteringInteractor¶
Applies the Wentzel single Coulomb scattering model.
This models incident high-energy electrons and positrons elastically scattering off of nuclei and atomic electrons. Scattering off of the nucleus versus electrons is randomly sampled based on the relative cross-sections. No secondaries are created in this process (in the future, with hadronic transport support, secondary ions may be emitted), however production cuts are used to determine the maximum scattering angle off of electrons.
Note
This performs the same sampling as in Geant4’s G4eCoulombScatteringModel, as documented in section 8.2 of the Geant4 Physics Reference Manual (release 11.1).
-
class WentzelDistribution¶
Helper class for
CoulombScatteringInteractor
.Samples the polar scattering angle for the Wentzel Coulomb scattering model.
This chooses between sampling scattering off an electron or nucleus based on the relative cross sections. Electron scattering angle imposes a maximum scattering angle (see WentzelHelper::cos_thetamax_electron), and nuclear sattering rejects an angular change based on the Mott cross section (see MottRatioCalculator). Nuclear scattering depends on the electronic and nuclear cross sections (calculated by
WentzelHelper
) and nuclear form factors (seeExpNuclearFormFactor
,GaussianNuclearFormFactor
, andUUNuclearFormFactor
) that describe an approximate spatial charge distribution of the nucleus.The polar angle distribution is given in [Fern] eqn 88 and is normalized on the interval \( cos\theta \in [\cos\theta_\mathrm{min}, \cos\theta_\mathrm{max}] \). The sampling function for \( \mu = \frac{1}{2}(1 - \cos\theta) \) is
\[ \mu = \mu_1 + \frac{(A + \mu_1) \xi (\mu_2 - \mu_1)}{A + \mu_2 - \xi (\mu_2 - \mu_1)}, \]References: [Fern] J.M. Fernandez-Varea, R. Mayol and F. Salvat. On the theory and simulation of multiple elastic scattering of electrons. Nucl. Instrum. and Method. in Phys. Research B, 73:447-473, Apr 1993. doi:10.1016/0168-583X(93)95827-R [PRM] Geant4 Physics Reference Manual (Release 11.1) sections 8.2 and 8.5
-
class MottRatioCalculator¶
Calculates the ratio of Mott cross section to the Rutherford cross section.
This ratio is an adjustment of the cross section from a purely classical treatment of a point nucleus in an electronic cloud (Rutherford scattering) to a quantum mechanical treatment. The implementation is an interpolated approximation developed in [LQZ95] and described in the Geant Physics Reference Manual PRM section 8.4
[LQZ95] T. Lijian, H. Quing and L. Zhengming, Radiat. Phys. Chem. 45 (1995), 235-245
The parameter
cos_theta
is the cosine of the scattered angle in the z-aligned momentum frame.
-
class ExpNuclearFormFactor : public celeritas::NuclearFormFactorTraits¶
Exponential nuclear form factor.
This nuclear form factor corresponds
NuclearFormFactorType::exponential
and assumes the nuclear charge decays exponentially from its center. This assumes a parameterization of the atomic nucleus valid for light and medium atomic nuclei (Eq. 7 of [BKM2002]):\[ R_N = 1.27A^{0.27} \,\mathrm{fm} \][BKM2002] A.V. Butkevich, R.P. Kokoulin, G.V. Matushko, S.P. Mikheyev, Comments on multiple scattering of high-energy muons in thick layers, Nuclear Instruments and Methods in Physics Research Section A: Accelerators, Spectrometers, Detectors and Associated Equipment 488 (2002) 282–294. https://doi.org/10.1016/S0168-9002(02)00478-3.
- Todo:
Instead of using this coarse parameterization, we should add nuclear radius to the isotope properties for a more accurate treatment, and construct these classes directly from the atomic radius.
Subclassed by celeritas::GaussianNuclearFormFactor
-
class GaussianNuclearFormFactor : public celeritas::ExpNuclearFormFactor¶
Gaussian nuclear form factor.
This nuclear form factor corresponds
NuclearFormFactorType::gaussian
and assumes a Gaussian distribution of nuclear charge. Its prefactor has the same value as the exponential.
-
class UUNuclearFormFactor : public celeritas::NuclearFormFactorTraits¶
Uniform-uniform folded nuclear form factor.
This nuclear form factor corresponds
NuclearFormFactorType::flat
and assumes a uniform nuclear charge at the center with a smoothly decreasing charge at the surface. This leads to a form factor:\[ F(q) = F'(x(R_0, q)) F'(x(R_1, q)) \]\[ F'(x) = \frac{3}{x^3} ( \sin x - x \cos x) \][LR16] C. Leroy and P.G. Rancoita. Principles of Radiation Interaction in Matter and Detection. World Scientific (Singapore), 4th edition, 2016.
[H56] R.H. Helm, Inelastic and Elastic Scattering of 187-Mev Electrons from Selected Even-Even Nuclei, Phys. Rev. 104 (1956) 1466–1475. https://doi.org/10.1103/PhysRev.104.1466.
[FMS93] J.M. Fernández-Varea, R. Mayol, F. Salvat, Cross sections for elastic scattering of fast electrons and positrons by atoms, Nuclear Instruments and Methods in Physics Research Section B: Beam Interactions with Materials and Atoms 82 (1993) 39–45. https://doi.org/10.1016/0168-583X(93)95079-K.
Warning
This form factor suffers from catastrophic numerical cancellation for small radii and momenta so should only be used for large nuclei or large momentum transfers.
-
class UrbanMscSafetyStepLimit¶
Sample a step limit for the Urban MSC model using the “safety” algorithm.
This distribution is to be used for tracks that have non-negligible steps and are near the boundary. Otherwise, no displacement or step limiting is needed.
Note
This code performs the same method as in ComputeTruePathLengthLimit of G4UrbanMscModel, as documented in section 8.1.6 of the Geant4 10.7 Physics Reference Manual or CERN-OPEN-2006-077 by L. Urban.
-
class UrbanMscScatter¶
Sample angular change and lateral displacement with the Urban multiple scattering model.
Note
This code performs the same method as in G4VMultipleScattering::AlongStepDoIt and G4UrbanMscModel::SampleScattering of the Geant4 10.7 release.
Discrete cross sections¶
Most physics processes use pre-calculated cross sections that are tabulated and interpolated.
-
class XsCalculator¶
Find and interpolate scaled cross sections on a uniform log grid.
This cross section calculator uses the same representation and interpolation as Geant4’s physics tables for EM physics:
The energy grid is uniformly spaced in log(E),
Values greater than or equal to an index i’ are scaled by E,
Linear interpolation between energy points is used to calculate the final value, and
If the energy is at or above the i’ index, the final result is scaled by 1/E.
This scaling and interpolation exactly reproduces functions \( f(E) \sim a E + b \) below the E’ threshold and \( f(E) \sim \frac{a'}{E} + b' \) above that threshold.
Note that linear interpolation is applied with energy points, not log-energy points.
XsCalculator calc_xs(xs_grid, xs_params.reals); real_type xs = calc_xs(particle.energy());
Cross sections for each process are evaluated at the beginning of the step along with range limiters.
-
inline StepLimit celeritas::calc_physics_step_limit(MaterialTrackView const &material, ParticleTrackView const &particle, PhysicsTrackView &physics, PhysicsStepView &pstep)¶
Calculate physics step limits based on cross sections and range limiters.
If undergoing an interaction, the process is sampled from the stored beginning-of-step cross sections.
-
template<class Engine>
ActionId celeritas::select_discrete_interaction(MaterialView const &material, ParticleTrackView const &particle, PhysicsTrackView const &physics, PhysicsStepView &pstep, Engine &rng)¶ Choose the physics model for a track’s pending interaction.
If the interaction MFP is zero, the particle is undergoing a discrete interaction. Otherwise, the result is a false ActionId.
Sample from the previously calculated per-process cross section/decay to determine the interacting process ID.
From the process ID and (post-slowing-down) particle energy, we obtain the applicable model ID.
For energy loss (continuous-discrete) processes, the post-step energy will be different from the pre-step energy, so the assumption that the cross section is constant along the step is no longer valid. Use the “integral
approach” to sample the discrete interaction from the correct probability distribution (section 7.4 of the Geant4 Physics Reference release 10.6).
Continuous slowing down¶
Most charged interactions emit one or more low-energy particles during their interaction. Instead of creating explicit daughter tracks that are immediately killed due to low energy, part of the interaction cross section is lumped into a “slowing down” term that continuously deposits energy locally over the step.
-
inline ParticleTrackView::Energy celeritas::calc_mean_energy_loss(ParticleTrackView const &particle, PhysicsTrackView const &physics, real_type step)¶
Calculate mean energy loss over the given “true” step length.
Stopping power is an integral over low-exiting-energy secondaries. Above some threshold energy T_c we treat exiting secondaries discretely; below it, we lump them into this continuous loss term that varies based on the energy, the atomic number density, and the element number:
\[ \frac{dE}{dx} = N_Z \int_0^{T_c} \frac{d \sigma_Z(E, T)}{dT} T dT \]The stopping range R due to these low-energy processes is an integral over the inverse of the stopping power: basically the distance that will take a particle (without discrete processes at play) from its current energy to an energy of zero.
\[ R = \int_0 ^{E_0} - \frac{dx}{dE} dE . \]Both Celeritas and Geant4 approximate the range limit as the minimum range over all processes, rather than the range as a result from integrating all energy loss processes over the allowed energy range. This is usually not a problem in practice because the range will get automatically decreased by
range_to_step
, and the above range calculation neglects energy loss by discrete processes.Both Geant4 and Celeritas integrate the energy loss terms across processes to get a single energy loss vector per particle type. The range table is an integral of the mean stopping power: the total distance for the particle’s energy to reach zero.
- Todo:
The geant3 manual makes the point that linear interpolation of energy loss rate results in a piecewise constant energy deposition curve, which is why they use spline interpolation. Investigate higher-order reconstruction of energy loss curve, e.g. through spline-based interpolation or log-log interpolation.
Zero energy loss can occur in the following cases:
The energy loss value at the given energy is zero (e.g. high energy particles)
The urban model is selected and samples zero collisions (possible in thin materials and/or small steps)
Note
The inverse range correction assumes range is always the integral of the stopping power/energy loss.
Note
See section 7.2.4 Run Time Energy Loss Computation of the Geant4 physics manual. See also the longer discussions in section 8 of PHYS010 of the Geant3 manual (1993).
Since true energy loss is a stochastic function of many small collisions, the mean energy loss term is an approximation. Additional models are implemented to adjust the loss per step with stochastic sampling for improved accuracy.
-
class EnergyLossHelper¶
Precalculate energy loss fluctuation properties.
Fluctuations in the energy loss of charged particles over a given thickness of material arise from statistical variation in both the number of collisions and the energy lost in each collision. Above a given energy threshold, fluctuations are simulated through the explicit sampling of secondaries. However, the continuous energy loss below the cutoff energy also has fluctuations, and these are not taken into account in the calculation of the mean loss. For continuous energy loss, fluctuation models are used to sample the actual restricted energy loss given the mean loss.
Different models are used depending on the value of the parameter \( \kappa = \xi / T_{max} \), the ratio of the mean energy loss to the maximum allowed energy transfer in a single collision.
For large \( \kappa \), when the particle loses all or most of its energy along the step, the number of collisions is large and the straggling function can be approximated by a Gaussian distribution.
Otherwise, the Urban model for energy loss fluctuations in thin layers is used.
Note
This performs the same sampling routine as in Geant4’s G4UniversalFluctuation, as documented in section 7.3 in the Geant4 Physics Reference Manual and in PHYS332 and PHYS333 in GEANT3, CERN Program Library Long Writeup, W5013 (1993).
-
class EnergyLossGammaDistribution¶
Sample energy loss from a gamma distribution.
This model is valid for heavy particles with small mean losses (less than 2 Bohr standard deviations). It is a special case of
EnergyLossGaussianDistribution
(see that class for more documentation).Note that while this appears in G4UniversalFluctuation, the Geant4 documentation does not explain why the loss is sampled from a gamma distribution in this case.
-
class EnergyLossGaussianDistribution¶
Sample energy loss from a Gaussian distribution.
In a thick absorber, the total energy transfer is a result of many small energy losses from a large number of collisions. The central limit theorem applies, and the energy loss fluctuations can be described by a Gaussian distribution. See section 7.3.1 of the Geant4 Physics Reference Manual and GEANT3 PHYS332 section 2.3.
The Gaussian approximation is valid for heavy particles and in the regime \( \kappa = \xi / T_\textrm{max} > 10 \). Fluctuations of the unrestricted energy loss follow a Gaussian distribution if \( \Delta E > \kappa T_{max} \), where \( T_{max} \) is the maximum energy transfer (PHYS332 section 2). For fluctuations of the restricted energy loss, the condition is modified to \( \Delta E > \kappa T_{c} \) and \( T_{max} \le 2 T_c \), where \( T_c \) is the delta ray cutoff energy (PRM Eq. 7.6-7.7).
-
class EnergyLossUrbanDistribution¶
Sample from the Urban model of energy loss fluctuations in thin layers.
The Urban model is used to compute the energy loss fluctuation when \( \kappa \) is small. It assumes atoms have only two energy levels with binding energies \( E_1 \) and \( E_2 \) and that the particle-atom interaction will either be an excitation with energy loss \( E_1 \) or \( E_2 \) or an ionization with energy loss proportional to \( 1 / E^2 \). The number of collisions for each interaction type has a Poisson distribution with mean proportional to the interaction cross section. The energy loss in a step will be the sum of the energy loss contributions from excitation and ionization.
When the number of ionizations is larger than a given threshold, a fast sampling method can be used instead. The possible energy loss interval is divided into two parts: a lower range in which the number of collisions is large and the energy loss can be sampled from a Gaussian distribution, and an upper range in which the energy loss is sampled for each collision.
See section 7.3.2 of the Geant4 Physics Reference Manual and GEANT3 PHYS332 section 2.4 for details.
Imported data¶
In addition to the core Imported data, these import parameters are used to provide cross sections, setup options, and other data to the EM physics.
-
struct ImportEmParameters¶
Common electromagnetic physics parameters (see G4EmParameters.hh).
Note
Geant4 v11 removed the Spline() option from G4EmParameters.hh.
Note
The Geant4 MSC models use the values in
G4EmParameters
as the defaults; however, the MSC parameters can also be set directly using the model setter methods (there is no way to retrieve the values from the model in that case).
-
struct ImportAtomicTransition¶
EADL transition data for atomic relaxation for a single element.
-
struct ImportAtomicSubshell¶
-
struct ImportAtomicRelaxation¶
-
struct ImportLivermoreSubshell¶
Livermore EPICS2014 photoelectric cross section data for a single element.
-
struct ImportLivermorePE¶
-
struct ImportMuPairProductionTable¶
Sampling table for electron-positron pair production by muons.
This 3-dimensional table is used to sample the energy transfer to the electron-positron pair, \( \epsilon_p \). The outer grid stores the atomic number using 5 equally spaced points in \( \log Z \); the x grid tabulates the ratio \( \log \epsilon_p / T \), where \( T \) is the incident muon energy; the y grid stores the incident muon energy using equal spacing in \( \log T \).
-
using celeritas::ImportSBTable = ImportPhysics2DVector¶
Seltzer Berger differential cross sections for a single element.
This 2-dimensional table stores the scaled bremsstrahlung differential cross section [mb]. The x grid is the log energy of the incident particle [MeV], and the y grid is the ratio of the gamma energy to the incident energy.