Celeritas
0.5.0-86+4a8eea4
|
Sample the energy of the delta ray for muon or hadron ionization. More...
#include <BetheBlochEnergyDistribution.hh>
Public Types | |
Type aliases | |
using | Energy = units::MevEnergy |
using | Mass = units::MevMass |
Public Member Functions | |
CELER_FUNCTION | BetheBlochEnergyDistribution (ParticleTrackView const &particle, Energy electron_cutoff, Mass electron_mass) |
Construct with incident and exiting particle data. | |
template<class Engine > | |
CELER_FUNCTION Energy | operator() (Engine &rng) |
CELER_FUNCTION Energy | min_secondary_energy () const |
Minimum energy of the secondary electron [MeV]. | |
CELER_FUNCTION Energy | max_secondary_energy () const |
Maximum energy of the secondary electron [MeV]. | |
template<class Engine > | |
CELER_FUNCTION auto | operator() (Engine &rng) -> Energy |
Sample secondary electron energy. More... | |
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] \]
and factorized as
\[ \difd{\sigma}{T} = C f(T) g(T) \]
with \( T \in [T_{cut}, T_{max}] \), where \( f(T) = \frac{1}{T^2} \), \( g(T) = 1 - \beta^2 \frac{T}{T_max} + s \frac{T^2}{2 E^2} \), \( T \) is the kinetic energy of the electron, \( E \) is the total energy of the incident particle, and \( s \) is 0 for spinless particles and 1 otherwise. The energy is sampled from \( f(T) \) and accepted with probability \( g(T) \).
CELER_FUNCTION auto celeritas::BetheBlochEnergyDistribution::operator() | ( | Engine & | rng | ) | -> Energy |