Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MicroElecLOPhononModel.cc
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// G4MicroElecLOPhononModel.cc,
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//
39// - J. Pierron, C. Inguimbert, M. Belhaj, T. Gineste, J. Puech, M. Raine
40// Electron emission yield for low energy electrons:
41// Monte Carlo simulation and experimental comparison for Al, Ag, and Si
42// Journal of Applied Physics 121 (2017) 215107.
43// https://doi.org/10.1063/1.4984761
44//
45// - P. Caron,
46// Study of Electron-Induced Single-Event Upset in Integrated Memory Devices
47// PHD, 16th October 2019
48//
49// - Q.Gibaru, C.Inguimbert, P.Caron, M.Raine, D.Lambert, J.Puech,
50// Geant4 physics processes for microdosimetry and secondary electron emission simulation :
51// Extension of MicroElec to very low energies and new materials
52// NIM B, 2020, in review.
53//
54////////////////////////////////////////////////////////////////////////
55
57#include "G4SystemOfUnits.hh"
59
61 const G4String& nam)
62 : G4VEmModel(nam),isInitialised(false)
63{
64 G4cout << "Phonon model is constructed " << G4endl
65 << "Phonon Energy = " << phononEnergy / eV << " eV "<< G4endl;
66}
67
68//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
69
71 const G4DataVector& /*cuts*/)
72{
73 fParticleChangeForGamma = GetParticleChangeForGamma();
74 isInitialised = true;
75}
76
77//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
78
80CrossSectionPerVolume(const G4Material* material,
81 const G4ParticleDefinition* p,
82 G4double ekin,
84{
85 const G4double e = CLHEP::eplus / CLHEP::coulomb;
86 const G4double m0 = CLHEP::electron_mass_c2 / (CLHEP::c_squared*CLHEP::kg);
87 const G4double h = CLHEP::hbar_Planck * CLHEP::s/ (CLHEP::m2*CLHEP::kg);
88 const G4double eps0 = CLHEP::epsilon0 * CLHEP::m/ (CLHEP::farad);
89 const G4double kb = CLHEP::k_Boltzmann * CLHEP::kelvin/ CLHEP::joule;
90 const G4double T = 300;
91 G4double eps = 9;
92 G4double einf = 3;
93
94 const G4DataVector cuts;
95 Initialise(p, cuts);
96
97 if (material->GetName() != "G4_SILICON_DIOXIDE"
98 && material->GetName() != "G4_ALUMINUM_OXIDE"
99 && material->GetName() != "G4_BORON_NITRIDE")
100 {
101 return 1 / DBL_MAX;
102 }
103
104 G4double E =(ekin/eV)*e;
105
106 if (material->GetName() == "G4_ALUMINUM_OXIDE")
107 {
108 eps = 9;
109 einf = 3;
110 phononEnergy = 0.1*eV;
111 }
112 if (material->GetName() == "G4_SILICON_DIOXIDE")
113 {
114 eps = 3.84;
115 einf = 2.25;
116 phononEnergy = (0.75*0.153+0.25*0.063 )* eV;
117 }
118
119 // Nuclear Instruments and Methods in Physics Research Section B:
120 // Beam Interactions with Materials and Atoms
121 // Volume 454, 1 September 2019, Pages 14 - 22
122 // Nuclear Instruments and Methods in Physics Research Section B:
123 // Beam Interactions with Materials and Atoms
124 // Monte Carlo modeling of low - energy electron - induced secondary
125 // electron emission yields in micro - architected boron nitride surfaces
126
127 if (material->GetName() == "G4_BORON_NITRIDE")
128 {
129 eps = 7.1;
130 einf = 4.5;
131 phononEnergy = 0.17 * eV;
132 }
133
134 G4double hw = (phononEnergy / eV) * e;
135 G4double n = 1.0 / (std::exp(hw / (kb*T)) - 1); //Phonon distribution
136
137 if (E<=hw)
138 {
139 return 1 / DBL_MAX;
140 }
141
142 if (absor) // Absorption
143 {
144 Eprim = E + hw;
145 signe = -1;
146 }
147 else // Emission
148 {
149 Eprim = E - hw;
150 signe = +1;
151 }
152
153 G4double racine = std::sqrt(1 + ((-signe*hw) / E));
154 G4double P = (std::pow(e, 2) / (4 * pi*eps0*h*h)) * (n + 0.5 + signe*0.5) * ((1 / einf) - (1 / eps)) * std::sqrt(m0 / (2 * E)) *hw* std::log((1 + racine) / (signe * 1 + ((-signe)*racine)));
155 G4double MFP = (std::sqrt(2 * E / m0) / P)*m;
156
157 if (material->GetName() == "G4_SILICON_DIOXIDE") { return 2 / MFP; }
158 return 1/(MFP);
159}
160
161//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
162
164SampleSecondaries(std::vector<G4DynamicParticle*>*,
166 const G4DynamicParticle* aDynamicElectron,
168{
169 G4double E = aDynamicElectron->GetKineticEnergy();
170 Eprim = (absor) ? E + phononEnergy : E - phononEnergy;
171
172 G4double rand = G4UniformRand();
173 G4double B = (E + Eprim + 2 * std::sqrt(E*Eprim))
174 / (E + Eprim - 2 * std::sqrt(E*Eprim));
175 G4double cosTheta = ((E + Eprim) / (2 * std::sqrt(E*Eprim)))
176 * (1 - std::pow(B, rand)) + std::pow(B, rand);
177 if(Interband)
178 {
179 cosTheta = 1 - 2 * G4UniformRand(); //Isotrope
180 }
181 G4double phi = twopi * G4UniformRand();
182 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
183 G4ThreeVector xVers = zVers.orthogonal();
184 G4ThreeVector yVers = zVers.cross(xVers);
185
186 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
187 G4double yDir = xDir;
188 xDir *= std::cos(phi);
189 yDir *= std::sin(phi);
190
191 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
192
193 fParticleChangeForGamma->ProposeMomentumDirection(zPrimeVers.unit());
194 fParticleChangeForGamma->SetProposedKineticEnergy(Eprim);
195}
196
197//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4double B(G4double temperature)
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector unit() const
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
const G4String & GetName() const
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin=0.0, G4double emax=DBL_MAX) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double maxEnergy=DBL_MAX) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4MicroElecLOPhononModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="G4MicroElecLOPhononModel")
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
#define DBL_MAX
Definition templates.hh:62