Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MicroElecInelasticModel_new.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_new.hh, 2011/08/29 A.Valentin, M. Raine are with CEA [a]
28// 2020/05/20 P. Caron, C. Inguimbert are with ONERA [b]
29// Q. Gibaru is with CEA [a], ONERA [b] and CNES [c]
30// M. Raine and D. Lambert are with CEA [a]
31//
32// A part of this work has been funded by the French space agency(CNES[c])
33// [a] CEA, DAM, DIF - 91297 ARPAJON, France
34// [b] ONERA - DPHY, 2 avenue E.Belin, 31055 Toulouse, France
35// [c] CNES, 18 av.E.Belin, 31401 Toulouse CEDEX, France
36//
37// Based on the following publications
38// - A.Valentin, M. Raine,
39// Inelastic cross-sections of low energy electrons in silicon
40// for the simulation of heavy ion tracks with the Geant4-DNA toolkit,
41// NSS Conf. Record 2010, pp. 80-85
42// https://doi.org/10.1109/NSSMIC.2010.5873720
43//
44// - A.Valentin, M. Raine, M.Gaillardin, P.Paillet
45// Geant4 physics processes for microdosimetry simulation:
46// very low energy electromagnetic models for electrons in Silicon,
47// https://doi.org/10.1016/j.nimb.2012.06.007
48// NIM B, vol. 288, pp. 66-73, 2012, part A
49// heavy ions in Si, NIM B, vol. 287, pp. 124-129, 2012, part B
50// https://doi.org/10.1016/j.nimb.2012.07.028
51//
52// - M. Raine, M. Gaillardin, P. Paillet
53// Geant4 physics processes for silicon microdosimetry simulation:
54// Improvements and extension of the energy-range validity up to 10 GeV/nucleon
55// NIM B, vol. 325, pp. 97-100, 2014
56// https://doi.org/10.1016/j.nimb.2014.01.014
57//
58// - J. Pierron, C. Inguimbert, M. Belhaj, T. Gineste, J. Puech, M. Raine
59// Electron emission yield for low energy electrons:
60// Monte Carlo simulation and experimental comparison for Al, Ag, and Si
61// Journal of Applied Physics 121 (2017) 215107.
62// https://doi.org/10.1063/1.4984761
63//
64// - P. Caron,
65// Study of Electron-Induced Single-Event Upset in Integrated Memory Devices
66// PHD, 16th October 2019
67//
68// - Q.Gibaru, C.Inguimbert, P.Caron, M.Raine, D.Lambert, J.Puech,
69// Geant4 physics processes for microdosimetry and secondary electron emission simulation :
70// Extension of MicroElec to very low energies and new materials
71// NIM B, 2020, in review.
72//
73//
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
75
76#ifndef G4MICROELECINELASTICMODEL_NEW_HH
77#define G4MICROELECINELASTICMODEL_NEW_HH 1
78
79#include "globals.hh"
80#include "G4VEmModel.hh"
85#include "G4Electron.hh"
86#include "G4Proton.hh"
87#include "G4GenericIon.hh"
91#include "G4NistManager.hh"
92
93
95{
96
97public:
98
100 const G4String& nam = "MicroElecInelasticModel");
102
103 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
104
106 const G4ParticleDefinition* p,
107 G4double ekin,
108 G4double emin,
109 G4double emax) override;
110
111 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
113 const G4DynamicParticle*,
114 G4double tmin,
115 G4double maxEnergy) override;
116
117 G4double DifferentialCrossSection(const G4ParticleDefinition * aParticleDefinition,
118 G4double k, G4double energyTransfer, G4int shell);
119
121
123
124 G4double BKZ(G4double Ep, G4double mp, G4int Zp, G4double EF);
125 // compute the effective charge according Brandt et Kitagawa theory
126
128
130
131protected:
133
134private:
135 //
136 // private methods
137 //
140
141 G4int RandomSelect(G4double energy,const G4String& particle, G4double originalMass, G4int originalZ );
142
143 G4double RandomizeCreatedElectronEnergy(G4double secondaryKinetic);
144
145 G4double RandomizeEjectedElectronEnergy(const G4ParticleDefinition * aParticleDefinition,
146 G4double incomingParticleEnergy, G4int shell,
147 G4double originalMass, G4int originalZ) ;
148
149 G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(const G4ParticleDefinition*,
150 G4double k, G4int shell);
151
152 G4double TransferedEnergy(const G4ParticleDefinition*, G4double k,
153 G4int ionizationLevelIndex, G4double random);
154
155 G4double Interpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
156
157 G4double QuadInterpolator( G4double e11, G4double e12, G4double e21, G4double e22,
158 G4double x11, G4double x12, G4double x21, G4double x22,
159 G4double t1, G4double t2, G4double t, G4double e);
160 //
161 // private elements
162 //
163 G4String currentMaterial = "";
164 G4bool fasterCode = false;
165 //deexcitation manager to produce fluo photns and e-
166 G4VAtomDeexcitation* fAtomDeexcitation = nullptr;
167 G4Material* nistSi = nullptr;
168 std::map<G4String,G4double,std::less<G4String> > lowEnergyLimit;
169 std::map<G4String,G4double,std::less<G4String> > highEnergyLimit;
170 G4bool isInitialised = false;
171 G4int verboseLevel = 0;
172 typedef std::map<G4String,G4String,std::less<G4String> > MapFile;
173 typedef std::map<G4String,G4MicroElecCrossSectionDataSet_new*,std::less<G4String> > MapData;
174 typedef std::map<G4double, std::map<G4double, G4double> > TriDimensionMap;
175 typedef std::map<G4double, std::vector<G4double> > VecMap;
176 //Tables for multilayers
177 typedef std::map<G4String, MapData*, std::less<G4String> > TCSMap;
178 TCSMap tableTCS; //TCS tables by particle
179 typedef std::map<G4String, std::vector<TriDimensionMap>* > dataDiffCSMap;
180 dataDiffCSMap eDiffDatatable, pDiffDatatable; //Transfer probabilities (for slower code)
181 dataDiffCSMap eNrjTransStorage, pNrjTransStorage; //Transfered energies and corresponding probability (faster code)
182 typedef std::map<G4String, std::vector<VecMap>* > dataProbaShellMap;
183 dataProbaShellMap eProbaShellStorage, pProbaShellStorage; //Cumulated Transfer probabilities (faster code)
184 typedef std::map<G4String, std::vector<G4double>* > incidentEnergyMap;
185 incidentEnergyMap eIncidentEnergyStorage, pIncidentEnergyStorage; //Incident energies for interpolation (faster code)
186 typedef std::map<G4String, VecMap* > TranfEnergyMap;
187 TranfEnergyMap eVecmStorage, pVecmStorage; //Transfered energy for interpolation (slower code)
188 typedef std::map<G4String, G4MicroElecMaterialStructure*, std::less<G4String> > MapStructure;
189 MapStructure tableMaterialsStructures; //Structures of all materials simulated
190 G4MicroElecMaterialStructure* currentMaterialStructure = nullptr;
191
192};
193
194#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4double DifferentialCrossSection(const G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
G4double ComputeElasticQmax(G4double T1i, G4double T2i, G4double m1, G4double m2)
G4ParticleChangeForGamma * fParticleChangeForGamma
G4double BKZ(G4double Ep, G4double mp, G4int Zp, G4double EF)
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double ComputeRelativistVelocity(G4double E, G4double mass)
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double vrkreussler(G4double v, G4double vF)