Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MicroElecInelasticModel.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.hh, 2011/08/29 A.Valentin, M. Raine
28//
29// Based on the following publications
30//
31// - Inelastic cross-sections of low energy electrons in silicon
32// for the simulation of heavy ion tracks with theGeant4-DNA toolkit,
33// NSS Conf. Record 2010, pp. 80-85
34// - Geant4 physics processes for microdosimetry simulation:
35// very low energy electromagnetic models for electrons in Si,
36// NIM B, vol. 288, pp. 66-73, 2012.
37// - Geant4 physics processes for microdosimetry simulation:
38// very low energy electromagnetic models for protons and
39// heavy ions in Si, NIM B, vol. 287, pp. 124-129, 2012.
40//
41//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
42
43#ifndef G4MicroElecInelasticModel_h
44#define G4MicroElecInelasticModel_h 1
45
46
47#include "globals.hh"
48#include "G4VEmModel.hh"
53
55class G4NistManager;
57
59{
60
61public:
62
64 const G4String& nam = "MicroElecInelasticModel");
66
67 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
68
70 const G4ParticleDefinition* p,
71 G4double ekin,
72 G4double emin,
73 G4double emax) override;
74 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
76 const G4DynamicParticle*,
77 G4double tmin,
78 G4double maxEnergy) override;
80 G4double k, G4double energyTransfer,
81 G4int shell);
83 G4double incomingParticleEnergy, G4int shell, G4double random) ;
84
85 inline void SelectFasterComputation(G4bool input);
86
89
90protected:
92
93private:
94 G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition * aParticleDefinition,
95 G4double incomingParticleEnergy, G4int shell) ;
96 G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(G4ParticleDefinition * aParticleDefinition,
97 G4double incomingParticleEnergy, G4int shell) ;
98 G4double Interpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
99 G4double QuadInterpolator( G4double e11,
100 G4double e12,
101 G4double e21,
102 G4double e22,
103 G4double x11,
104 G4double x12,
105 G4double x21,
106 G4double x22,
107 G4double t1,
108 G4double t2,
109 G4double t,
110 G4double e);
111 // Partial cross section
112 G4int RandomSelect(G4double energy,const G4String& particle );
113
114 //deexcitation manager to produce fluo photns and e-
115 G4VAtomDeexcitation* fAtomDeexcitation;
116 G4Material* nistSi;
117 std::map<G4String,G4double,std::less<G4String> > lowEnergyLimit;
118 std::map<G4String,G4double,std::less<G4String> > highEnergyLimit;
119
120 // Cross section
121 typedef std::map<G4String,G4String,std::less<G4String> > MapFile;
122 MapFile tableFile;
123
124 typedef std::map<G4String,G4MicroElecCrossSectionDataSet*,std::less<G4String> > MapData;
125 MapData tableData;
126
127 typedef std::map<G4double, std::map<G4double, G4double> > TriDimensionMap;
128 TriDimensionMap eDiffCrossSectionData[7];
129 TriDimensionMap eNrjTransfData[7]; // for cumulated dcs
130 TriDimensionMap pDiffCrossSectionData[7];
131 TriDimensionMap pNrjTransfData[7]; // for cumulated dcs
132 std::vector<G4double> eTdummyVec;
133 std::vector<G4double> pTdummyVec;
134
135 typedef std::map<G4double, std::vector<G4double> > VecMap;
136 VecMap eVecm;
137 VecMap pVecm;
138 VecMap eProbaShellMap[7]; // for cumulated dcs
139 VecMap pProbaShellMap[7]; // for cumulated dcs
140
141 // Final state
142 G4MicroElecSiStructure SiStructure;
143
144 G4int verboseLevel;
145 G4bool isInitialised;
146 G4bool fasterCode;
147 //
148};
149
150//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
151
153{
154 fasterCode = input;
155}
156
157//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
158
159#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4MicroElecInelasticModel(const G4MicroElecInelasticModel &)=delete
G4double TransferedEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell, G4double random)
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4ParticleChangeForGamma * fParticleChangeForGamma
G4MicroElecInelasticModel & operator=(const G4MicroElecInelasticModel &right)=delete
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
G4double DifferentialCrossSection(G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)