OMSim
Geant4 for IceCube optical module studies
Loading...
Searching...
No Matches
OMSimSNParticleGenerators.hh
Go to the documentation of this file.
1
7#pragma once
8
9#include "G4VUserPrimaryGeneratorAction.hh"
10#include "OMSimPrimaryGeneratorMessenger.hh"
11#include "OMSimSNTools.hh"
12
13#include <globals.hh>
14#include <vector>
15#include <G4PhysicalConstants.hh>
16#include <G4SystemOfUnits.hh>
17
18class G4ParticleGun;
19class G4Event;
20
26class SNBaseParticleGenerator : public G4VUserPrimaryGeneratorAction
27{
28public:
29
30 SNBaseParticleGenerator(G4ParticleGun*);
32
33 void GeneratePrimaries(G4Event* anEvent);
34
35protected:
36 G4ParticleGun* m_particleGun;
37 OMSimSNTools m_SNToolBox;
38
39 void initialiseDistribution(int column, int targets);
40 G4double calculateNeutrinoEnergy();
41 G4ThreeVector calculateMomentumDirection(G4double pNeutrinoEnergy);
42 G4double calculateWeight(G4double pNuEnergy);
43
44 const G4double m_Gf = 1.166e-5 * 1e-6 / (MeV * MeV);
45 const G4double m_Consg = 1.26;
46 const G4double m_deltaMass = neutron_mass_c2 - proton_mass_c2;
47 const G4double m_massSquaredDifference = (pow(m_deltaMass, 2) - pow(electron_mass_c2, 2)) / 2.;
48
49 G4double m_NrTargets;
50 DistributionSampler m_timeDistribution;
51 DistributionSampler m_meanEnergyDistribution;
52 DistributionSampler m_meanEnergySquaredDistribution;
53
54 virtual void initialiseParticle() = 0;
55 virtual DistributionSampler angularDistribution(G4double Enu) = 0;
56 virtual G4double calculateSecondaryParticleEnergy(G4double nu_energy, G4double costheta) = 0;
57 virtual G4double totalCrossSection(G4double energy) = 0;
58};
59
60
61
73{
74public:
75
76 OMSimENES(G4ParticleGun*);
77 ~OMSimENES(){};
78
79private:
80 void initialiseParticle();
82 G4double calculateSecondaryParticleEnergy(G4double nu_energy, G4double costheta);
83 G4double totalCrossSection(G4double energy);
84};
85
86
87
98{
99public:
100
101 OMSimIBD(G4ParticleGun*);
102 ~OMSimIBD(){};
103
104private:
105 void initialiseParticle();
107 G4double calculateSecondaryParticleEnergy(G4double nu_energy, G4double costheta);
108 G4double totalCrossSection(G4double energy);
109};
110
111
112#endif
113
114
Provides the definition for the OMSimSNTools class, centralizing utility methods for Supernova (SN) n...
Utility class for sampling from a given distribution using the inverse cumulative function and interp...
Definition OMSimSNTools.hh:24
Class in charge of generating electrons from Electron-electron_Neutrino Elastic Scattering (ENES) nu_...
Definition OMSimSNParticleGenerators.hh:73
DistributionSampler angularDistribution(G4double Enu)
Computes the angular distribution for electron scattering based on incident neutrino energy.
Definition OMSimSNParticleGenerators.cc:177
G4double totalCrossSection(G4double energy)
Computes the total cross section for ENES based on incident neutrino energy.
Definition OMSimSNParticleGenerators.cc:240
G4double calculateSecondaryParticleEnergy(G4double nu_energy, G4double costheta)
Computes the electron energy from elastic scattering based on incident neutrino energy and scatter an...
Definition OMSimSNParticleGenerators.cc:223
Class in charge of generating positrons from IBD of supernova antineutrino interactions.
Definition OMSimSNParticleGenerators.hh:98
DistributionSampler angularDistribution(G4double Enu)
Computes the angular distribution of the reaction ( \bar{\nu}_e + p \rightarrow e^+ + n ).
Definition OMSimSNParticleGenerators.cc:285
G4double totalCrossSection(G4double energy)
Calculates the total cross-section of the inverse beta decay reaction for a given energy.
Definition OMSimSNParticleGenerators.cc:367
G4double calculateSecondaryParticleEnergy(G4double nu_energy, G4double costheta)
Computes the energy of a positron resulting from the inverse beta decay.
Definition OMSimSNParticleGenerators.cc:343
Provides utility methods for generating ENES and IBD interactions.
Definition OMSimSNTools.hh:56
Definition OMSimSNParticleGenerators.hh:27