Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNABornIonisationModel2.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 G4DNABornIonisationModel2_h
29#define G4DNABornIonisationModel2_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
81 G4double TransferedEnergy(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell, G4double random) ;
82
83 inline void SelectFasterComputation(G4bool input);
84
85 inline void SelectStationary(G4bool input);
86
87 inline void SelectSPScaling(G4bool input);
88
89protected:
90
92
93private:
94
95 G4bool fasterCode;
96 G4bool statCode;
97 G4bool spScaling;
98
99 // Water density table
100 const std::vector<G4double>* fpMolWaterDensity;
101
102 // Deexcitation manager to produce fluo photons and e-
103 G4VAtomDeexcitation* fAtomDeexcitation;
104
105 G4double fLowEnergyLimit;
106 G4double fHighEnergyLimit;
107
108 const G4ParticleDefinition* fParticleDef;
109
110 G4bool isInitialised{false};
111 G4int verboseLevel;
112
113 // Cross section
114 G4String fTableFile;
115 G4DNACrossSectionDataSet* fTableData;
116
117 // Final state
118
119 G4DNAWaterIonisationStructure waterStructure;
120
121 G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
122
123 G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
124
125 G4double Interpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
126
127 G4double QuadInterpolator( G4double e11,
128 G4double e12,
129 G4double e21,
130 G4double e22,
131 G4double x11,
132 G4double x12,
133 G4double x21,
134 G4double x22,
135 G4double t1,
136 G4double t2,
137 G4double t,
138 G4double e);
139
140 using TriDimensionMap = std::map<G4double, std::map<G4double, G4double>>;
141
142 TriDimensionMap fDiffCrossSectionData[6];
143 TriDimensionMap fNrjTransfData[6]; // for cumulated dcs
144
145 std::vector<G4double> fTdummyVec;
146
147 using VecMap = std::map<G4double, std::vector<G4double>>;
148
149 VecMap fVecm;
150 VecMap fProbaShellMap[6]; // for cumulated dcs
151
152 // Partial cross section
153
154 G4int RandomSelect(G4double energy);
155
156};
157
158//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
159
161{
162 fasterCode = input;
163}
164
165//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
166
168{
169 statCode = input;
170}
171
172//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
173
175{
176 spScaling = input;
177}
178
179//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
180
181#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4DNABornIonisationModel2(const G4DNABornIonisationModel2 &)=delete
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector())) override
G4double GetPartialCrossSection(const G4Material *, G4int, const G4ParticleDefinition *, G4double) override
G4ParticleChangeForGamma * fParticleChangeForGamma
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double DifferentialCrossSection(G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
G4DNABornIonisationModel2(const G4ParticleDefinition *p=nullptr, const G4String &nam="DNABornIonisationModel")
G4DNABornIonisationModel2 & operator=(const G4DNABornIonisationModel2 &right)=delete
G4double TransferedEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell, G4double random)