Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ElectroNuclearCrossSection.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 // GEANT4 tag $Name: not supported by cvs2svn $
27 //
28 // GEANT4 physics class: G4ElectroNuclearCrossSection -- header file
29 // M.V. Kossov, ITEP(Moscow), 24-OCT-01
30 // The last update: M.V. Kossov, CERN/ITEP (Moscow) 25-Sept-03
31 
32 #ifndef G4ElectroNuclearCrossSection_h
33 #define G4ElectroNuclearCrossSection_h 1
34 
36 #include "G4DynamicParticle.hh"
37 #include "G4Element.hh"
38 #include "G4ParticleTable.hh"
39 #include "G4NucleiProperties.hh"
40 #include <vector>
41 #include "Randomize.hh"
42 #include "G4Electron.hh"
43 #include "G4Positron.hh"
44 
46 {
47 public:
48 
49  G4ElectroNuclearCrossSection(const G4String& name = "ElectroNuclearXS");
51 
52  virtual void CrossSectionDescription(std::ostream&) const;
53 
54  virtual G4bool
55  IsIsoApplicable(const G4DynamicParticle* aParticle, G4int /*Z*/,
56  G4int /*A*/, const G4Element*, const G4Material*);
57 
58  virtual G4double
59  GetIsoCrossSection(const G4DynamicParticle* aParticle,
60  G4int /*Z*/, G4int /*A*/,
61  const G4Isotope*, const G4Element*, const G4Material*);
62 
64 
66 
68 
69 private:
70  G4int GetFunctions(G4double a, G4double* x, G4double* y, G4double* z);
71 
72  G4double ThresholdEnergy(G4int Z, G4int N);
73  G4double HighEnergyJ1(G4double lE);
74  G4double HighEnergyJ2(G4double lE);
75  G4double HighEnergyJ3(G4double lE);
76  G4double SolveTheEquation(G4double f);
77  G4double Fun(G4double x);
78  G4double DFun(G4double x);
79 
80 // Body
81 private:
82  static G4int lastN; // The last N of calculated nucleus
83  static G4int lastZ; // The last Z of calculated nucleus
84  static G4int lastF; // Last used in the cross section TheFirstBin
85  static G4double* lastJ1; // Pointer to the last array of the J1 function
86  static G4double* lastJ2; // Pointer to the last array of the J2 function
87  static G4double* lastJ3; // Pointer to the last array of the J3 function
88  static G4int lastL; // Last used in the cross section TheLastBin
89  static G4double lastE; // Last used in the cross section Energy
90  static G4double lastTH; // Last value of the Energy Threshold
91  static G4double lastSig; // Last value of the Cross Section
92  static G4double lastG; // Last value of gamma=lnE-ln(me)
93  static G4double lastH; // Last value of the High energy A-dependence
94 
95  // Vector of pointers to the J1 tabulated functions
96  static std::vector <G4double*> J1;
97 
98  // Vector of pointers to the J2 tabulated functions
99  static std::vector <G4double*> J2;
100 
101  // Vector of pointers to the J3 tabulated functions
102  static std::vector <G4double*> J3;
103 };
104 
105 
106 inline G4double
107 G4ElectroNuclearCrossSection::DFun(G4double x)
108 {
109  // Parametrization of the PhotoNucCS
110  static const G4double shd=1.0734; // HE PomShadowing(D)
111  static const G4double poc=0.0375; // HE Pomeron coefficient
112  static const G4double pos=16.5; // HE Pomeron shift
113  static const G4double reg=.11; // HE Reggeon slope
114  static const G4double mel=0.5109989; // Mass of an electron in MeV
115  static const G4double lmel=std::log(mel); // Log of an electron mass
116  G4double y=std::exp(x-lastG-lmel); // y for the x
117  G4double flux=lastG*(2.-y*(2.-y))-1.; // flux factor
118  return (poc*(x-pos)+shd*std::exp(-reg*x))*flux;
119 }
120 
121 
122 inline G4double
123 G4ElectroNuclearCrossSection::Fun(G4double x)
124 {
125  // Integrated PhoNuc cross section
126  G4double dlg1=lastG+lastG-1.;
127  G4double lgoe=lastG/lastE;
128  G4double HE2=HighEnergyJ2(x);
129  return dlg1*HighEnergyJ1(x)-lgoe*(HE2+HE2-HighEnergyJ3(x)/lastE);
130 }
131 
132 
133 inline G4double
134 G4ElectroNuclearCrossSection::HighEnergyJ1(G4double lEn)
135 {
136  static const G4double le=std::log(50000.); // std::log(E0)
137  static const G4double le2=le*le; // std::log(E0)^2
138  static const G4double a=.0375; // a
139  static const G4double ha=a*.5; // a/2
140  static const G4double ab=a*16.5; // a*b
141  static const G4double d=0.11; // d
142  static const G4double cd=1.0734/d; // c/d
143  static const G4double ele=std::exp(-d*le); // E0^(-d)
144  return ha*(lEn*lEn-le2)-ab*(lEn-le)-cd*(std::exp(-d*lEn)-ele);
145 }
146 
147 
148 inline G4double
149 G4ElectroNuclearCrossSection::HighEnergyJ2(G4double lEn)
150 {
151  static const G4double e=50000.; // E0
152  static const G4double le=std::log(e); // std::log(E0)
153  static const G4double le1=(le-1.)*e; // (std::log(E0)-1)*E0
154  static const G4double a=.0375; // a
155  static const G4double ab=a*16.5; // a*b
156  static const G4double d=1.-0.11; // 1-d
157  static const G4double cd=1.0734/d; // c/(1-d)
158  static const G4double ele=std::exp(d*le); // E0^(1-d)
159  G4double En=std::exp(lEn);
160  return a*((lEn-1.)*En-le1)-ab*(En-e)+cd*(std::exp(d*lEn)-ele);
161 }
162 
163 
164 inline G4double
165 G4ElectroNuclearCrossSection::HighEnergyJ3(G4double lEn)
166 {
167  static const G4double e=50000.; // E0
168  static const G4double le=std::log(e); // std::log(E0)
169  static const G4double e2=e*e; // E0^2
170  static const G4double leh=(le-.5)*e2; // (std::log(E0)-.5)*E0^2
171  static const G4double ha=.0375*.5; // a/2
172  static const G4double hab=ha*16.5; // a*b/2
173  static const G4double d=2.-.11; // 2-d
174  static const G4double cd=1.0734/d; // c/(2-d)
175  static const G4double ele=std::exp(d*le); // E0^(2-d)
176  G4double lastE2=std::exp(lEn+lEn);
177  return ha*((lEn-.5)*lastE2-leh)-hab*(lastE2-e2)+cd*(std::exp(d*lEn)-ele);
178 }
179 
180 #endif