Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNABornIonisationModel1.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 G4DNABornIonisationModel1_h
29#define G4DNABornIonisationModel1_h 1
30
31#include "G4VEmModel.hh"
34
36#include "G4Electron.hh"
37#include "G4Proton.hh"
39
41
44#include "G4NistManager.hh"
45
46
48{
49
50public:
51
53 const G4String& nam = "DNABornIonisationModel");
54
56
59
60 void Initialise(const G4ParticleDefinition*, const G4DataVector& = *(new G4DataVector())) override;
61
63 const G4ParticleDefinition* p,
64 G4double ekin,
65 G4double emin,
66 G4double emax) override;
67
68 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
70 const G4DynamicParticle*,
71 G4double tmin,
72 G4double maxEnergy) override;
73
75 G4int /*level*/,
77 G4double /*kineticEnergy*/) override;
78
79 G4double DifferentialCrossSection(G4ParticleDefinition * aParticleDefinition, G4double k, G4double energyTransfer, G4int shell);
80
82 G4double incomingParticleEnergy, G4int shell, G4double random) ;
83
84 inline void SelectFasterComputation(G4bool input);
85
86 inline void SelectStationary(G4bool input);
87
88 inline void SelectSPScaling(G4bool input);
89
90protected:
91
93
94private:
95
96 G4bool fasterCode;
97 G4bool statCode;
98 G4bool spScaling;
99
100 // Water density table
101 const std::vector<G4double>* fpMolWaterDensity;
102
103 // Deexcitation manager to produce fluo photons and e-
104 G4VAtomDeexcitation* fAtomDeexcitation;
105
106 std::map<G4String,G4double,std::less<G4String> > lowEnergyLimit;
107 std::map<G4String,G4double,std::less<G4String> > highEnergyLimit;
108
109 // TODO :
110// std::map<const G4ParticleDefinition*,std::pair<G4double,G4double> > fEnergyLimits;
111
112
113 G4bool isInitialised{false};
114 G4int verboseLevel;
115
116 // Cross section
117
118 using MapFile = std::map<G4String, G4String, std::less<G4String>>;
119 MapFile tableFile; // useful ?
120
121 using MapData = std::map<G4String, G4DNACrossSectionDataSet *, std::less<G4String>>;
122 MapData tableData;
123
124 // Final state
125
126 G4DNAWaterIonisationStructure waterStructure;
127
128 G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
129
130 G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
131
132 G4double Interpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
133
134 G4double QuadInterpolator( G4double e11,
135 G4double e12,
136 G4double e21,
137 G4double e22,
138 G4double x11,
139 G4double x12,
140 G4double x21,
141 G4double x22,
142 G4double t1,
143 G4double t2,
144 G4double t,
145 G4double e);
146
147 using TriDimensionMap = std::map<G4double, std::map<G4double, G4double>>;
148
149 TriDimensionMap eDiffCrossSectionData[6];
150 TriDimensionMap eNrjTransfData[6]; // for cumulated dcs
151
152 TriDimensionMap pDiffCrossSectionData[6];
153 TriDimensionMap pNrjTransfData[6]; // for cumulated dcs
154
155 std::vector<G4double> eTdummyVec;
156 std::vector<G4double> pTdummyVec;
157
158 using VecMap = std::map<G4double, std::vector<G4double>>;
159
160 VecMap eVecm;
161 VecMap pVecm;
162
163 VecMap eProbaShellMap[6]; // for cumulated dcs
164 VecMap pProbaShellMap[6]; // for cumulated dcs
165
166 // Partial cross section
167
168 G4int RandomSelect(G4double energy,const G4String& particle );
169
170};
171
173{
174 fasterCode = input;
175}
176
177//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
178
180{
181 statCode = input;
182}
183
184//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
185
187{
188 spScaling = input;
189}
190
191//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
192
193#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 &= *(new G4DataVector())) override
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
G4DNABornIonisationModel1(const G4DNABornIonisationModel1 &)=delete
G4double GetPartialCrossSection(const G4Material *, G4int, const G4ParticleDefinition *, G4double) override
G4double DifferentialCrossSection(G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
G4DNABornIonisationModel1 & operator=(const G4DNABornIonisationModel1 &right)=delete
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double TransferedEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell, G4double random)
G4ParticleChangeForGamma * fParticleChangeForGamma
G4DNABornIonisationModel1(const G4ParticleDefinition *p=nullptr, const G4String &nam="DNABornIonisationModel")