Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MicroElecElasticModel.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.hh, 2011/08/29 A.Valentin, M. Raine
28//
29// Based on the following publications
30// - Geant4 physics processes for microdosimetry simulation:
31// very low energy electromagnetic models for electrons in Si,
32// NIM B, vol. 288, pp. 66 - 73, 2012.
33//
34//
35//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
36
37#ifndef G4MicroElecElasticModel_h
38#define G4MicroElecElasticModel_h 1
39
40#include <map>
42
44#include "G4VEmModel.hh"
45#include "G4Electron.hh"
49#include "G4NistManager.hh"
50
52{
53
54public:
56 const G4String& nam = "MicroElecElasticModel");
58
59 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
60
62 const G4ParticleDefinition* p,
63 G4double ekin,
64 G4double emin,
65 G4double emax) override;
66
67 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
69 const G4DynamicParticle*,
70 G4double tmin,
71 G4double maxEnergy) override;
72
73 inline void SetKillBelowThreshold (G4double threshold);
74 G4double GetKillBelowThreshold () { return killBelowEnergy; }
75
78
79protected:
81
82private:
83 // Final state
84 G4double Theta(G4ParticleDefinition * aParticleDefinition, G4double k, G4double integrDiff);
85 G4double LinLinInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
86 G4double LogLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
87 G4double LinLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
88 G4double QuadInterpolator(G4double e11,
89 G4double e12,
90 G4double e21,
91 G4double e22,
92 G4double x11,
93 G4double x12,
94 G4double x21,
95 G4double x22,
96 G4double t1,
97 G4double t2,
98 G4double t,
99 G4double e);
100 G4double RandomizeCosTheta(G4double k);
101
102 // Cross section
103 typedef std::map<G4String,G4String,std::less<G4String> > MapFile;
104 MapFile tableFile;
105
106 typedef std::map<G4String,G4MicroElecCrossSectionDataSet*,std::less<G4String> > MapData;
107 MapData tableData;
108
109 typedef std::map<G4double, std::map<G4double, G4double> > TriDimensionMap;
110 TriDimensionMap eDiffCrossSectionData;
111 std::vector<G4double> eTdummyVec;
112
113 typedef std::map<G4double, std::vector<G4double> > VecMap;
114 VecMap eVecm;
115
116 G4Material* nistSi;
117 G4double killBelowEnergy;
118 G4double lowEnergyLimit;
119 G4double lowEnergyLimitOfModel;
120 G4double highEnergyLimit;
121 G4int verboseLevel;
122 G4bool isInitialised;
123};
124
126{
127 killBelowEnergy = threshold;
128
129 if (threshold < 5*CLHEP::eV)
130 {
131 G4Exception ("*** WARNING : the G4MicroElecElasticModel class is not validated below 5 eV !","",JustWarning,"") ;
132 threshold = 5*CLHEP::eV;
133 }
134
135}
136
137
138//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
139
140#endif
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4MicroElecElasticModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="MicroElecElasticModel")
void SetKillBelowThreshold(G4double threshold)
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
G4ParticleChangeForGamma * fParticleChangeForGamma
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4MicroElecElasticModel(const G4MicroElecElasticModel &)=delete
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4MicroElecElasticModel & operator=(const G4MicroElecElasticModel &right)=delete