Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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"
47 #include "G4LogLogInterpolation.hh"
48 #include "G4ProductionCutsTable.hh"
49 #include "G4NistManager.hh"
50 
52 {
53 
54 public:
55 
57  const G4String& nam = "MicroElecElasticModel");
58 
59  virtual ~G4MicroElecElasticModel();
60 
61  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
62 
63  virtual G4double CrossSectionPerVolume(const G4Material* material,
64  const G4ParticleDefinition* p,
65  G4double ekin,
66  G4double emin,
67  G4double emax);
68 
69  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
70  const G4MaterialCutsCouple*,
71  const G4DynamicParticle*,
72  G4double tmin,
73  G4double maxEnergy);
74 
75  inline void SetKillBelowThreshold (G4double threshold);
76  G4double GetKillBelowThreshold () { return killBelowEnergy; }
77 
78 protected:
79 
81 
82 private:
83 
84  G4Material* nistSi;
85  G4double killBelowEnergy;
86  G4double lowEnergyLimit;
87  G4double lowEnergyLimitOfModel;
88  G4double highEnergyLimit;
89  G4bool isInitialised;
90  G4int verboseLevel;
91 
92  // Cross section
93 
94  typedef std::map<G4String,G4String,std::less<G4String> > MapFile;
95  MapFile tableFile;
96 
97  typedef std::map<G4String,G4MicroElecCrossSectionDataSet*,std::less<G4String> > MapData;
98  MapData tableData;
99 
100  // Final state
101 
102  G4double Theta(G4ParticleDefinition * aParticleDefinition, G4double k, G4double integrDiff);
103 
104  G4double LinLinInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
105 
106  G4double LogLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
107 
108  G4double LinLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
109 
110  G4double QuadInterpolator(G4double e11,
111  G4double e12,
112  G4double e21,
113  G4double e22,
114  G4double x11,
115  G4double x12,
116  G4double x21,
117  G4double x22,
118  G4double t1,
119  G4double t2,
120  G4double t,
121  G4double e);
122 
123  typedef std::map<double, std::map<double, double> > TriDimensionMap;
124 
125  TriDimensionMap eDiffCrossSectionData;
126  std::vector<double> eTdummyVec;
127 
128  typedef std::map<double, std::vector<double> > VecMap;
129  VecMap eVecm;
130 
131  G4double RandomizeCosTheta(G4double k);
132 
133  //
134 
135  G4MicroElecElasticModel & operator=(const G4MicroElecElasticModel &right);
137 
138 };
139 
141 {
142  killBelowEnergy = threshold;
143 
144  if (threshold < 5*CLHEP::eV)
145  {
146  G4Exception ("*** WARNING : the G4MicroElecElasticModel class is not validated below 5 eV !","",JustWarning,"") ;
147  threshold = 5*CLHEP::eV;
148  }
149 
150 }
151 
152 
153 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
154 
155 #endif
G4MicroElecElasticModel(const G4ParticleDefinition *p=0, const G4String &nam="MicroElecElasticModel")
const char * p
Definition: xmltok.h:285
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
int G4int
Definition: G4Types.hh:78
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
bool G4bool
Definition: G4Types.hh:79
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
static constexpr double eV
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const G4double emax
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetKillBelowThreshold(G4double threshold)
double G4double
Definition: G4Types.hh:76