Geant4-11
G4MicroElecSurface.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// G4MicroElecSurface.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// - Q.Gibaru, C.Inguimbert, P.Caron, M.Raine, D.Lambert, J.Puech,
40// Geant4 physics processes for microdosimetry and secondary electron emission simulation :
41// Extension of MicroElec to very low energies and new materials
42// NIM B, 2020, in review.
43//
44//
45// - Modèle de transport d'électrons à basse énergie (10 eV- 2 keV) pour
46// applications spatiales (OSMOSEE, GEANT4), PhD dissertation, 2017.
47//
48//
50
51#include "G4MicroElecSurface.hh"
52
53#include "G4ios.hh"
55#include "G4EmProcessSubType.hh"
57#include "G4EmProcessSubType.hh"
58
59#include "G4SystemOfUnits.hh"
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62
64 : G4VDiscreteProcess(processName, type),
65 oldMomentum(0.,0.,0.), previousMomentum(0.,0.,0.),
66 theGlobalNormal(0.,0.,0.), theFacetNormal(0.,0.,0.)
67{
68 if ( verboseLevel > 0)
69 {
70 G4cout << GetProcessName() << " is created " << G4endl;
71 }
72
73 isInitialised=false;
75
77 material1 = nullptr;
78 material2 = nullptr;
79
82
84 flag_normal = false;
85 flag_reflexion = false;
86 teleportToDo = teleportDone = false;
87
89}
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
92
94{}
95
96//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
97
98G4bool
100{
101 return ( aParticleType.GetPDGEncoding() == 11 );
102}
103
104//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
105
107{
108 if (isInitialised) { return; }
109
110 G4ProductionCutsTable* theCoupleTable =
112 G4int numOfCouples = theCoupleTable->GetTableSize();
113 G4cout << "G4MicroElecSurface::Initialise: Ncouples= "
114 << numOfCouples << G4endl;
115
116 for (G4int i = 0; i < numOfCouples; ++i) {
117 const G4Material* material =
118 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
119
120 G4cout << "G4Surface, Material " << i + 1 << " / " << numOfCouples << " : " << material->GetName() << G4endl;
121 if (material->GetName() == "Vacuum") { tableWF[material->GetName()] = 0; continue; }
122 G4String mat = material->GetName();
124 tableWF[mat] = str.GetWorkFunction();
125 }
126 isInitialised = true;
127}
128
129//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
130
132{
134
135 //Definition of the parameters for the particle
138
139 G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
140 G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
141
142 material1 = pPreStepPoint -> GetMaterial();
143 material2 = pPostStepPoint -> GetMaterial();
144
145 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
146
149 oldMomentum = aParticle->GetMomentumDirection();
150
151 //First case: not a boundary
152 if (pPostStepPoint->GetStepStatus() != fGeomBoundary ||
153 pPostStepPoint->GetPhysicalVolume() == pPreStepPoint->GetPhysicalVolume())
154 {
157 flag_reflexion = false;
158 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
159
160 }
162
163 //Third case: same material
164 if (material1 == material2)
165 {
167 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
168
169 }
170 if (verboseLevel > 0)
171 {
172 G4cout << G4endl << " Electron at Boundary! " << G4endl;
173 G4VPhysicalVolume* thePrePV = pPreStepPoint->GetPhysicalVolume();
174 G4VPhysicalVolume* thePostPV = pPostStepPoint->GetPhysicalVolume();
175 if (thePrePV) G4cout << " thePrePV: " << thePrePV->GetName() << G4endl;
176 if (thePostPV) G4cout << " thePostPV: " << thePostPV->GetName() << G4endl;
177 G4cout << " Old Momentum Direction: " << oldMomentum << G4endl;
178 }
179
180 //Definition of the parameters for the surface
181 G4ThreeVector theGlobalPoint = pPostStepPoint->GetPosition();
182
183 G4Navigator* theNavigator =
185 GetNavigatorForTracking();
186
187 G4bool valid;
188
189 theGlobalNormal = theNavigator->GetGlobalExitNormal(theGlobalPoint, &valid);
190 // G4cout << "Global exit normal = " << theGlobalNormal << " valid = " << valid << G4endl;
191
192 if (valid)
193 {
195 }
196 else
197 {
199 ed << " G4MicroElecSurface/PostStepDoIt(): "
200 << " The Navigator reports that it returned an invalid normal.\n"
201 << "PV: " << pPreStepPoint->GetPhysicalVolume()->GetName()
202 << " TrackID= " << aTrack.GetTrackID()
203 << " Ekin(MeV)= " << aTrack.GetKineticEnergy()
204 << " position: " << theGlobalPoint
205 << " direction: " << oldMomentum
206 << G4endl;
207 G4Exception("G4MuElecSurf::PostStepDoIt", "OpBoun01",
208 FatalException, ed,
209 "Invalid Surface Normal - Geometry must return valid surface normal");
210 return 0;
211 }
212
213 //Exception: the particle is not in the right direction
214 if (oldMomentum * theGlobalNormal > 0.0)
215 {
217 }
218
219 //Second case: step too small
220 //Corrections bug rotation + réflexion
221 if (aTrack.GetStepLength()<=kCarTolerance)
222 {
224
225 WorkFunctionTable::iterator postStepWF;
226 postStepWF = tableWF.find(pPostStepPoint->GetMaterial()->GetName());
227 WorkFunctionTable::iterator preStepWF;
228 preStepWF = tableWF.find(pPreStepPoint->GetMaterial()->GetName());
229
230 if (postStepWF == tableWF.end()) {
231 G4String str = "Material ";
232 str += pPostStepPoint->GetMaterial()->GetName() + " not found!";
233 G4Exception("G4Surface::G4Surface", "em0002", FatalException, str);
234 return 0;
235 }
236
237 else if (preStepWF == tableWF.end()) {
238 G4String str = "Material ";
239 str += pPreStepPoint->GetMaterial()->GetName() + " not found!";
240 G4Exception("G4Surface::G4Surface", "em0002", FatalException, str);
241 return 0;
242 }
243
244 if (pPreStepPoint->GetMaterial() != pPostStepPoint->GetMaterial()) {
245
247
248 if (flag_reflexion == true && flag_normal == true) {
250 flag_reflexion = false;
251 flag_normal = false;
252 }
253 }
254 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
255 }
256
257 flag_normal = (theGlobalNormal.x() == 0.0 && theGlobalNormal.y() == 0.0);
258
259 G4LogicalSurface* Surface = nullptr;
260
262 (pPreStepPoint ->GetPhysicalVolume(),
263 pPostStepPoint->GetPhysicalVolume());
264
265 if (Surface == nullptr)
266 {
267 G4bool enteredDaughter=(pPostStepPoint->GetPhysicalVolume()
268 ->GetMotherLogical() ==
269 pPreStepPoint->GetPhysicalVolume()
270 ->GetLogicalVolume());
271 if(enteredDaughter)
272 {
274 (pPostStepPoint->GetPhysicalVolume()->
275 GetLogicalVolume());
276
277 if(Surface == nullptr)
279 (pPreStepPoint->GetPhysicalVolume()->
280 GetLogicalVolume());
281 }
282 else
283 {
285 (pPreStepPoint->GetPhysicalVolume()->
286 GetLogicalVolume());
287
288 if(Surface == nullptr)
290 (pPostStepPoint->GetPhysicalVolume()->
291 GetLogicalVolume());
292 }
293 }
294
295 G4VPhysicalVolume* thePrePV = pPreStepPoint->GetPhysicalVolume();
296 G4VPhysicalVolume* thePostPV = pPostStepPoint->GetPhysicalVolume();
297
298 if (thePostPV)
299 {
300 WorkFunctionTable::iterator postStepWF;
301 postStepWF = tableWF.find(thePostPV->GetLogicalVolume()->GetMaterial()->GetName());
302 WorkFunctionTable::iterator preStepWF;
303 preStepWF = tableWF.find(thePrePV->GetLogicalVolume()->GetMaterial()->GetName());
304
305 if (postStepWF == tableWF.end()) {
306 G4String str = "Material ";
307 str += thePostPV->GetLogicalVolume()->GetMaterial()->GetName() + " not found!";
308 G4Exception("G4Surface::G4Surface", "em0002", FatalException, str);
309 return 0;
310 }
311 else if (preStepWF == tableWF.end()) {
312 G4String str = "Material ";
313 str += thePrePV->GetLogicalVolume()->GetMaterial()->GetName() + " not found!";
314 G4Exception("G4Surface::G4Surface", "em0002", FatalException, str);
315 return 0;
316 }
317 else
318 {
319 G4double thresholdNew = postStepWF->second;
320 G4double thresholdOld = preStepWF->second;
321 energyThreshold = thresholdNew - thresholdOld;
322 }
323 }
324
325 ekint = pPreStepPoint->GetKineticEnergy();
326 thetat= GetIncidentAngle(); //angle d'incidence
327 G4double ekinNormalt=ekint*std::cos(thetat)*std::cos(thetat);
328
329 G4double atet = std::sqrt(ekint/(ekint+energyThreshold))*std::sin(thetat);
330 thetaft = (atet > 1.0) ? pi*0.5 : std::asin(atet);//Angle de réfraction
331
332 G4double aleat=G4UniformRand();
333
334 const G4double waveVectort=std::sqrt(2*9.1093826E-31*1.602176487E-19)/(6.6260755E-34/(2.0*pi));
335
336 //Parameter for an exponential barrier of potential (Thèse P68)
337 const G4double at=0.5E-10;
338
339 //G4double modif already declared in .hh
341
342 G4double kft=waveVectort*std::sqrt(ekint+energyThreshold)*std::cos(thetaft);
343 G4double kit=waveVectort*std::sqrt(ekinNormalt);
344
345 G4double yy = std::sinh(pi*at*(kit-kft))/std::sinh(pi*at*(kit+kft));
346 crossingProbability = 1 - yy*yy;
347
348 //First case: the electron crosses the surface
349 if((aleat<=crossingProbability)&&(ekint>std::abs(energyThreshold)))
350 {
351 if (pPreStepPoint->GetMaterial() != pPostStepPoint->GetMaterial()) {
353 }
354
355 thetaft=std::abs(thetaft-thetat);
356
358 G4ThreeVector xVerst = zVerst.orthogonal();
359 G4ThreeVector yVerst = zVerst.cross(xVerst);
360
361 G4double cost = std::cos(thetaft);
362 G4double xDirt = std::sqrt(1. - cost*cost);
363 G4double yDirt = xDirt;
364
365 G4ThreeVector zPrimeVerst = xDirt*xVerst + yDirt*yVerst + cost*zVerst;
366
368 }
369 else if ((aleat > crossingProbability) && (ekint>std::abs(energyThreshold)))
370 {
371 flag_reflexion = true;
374
375 }
376 else {
377
380 flag_reflexion = true;
381 }
382
383 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
384}
385
386//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
387
390{
391 *condition = Forced;
392 return DBL_MAX;
393}
394
395//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
396
398{
400
402 G4double magP= oldMomentum.mag();
404
405 G4double incidentangle = pi - std::acos(PdotN/(magP*magN));
406
407 return incidentangle;
408}
409
410//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
411
413{
414 //Normale
418
419 //PostStepPoint
420 G4double PSx = PostStepPoint->GetPosition().x();
421 G4double PSy = PostStepPoint->GetPosition().y();
422 G4double PSz = PostStepPoint->GetPosition().z();
423
424 //P(alpha,beta,gamma) - PostStep avec translation momentum
425 G4double alpha = PSx + oldMomentum.x();
426 G4double beta = PSy + oldMomentum.y();
427 G4double gamma = PSz + oldMomentum.z();
429 G4double x, y, z, d, A, B, PM2x, PM2y, PM2z;
430 d = -(Nx*PSx + Ny*PSy + Nz*PSz);
431
432 if (Ny == 0 && Nx == 0) {
433 gamma = -gamma;
434 }
435 else {
436 if (Ny == 0) {
437
438 A = (Nz*Nz*alpha) + (Nx*Nx*PSx) + (Nx*Nz*(PSz - gamma));
439 B = r*r;
440
441 //M(x,y,z) - Projection de P sur la surface
442 x = A / B;
443 y = beta;
444 z = (x - alpha)*(Nz / Nx) + gamma;
445 }
446
447 else {
448
449 A = (r*r) / Ny;
450 B = (beta / Ny)*(Nx*Nx + Nz*Nz) - (Nx*alpha + Nz*gamma + d);
451
452 //M(x,y,z) - Projection de P sur la surface
453 y = B / A;
454 x = (y - beta)*(Nx / Ny) + alpha;
455 z = (y - beta)*(Nz / Ny) + gamma;
456 }
457
458 //Vecteur 2*PM
459 PM2x = 2 * (x - alpha); PM2y = 2 * (y - beta); PM2z = 2 * (z - gamma);
460
461 //Nouveau point P
462 alpha += PM2x; beta += PM2y; gamma += PM2z;
463
464 }
465
466 G4ThreeVector newMomentum = G4ThreeVector(alpha-PSx,beta-PSy,gamma-PSz);
467 return newMomentum.unit();
468}
469
470//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
471
473{
474 return theStatus;
475}
476
477//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
478
479
481{
483}
484
485//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
486
G4double B(G4double temperature)
@ fSurfaceReflection
G4double condition(const G4ErrorSymMatrix &m)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ForceCondition
@ Forced
static const G4double alpha
G4MicroElecSurfaceStatus
@ StepTooSmallSurf
@ SameMaterialSurf
@ UndefinedSurf
@ NotAtBoundarySurf
G4ProcessType
static constexpr double pi
Definition: G4SIunits.hh:55
@ fGeomBoundary
Definition: G4StepStatus.hh:43
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
double z() const
Hep3Vector unit() const
Hep3Vector orthogonal() const
double x() const
double y() const
Hep3Vector cross(const Hep3Vector &) const
double mag() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetTotalMomentum() const
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
static G4LogicalBorderSurface * GetSurface(const G4VPhysicalVolume *vol1, const G4VPhysicalVolume *vol2)
static G4LogicalSkinSurface * GetSurface(const G4LogicalVolume *vol)
G4Material * GetMaterial() const
const G4Material * GetMaterial() const
const G4String & GetName() const
Definition: G4Material.hh:173
G4ThreeVector theFacetNormal
WorkFunctionTable tableWF
G4MicroElecSurfaceStatus GetStatus() const
G4MicroElecSurfaceStatus theStatus
G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *condition) override
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
~G4MicroElecSurface() override
G4ThreeVector theGlobalNormal
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4ThreeVector Reflexion(const G4StepPoint *PostStepPoint)
G4MicroElecSurface(const G4String &processName="MicroElecSurface", G4ProcessType type=fElectromagnetic)
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
const G4Material * material2
const G4Material * material1
G4ThreeVector previousMomentum
virtual G4ThreeVector GetGlobalExitNormal(const G4ThreeVector &point, G4bool *valid)
void ProposeVelocity(G4double finalVelocity)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
virtual void Initialize(const G4Track &)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4StepStatus GetStepStatus() const
G4Material * GetMaterial() const
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
G4VPhysicalVolume * GetPhysicalVolume() const
G4double GetKineticEnergy() const
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
G4double GetVelocity() const
G4int GetTrackID() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
G4double GetStepLength() const
static G4TransportationManager * GetTransportationManager()
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4LogicalVolume * GetMotherLogical() const
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:327
G4int verboseLevel
Definition: G4VProcess.hh:356
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
string material
Definition: eplot.py:19
#define DBL_MAX
Definition: templates.hh:62