Geant4 11.4.0
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
47
49{
50
51public:
52
54 const G4String& nam = "DNABornIonisationModel");
55
57
60
61 void Initialise(const G4ParticleDefinition*, const G4DataVector& = *(new G4DataVector())) override;
62
64 const G4ParticleDefinition* p,
65 G4double ekin,
66 G4double emin,
67 G4double emax) override;
68
69 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
71 const G4DynamicParticle*,
72 G4double tmin,
73 G4double maxEnergy) override;
74
76 G4int /*level*/,
78 G4double /*kineticEnergy*/) override;
79 void StartTracking(G4Track*) override;
80
81 G4double DifferentialCrossSection(G4ParticleDefinition * aParticleDefinition, G4double k, G4double energyTransfer, G4int shell);
82
84 G4double incomingParticleEnergy, G4int shell, G4double random);
85
86 inline void SelectFasterComputation(G4bool input);
87
88 inline void SelectStationary(G4bool input);
89
90 inline void SelectSPScaling(G4bool input);
91
92protected:
93
95
96private:
97
98 G4bool fasterCode{false};
99 G4bool statCode{false};
100 G4bool spScaling{true};
101
102 // Water density table
103 const std::vector<G4double>* fpMolWaterDensity;
104
105 // Deexcitation manager to produce fluo photons and e-
106 G4VAtomDeexcitation* fAtomDeexcitation;
107
108 const G4Track* fTrack{nullptr};
109 G4DNAChemistryManager* fChemistry{nullptr};
110
111
112 std::map<G4String,G4double,std::less<G4String> > lowEnergyLimit;
113 std::map<G4String,G4double,std::less<G4String> > highEnergyLimit;
114
115 G4bool isInitialised{false};
116 G4int verboseLevel;
117
118 // Cross section
119 using MapFile = std::map<G4String, G4String, std::less<G4String>>;
120 MapFile tableFile; // useful ?
121
122 using MapData = std::map<G4String, G4DNACrossSectionDataSet *, std::less<G4String>>;
123 MapData tableData;
124
125 // Final state
126
127 G4DNAWaterIonisationStructure waterStructure;
128
129 G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
130
131 G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
132
133 G4double Interpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
134
135 G4double QuadInterpolator( G4double e11,
136 G4double e12,
137 G4double e21,
138 G4double e22,
139 G4double x11,
140 G4double x12,
141 G4double x21,
142 G4double x22,
143 G4double t1,
144 G4double t2,
145 G4double t,
146 G4double e);
147
148 using TriDimensionMap = std::map<G4double, std::map<G4double, G4double>>;
149
150 TriDimensionMap eDiffCrossSectionData[6];
151 TriDimensionMap eNrjTransfData[6]; // for cumulated dcs
152
153 TriDimensionMap pDiffCrossSectionData[6];
154 TriDimensionMap pNrjTransfData[6]; // for cumulated dcs
155
156 std::vector<G4double> eTdummyVec;
157 std::vector<G4double> pTdummyVec;
158
159 using VecMap = std::map<G4double, std::vector<G4double>>;
160
161 VecMap eVecm;
162 VecMap pVecm;
163
164 VecMap eProbaShellMap[6]; // for cumulated dcs
165 VecMap pProbaShellMap[6]; // for cumulated dcs
166
167 // Partial cross section
168 G4int RandomSelect(G4double energy,const G4String& particle );
169};
170
172{
173 fasterCode = input;
174}
175
176//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
177
179{
180 statCode = input;
181}
182
183//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
184
186{
187 spScaling = input;
188}
189
190//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
191
192#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
void StartTracking(G4Track *) override
G4double TransferedEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell, G4double random)
G4ParticleChangeForGamma * fParticleChangeForGamma
G4DNABornIonisationModel1(const G4ParticleDefinition *p=nullptr, const G4String &nam="DNABornIonisationModel")
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67