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#include "G4RandomTools.hh"
13
14#include <globals.hh>
15#include <vector>
16#include <G4PhysicalConstants.hh>
17#include <G4SystemOfUnits.hh>
18
19class G4ParticleGun;
20class G4Event;
21
27class SNBaseParticleGenerator : public G4VUserPrimaryGeneratorAction
28{
29public:
30
31 SNBaseParticleGenerator(G4ParticleGun*);
33
34 void GeneratePrimaries(G4Event* anEvent);
35
36protected:
37 G4ParticleGun* m_particleGun;
38 OMSimSNTools m_SNToolBox;
39
40 void initialiseDistribution(int column, int targets);
41 G4double calculateNeutrinoEnergy();
42 G4ThreeVector calculateMomentumDirection(G4double pNeutrinoEnergy);
43 G4double calculateWeight(G4double pNuEnergy);
44
45 const G4double m_Gf = 1.166e-5 * 1e-6 / (MeV * MeV);
46 const G4double m_Consg = 1.26;
47 const G4double m_deltaMass = neutron_mass_c2 - proton_mass_c2;
48 const G4double m_massSquaredDifference = (pow(m_deltaMass, 2) - pow(electron_mass_c2, 2)) / 2.;
49
50 G4double m_NrTargets;
51 DistributionSampler m_timeDistribution;
52 DistributionSampler m_meanEnergyDistribution;
53 DistributionSampler m_meanEnergySquaredDistribution;
54
55 virtual void initialiseParticle() = 0;
56 virtual DistributionSampler angularDistribution(G4double Enu) = 0;
57 virtual G4double calculateSecondaryParticleEnergy(G4double nu_energy, G4double costheta) = 0;
58 virtual G4double totalCrossSection(G4double energy) = 0;
59};
60
61
62
74{
75public:
76
77 OMSimENES(G4ParticleGun*);
78 ~OMSimENES(){};
79
80private:
81 void initialiseParticle();
83 G4double calculateSecondaryParticleEnergy(G4double nu_energy, G4double costheta);
84 G4double totalCrossSection(G4double energy);
85};
86
87
88
99{
100public:
101
102 OMSimIBD(G4ParticleGun*);
103 ~OMSimIBD(){};
104
105private:
106 void initialiseParticle();
108 G4double calculateSecondaryParticleEnergy(G4double nu_energy, G4double costheta);
109 G4double totalCrossSection(G4double energy);
110};
111
112
113
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:74
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:99
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:28