Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MicroElecInelasticModel.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// G4MicroElecInelasticModel.hh, 2011/08/29 A.Valentin, M. Raine
28//
29// Based on the following publications
30//
31// - Inelastic cross-sections of low energy electrons in silicon
32// for the simulation of heavy ion tracks with theGeant4-DNA toolkit,
33// NSS Conf. Record 2010, pp. 80-85
34// - Geant4 physics processes for microdosimetry simulation:
35// very low energy electromagnetic models for electrons in Si,
36// NIM B, vol. 288, pp. 66-73, 2012.
37// - Geant4 physics processes for microdosimetry simulation:
38// very low energy electromagnetic models for protons and
39// heavy ions in Si, NIM B, vol. 287, pp. 124-129, 2012.
40//
41//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
42
43#ifndef G4MicroElecInelasticModel_h
44#define G4MicroElecInelasticModel_h 1
45
46
47#include "globals.hh"
48#include "G4VEmModel.hh"
51
53#include "G4Electron.hh"
54#include "G4Proton.hh"
55#include "G4GenericIon.hh"
57
59
62#include "G4NistManager.hh"
63
65{
66
67public:
68
70 const G4String& nam = "MicroElecInelasticModel");
71
73
74 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
75
76 virtual G4double CrossSectionPerVolume( const G4Material* material,
77 const G4ParticleDefinition* p,
78 G4double ekin,
79 G4double emin,
80 G4double emax);
81
82 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
84 const G4DynamicParticle*,
85 G4double tmin,
86 G4double maxEnergy);
87
88 double DifferentialCrossSection(G4ParticleDefinition * aParticleDefinition, G4double k, G4double energyTransfer, G4int shell);
89
91 G4double incomingParticleEnergy, G4int shell, G4double random) ;
92
93 inline void SelectFasterComputation(G4bool input);
94
95protected:
96
98
99private:
100
101 G4bool fasterCode;
102
103 //deexcitation manager to produce fluo photns and e-
104 G4VAtomDeexcitation* fAtomDeexcitation;
105
106 G4Material* nistSi;
107
108 std::map<G4String,G4double,std::less<G4String> > lowEnergyLimit;
109 std::map<G4String,G4double,std::less<G4String> > highEnergyLimit;
110
111 G4bool isInitialised;
112 G4int verboseLevel;
113
114 // Cross section
115
116 typedef std::map<G4String,G4String,std::less<G4String> > MapFile;
117 MapFile tableFile;
118
119 typedef std::map<G4String,G4MicroElecCrossSectionDataSet*,std::less<G4String> > MapData;
120 MapData tableData;
121
122 // Final state
123
124 G4MicroElecSiStructure SiStructure;
125
126 G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
127
128 G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
129
130 G4double Interpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
131
132 G4double QuadInterpolator( G4double e11,
133 G4double e12,
134 G4double e21,
135 G4double e22,
136 G4double x11,
137 G4double x12,
138 G4double x21,
139 G4double x22,
140 G4double t1,
141 G4double t2,
142 G4double t,
143 G4double e);
144
145 typedef std::map<double, std::map<double, double> > TriDimensionMap;
146
147 TriDimensionMap eDiffCrossSectionData[7];
148 TriDimensionMap eNrjTransfData[7]; // for cumulated dcs
149
150 TriDimensionMap pDiffCrossSectionData[7];
151 TriDimensionMap pNrjTransfData[7]; // for cumulated dcs
152
153 std::vector<double> eTdummyVec;
154 std::vector<double> pTdummyVec;
155
156 typedef std::map<double, std::vector<double> > VecMap;
157
158 VecMap eVecm;
159 VecMap pVecm;
160
161 VecMap eProbaShellMap[7]; // for cumulated dcs
162 VecMap pProbaShellMap[7]; // for cumulated dcs
163
164 // Partial cross section
165
166 G4int RandomSelect(G4double energy,const G4String& particle );
167
168 //
169
172
173};
174
175//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
176
178{
179 fasterCode = input;
180}
181
182//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
183
184#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4double TransferedEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell, G4double random)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4ParticleChangeForGamma * fParticleChangeForGamma
double DifferentialCrossSection(G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)