Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNADingfelderChargeDecreaseModel.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27
28#ifndef G4DNADingfelderChargeDecreaseModel_h
29#define G4DNADingfelderChargeDecreaseModel_h 1
30
31#include "G4VEmModel.hh"
34
35#include "G4Proton.hh"
37#include "G4NistManager.hh"
38
40{
41
42public:
43
45 const G4String& nam = "DNADingfelderChargeDecreaseModel");
46
48
51
52 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
53
55 const G4ParticleDefinition* p,
56 G4double ekin,
57 G4double emin,
58 G4double emax) override;
59
60 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
62 const G4DynamicParticle*,
63 G4double tmin,
64 G4double maxEnergy) override;
65
66 inline void SelectStationary(G4bool input);
67
68protected:
69
71
72private:
73
74 G4bool statCode;
75
76 // Water density table
77 const std::vector<G4double>* fpMolWaterDensity = nullptr;
78
79 std::map<G4String,G4double,std::less<G4String> > lowEnergyLimit;
80 std::map<G4String,G4double,std::less<G4String> > highEnergyLimit;
81
82 G4bool isInitialised{false};
83 G4int verboseLevel;
84
85 // Partial cross section
86
87 G4double PartialCrossSection(G4double energy, G4int level, const G4ParticleDefinition* particle);
88
89 G4double Sum(G4double energy, const G4ParticleDefinition* particle);
90
91 G4int RandomSelect(G4double energy, const G4ParticleDefinition* particle);
92
93 G4int numberOfPartialCrossSections[3]; // 3 is the particle type index
94
95 G4double f0[2][3];
96 G4double a0[2][3];
97 G4double a1[2][3];
98 G4double b0[2][3];
99 G4double b1[2][3];
100 G4double c0[2][3];
101 G4double d0[2][3];
102 G4double x0[2][3];
103 G4double x1[2][3];
104
105 // Final state
106
107 G4int NumberOfFinalStates(G4ParticleDefinition* particleDefinition, G4int finalStateIndex);
108
109 G4ParticleDefinition* OutgoingParticleDefinition(G4ParticleDefinition* particleDefinition, G4int finalStateIndex);
110
111 G4double WaterBindingEnergyConstant(G4ParticleDefinition* particleDefinition, G4int finalStateIndex);
112
113 G4double OutgoingParticleBindingEnergyConstant(G4ParticleDefinition* particleDefinition, G4int finalStateIndex);
114
115 // Reusable particle definitions
116 G4ParticleDefinition* protonDef = nullptr;
117 G4ParticleDefinition* alphaPlusPlusDef = nullptr;
118 G4ParticleDefinition* alphaPlusDef = nullptr;
119 G4ParticleDefinition* hydrogenDef = nullptr;
120 G4ParticleDefinition* heliumDef = nullptr;
121
122};
123
124//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
125
127{
128 statCode = input;
129}
130
131//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
132
133#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4DNADingfelderChargeDecreaseModel(const G4DNADingfelderChargeDecreaseModel &)=delete
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
G4DNADingfelderChargeDecreaseModel & operator=(const G4DNADingfelderChargeDecreaseModel &right)=delete
~G4DNADingfelderChargeDecreaseModel() override=default
G4DNADingfelderChargeDecreaseModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="DNADingfelderChargeDecreaseModel")