Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNACPA100IonisationModel.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// CPA100 ionisation model class for electrons
27//
28// Based on the work of M. Terrissol and M. C. Bordage
29//
30// Users are requested to cite the following papers:
31// - M. Terrissol, A. Baudre, Radiat. Prot. Dosim. 31 (1990) 175-177
32// - M.C. Bordage, J. Bordes, S. Edel, M. Terrissol, X. Franceries,
33// M. Bardies, N. Lampe, S. Incerti, Phys. Med. 32 (2016) 1833-1840
34//
35// Authors of this class:
36// M.C. Bordage, M. Terrissol, S. Edel, J. Bordes, S. Incerti
37//
38// 15.01.2014: creation
39//
40
41#ifndef G4DNACPA100IonisationModel_h
42#define G4DNACPA100IonisationModel_h 1
43
44#include "G4VEmModel.hh"
47
49#include "G4Electron.hh"
50#include "G4Proton.hh"
51
53//#include "G4DNACPA100LogLogInterpolation.hh"
54
57#include "G4NistManager.hh"
58
59
61{
62
63public:
64
66 const G4String& nam = "DNACPA100IonisationModel");
67
69
70 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector& = *(new G4DataVector()));
71
72 virtual G4double CrossSectionPerVolume( const G4Material* material,
73 const G4ParticleDefinition* p,
74 G4double ekin,
75 G4double emin,
76 G4double emax);
77
78 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
80 const G4DynamicParticle*,
81 G4double tmin,
82 G4double maxEnergy);
83
84 G4double DifferentialCrossSection(G4ParticleDefinition * aParticleDefinition, G4double k, G4double energyTransfer, G4int shell);
85
86 inline void SelectFasterComputation(G4bool input);
87
88 inline void SelectUseDcs(G4bool input);
89
90 inline void SelectStationary(G4bool input);
91
92protected:
93
95
96private:
97
98 G4bool statCode;
99
100 G4bool fasterCode;
101 G4bool useDcs;
102
103 // Water density table
104 const std::vector<G4double>* fpMolWaterDensity;
105
106 // Deexcitation manager to produce fluo photons and e-
107 G4VAtomDeexcitation* fAtomDeexcitation;
108
109 std::map<G4String,G4double,std::less<G4String> > lowEnergyLimit;
110 std::map<G4String,G4double,std::less<G4String> > highEnergyLimit;
111
112 G4bool isInitialised;
113 G4int verboseLevel;
114
115 // Cross section
116
117 typedef std::map<G4String,G4String,std::less<G4String> > MapFile;
118 MapFile tableFile;
119
120 typedef std::map<G4String,G4DNACrossSectionDataSet*,std::less<G4String> > MapData;
121 MapData tableData;
122
123 // Final state
124
126
127 G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
128
129 G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
130
131 G4double RandomizeEjectedElectronEnergyFromCompositionSampling(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
132
133 G4double RandomTransferedEnergy(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
134
135 void RandomizeEjectedElectronDirection(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4double
136 outgoingParticleEnergy, G4double & cosTheta, G4double & phi );
137
138 G4double Interpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
139
140 G4double QuadInterpolator( G4double e11,
141 G4double e12,
142 G4double e21,
143 G4double e22,
144 G4double x11,
145 G4double x12,
146 G4double x21,
147 G4double x22,
148 G4double t1,
149 G4double t2,
150 G4double t,
151 G4double e);
152
153 typedef std::map<G4double, std::map<G4double, G4double> > TriDimensionMap;
154
155 TriDimensionMap eDiffCrossSectionData[6];
156 TriDimensionMap eNrjTransfData[6]; // for cumulated dcs
157
158 std::vector<G4double> eTdummyVec;
159
160 typedef std::map<G4double, std::vector<G4double> > VecMap;
161
162 VecMap eVecm;
163
164 VecMap eProbaShellMap[6]; // for cumulated dcs
165
166 // Partial cross section
167
168 G4int RandomSelect(G4double energy,const G4String& particle );
169
170 //
171
174
175};
176
177//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
178
180{
181 fasterCode = input;
182}
183
184//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
185
187{
188 useDcs = input;
189}
190
191//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
192
193//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
194
196{
197 statCode = input;
198}
199
200//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
201
202#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4double DifferentialCrossSection(G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
G4ParticleChangeForGamma * fParticleChangeForGamma
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector()))
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)