Class in charge of generating positrons from IBD of supernova antineutrino interactions.
More...
Supported by OMSimSNTools, this class generates the corresponding positrons after the inverse beta decay and weights the interactions corresponding to the cross section and the generated volume
|
| OMSimIBD (G4ParticleGun *) |
| Constructor of OMSimIBD class.
|
|
Public Member Functions inherited from SNBaseParticleGenerator |
| SNBaseParticleGenerator (G4ParticleGun *) |
|
void | GeneratePrimaries (G4Event *anEvent) |
|
|
Protected Member Functions inherited from SNBaseParticleGenerator |
void | initialiseDistribution (int column, int targets) |
|
G4double | calculateNeutrinoEnergy () |
|
G4ThreeVector | calculateMomentumDirection (G4double pNeutrinoEnergy) |
|
G4double | calculateWeight (G4double pNuEnergy) |
|
Protected Attributes inherited from SNBaseParticleGenerator |
G4ParticleGun * | m_particleGun |
|
OMSimSNTools | m_SNToolBox |
|
const G4double | m_Gf = 1.166e-5 * 1e-6 / (MeV * MeV) |
|
const G4double | m_Consg = 1.26 |
|
const G4double | m_deltaMass = neutron_mass_c2 - proton_mass_c2 |
|
const G4double | m_massSquaredDifference = (pow(m_deltaMass, 2) - pow(electron_mass_c2, 2)) / 2. |
|
G4double | m_NrTargets |
|
DistributionSampler | m_timeDistribution |
|
DistributionSampler | m_meanEnergyDistribution |
|
DistributionSampler | m_meanEnergySquaredDistribution |
|
◆ OMSimIBD()
OMSimIBD::OMSimIBD |
( |
G4ParticleGun * |
p_gun | ) |
|
- Parameters
-
p_gun | Pointer to the particle p_gun object used for simulation. |
◆ angularDistribution()
This function evaluates the angular distribution of the positron based on the energy of the incident electronic antineutrino. The calculation follows the theoretical approach outlined in "The angular distribution of the reaction \( \bar{\nu}_e + p \rightarrow e^+ + n \)" by P. Vogel and J. F. Beacom (1999), specifically referencing Eq. 14.
- Parameters
-
p_energy | Energy of the incoming electronic antineutrino. |
- Note
- This function populates the xAngularDistribution and yAngularDistribution vectors, which later can be used in other parts of the simulation.
Implements SNBaseParticleGenerator.
◆ calculateSecondaryParticleEnergy()
G4double OMSimIBD::calculateSecondaryParticleEnergy |
( |
G4double |
p_energy, |
|
|
G4double |
p_cosTheta |
|
) |
| |
|
privatevirtual |
Given the energy of the incident electronic antineutrino and the scatter angle, this function calculates the energy of the emitted positron following inverse beta decay. The theoretical basis for this calculation is given in "The angular distribution of the reaction \( \nu_e + p \rightarrow e^+ + n \)" by P. Vogel and J. F. Beacom (1999), specifically referencing Eq. 13.
- Parameters
-
p_energy | Energy of the incoming electronic antineutrino. |
p_cosTheta | Cosine of the scatter angle between the direction of the antineutrino's momentum and the positron's momentum. |
- Returns
- Energy of the emitted positron as a result of the inverse beta decay process. @reference P. Vogel, J. F. Beacom. (1999). The angular distribution of the reaction ( \nu_e + p \rightarrow e^+ + n ). Phys.Rev., D60, 053003
Implements SNBaseParticleGenerator.
◆ initialiseParticle()
void OMSimIBD::initialiseParticle |
( |
| ) |
|
|
privatevirtual |
◆ totalCrossSection()
G4double OMSimIBD::totalCrossSection |
( |
G4double |
p_energy | ) |
|
|
privatevirtual |
This function estimates the total cross-section of the inverse beta decay, which can be used to weigh each event. The theoretical basis for this calculation is presented in "Future detection
of supernova neutrino burst and explosion mechanism" by T. Totani, K. Sato, H. E. Dalhed, and J. R. Wilson (1998), specifically referencing Equation 9.
- Parameters
-
p_energy | Energy of the incoming electronic antineutrino. |
- Returns
- Total cross-section for the given energy. @reference T. Totani, K. Sato, H. E. Dalhed, J. R. Wilson, "Future detection of supernova neutrino burst
and explosion mechanism", Astrophys. J. 496, 1998, 216–225, Preprint astro-ph/9710203, 1998.
Implements SNBaseParticleGenerator.
The documentation for this class was generated from the following files: