Geant4-11
G4MicroElecElasticModel_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// G4MicroElecElasticModel_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 G4MICROELECELASTICMODEL_NEW_HH
77#define G4MICROELECELASTICMODEL_NEW_HH 1
78
79#include <map>
81
84#include "G4VEmModel.hh"
85#include "G4Electron.hh"
89#include "G4NistManager.hh"
90
92{
93
94public:
96 const G4String& nam = "MicroElecElasticModel");
98
99 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
100
102 const G4ParticleDefinition* p,
103 G4double ekin,
104 G4double emin,
105 G4double emax) override;
106
108 G4double cs, G4double Aac, G4double Eac,
109 G4double prefactor);
110
111 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
113 const G4DynamicParticle*,
114 G4double tmin,
115 G4double maxEnergy) override;
116
117 void SetKillBelowThreshold (G4double threshold);
118
120
122
123protected:
125
126private:
127
130
131 // Final state
132 G4double Theta(G4ParticleDefinition * aParticleDefinition, G4double k, G4double integrDiff);
137 G4double x11, G4double x12, G4double x21, G4double x22,
138 G4double t1, G4double t2, G4double t, G4double e);
139
141
142 G4Material* nistSi = nullptr;
149 // Cross section
150 typedef std::map<G4String,G4String,std::less<G4String> > MapFile;
152 typedef std::map<G4String,G4MicroElecCrossSectionDataSet_new*,std::less<G4String> > MapData;
153 //MapData tableData;
154
155 typedef std::map<G4String, MapData*, std::less<G4String> > TCSMap;
157
158 //Maps for multilayers
159 typedef std::map<G4double, std::map<G4double, G4double> > TriDimensionMap;
160
161 typedef std::map<G4String, TriDimensionMap* > ThetaMap;
162 ThetaMap thetaDataStorage; //Storage of angles (cumulated)
163
164 typedef std::map<G4String, std::vector<G4double>* > energyMap;
166
167 typedef std::map<G4double, std::vector<G4double> > VecMap;
168
169 typedef std::map<G4String, VecMap* > ProbaMap;
170 ProbaMap eProbaStorage; //Storage of probabilities for cumulated sections
171
172 typedef std::map<G4String, G4MicroElecMaterialStructure*, std::less<G4String> > MapStructure;
173
174 MapStructure tableMaterialsStructures; //Structures of all materials simulated
175
177 typedef std::map<G4String, G4double, std::less<G4String> > MapEnergy;
181
185
186};
187
188//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
189
190#endif
static const G4double e1[44]
static const G4double e2[44]
static const G4double emax
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
std::map< G4double, std::map< G4double, G4double > > TriDimensionMap
std::map< G4String, G4MicroElecMaterialStructure *, std::less< G4String > > MapStructure
std::map< G4String, G4double, std::less< G4String > > MapEnergy
G4double AcousticCrossSectionPerVolume(G4double ekin, G4double kbz, G4double rho, G4double cs, G4double Aac, G4double Eac, G4double prefactor)
void SetKillBelowThreshold(G4double threshold)
std::map< G4String, VecMap * > ProbaMap
G4double LinLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
G4MicroElecElasticModel_new(const G4ParticleDefinition *p=0, const G4String &nam="MicroElecElasticModel")
std::map< G4String, G4String, std::less< G4String > > MapFile
G4MicroElecElasticModel_new & operator=(const G4MicroElecElasticModel_new &right)
std::map< G4String, MapData *, std::less< G4String > > TCSMap
std::map< G4double, std::vector< G4double > > VecMap
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4MicroElecElasticModel_new(const G4MicroElecElasticModel_new &)
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double DamageEnergy(G4double T, G4double A, G4double Z)
std::map< G4String, G4MicroElecCrossSectionDataSet_new *, std::less< G4String > > MapData
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
std::map< G4String, TriDimensionMap * > ThetaMap
G4double LinLinInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
G4double QuadInterpolator(G4double e11, G4double e12, G4double e21, G4double e22, G4double x11, G4double x12, G4double x21, G4double x22, G4double t1, G4double t2, G4double t, G4double e)
std::map< G4String, std::vector< G4double > * > energyMap
G4MicroElecMaterialStructure * currentMaterialStructure
G4double Theta(G4ParticleDefinition *aParticleDefinition, G4double k, G4double integrDiff)
G4ParticleChangeForGamma * fParticleChangeForGamma
G4double LogLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
string material
Definition: eplot.py:19