Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4QElectronNuclearCrossSection.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// GEANT4 tag $Name: not supported by cvs2svn $
28//
29//
30// GEANT4 physics class: G4QElectronNuclearCrossSection -- header file
31// M.V. Kossov, ITEP(Moscow), 24-OCT-01
32// The last update: M.V. Kossov, CERN/ITEP (Moscow) 25-Sept-03
33//
34// --------------------------------------------------------------------------------
35// Short description: reaction cross-sections for electron-nuclear reactions, which
36// are integrals over virtual equivalent photons photons.
37// --------------------------------------------------------------------------------
38
39#ifndef G4QElectronNuclearCrossSection_h
40#define G4QElectronNuclearCrossSection_h 1
41
42#include <vector>
43#include "Randomize.hh"
44#include "G4VQCrossSection.hh"
45
47{
48protected:
49
51
52public:
53
55
56 static G4VQCrossSection* GetPointer(); // Gives a pointer to this singletone
57
59
60 // At present momentum (pMom) must be in GeV (@@ Units)
61 virtual G4double GetCrossSection(G4bool fCS, G4double pMom, G4int tgZ, G4int tgN,
62 G4int pPDG=0);
63
65 G4double Momentum);
66
68
70
72
74
75private:
76 G4int GetFunctions(G4double a, G4double* x, G4double* y, G4double* z);
77 G4double HighEnergyJ1(G4double lE);
78 G4double HighEnergyJ2(G4double lE);
79 G4double HighEnergyJ3(G4double lE);
80 G4double SolveTheEquation(G4double f);
81 G4double Fun(G4double x);
82 G4double DFun(G4double x);
83
84// Body
85private:
86 static G4bool onlyCS; // flag to calculate only CrossSection
87 static G4double lastSig; // Last calculated cross section
88 static G4int lastL; // Last used in the cross section TheLastBin
89 static G4double lastE; // Last used in the cross section Energy
90 static G4int lastF; // Last used in the cross section TheFirstBin
91 static G4double lastG; // Last value of gamma=lnE-ln(m)
92 static G4double lastH; // Last value of the High energy A-dependence
93 static G4double* lastJ1; // Pointer to the last array of the J1 function
94 static G4double* lastJ2; // Pointer to the last array of the J2 function
95 static G4double* lastJ3; // Pointer to the last array of the J3 function
96 static G4int lastPDG; // The last projectile PDG
97 static G4int lastN; // The last N of calculated nucleus
98 static G4int lastZ; // The last Z of calculated nucleus
99 static G4double lastP; // Last used in the cross section Momentum
100 static G4double lastTH; // Last value of the Momentum Threshold
101 static G4double lastCS; // Last value of the Cross Section
102 static G4int lastI; // The last position in the DAMDB
103 static std::vector <G4double*>* J1; // Vector of pointers to the J1 tabulated functions
104 static std::vector <G4double*>* J2; // Vector of pointers to the J2 tabulated functions
105 static std::vector <G4double*>* J3; // Vector of pointers to the J3 tabulated functions
106};
107
108inline G4double G4QElectronNuclearCrossSection::DFun(G4double x)//PhotoNucCSParametrization
109{
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
121inline G4double G4QElectronNuclearCrossSection::Fun(G4double x) // Integrated PhoNucCS
122{
123 G4double dlg1=lastG+lastG-1.;
124 G4double lgoe=lastG/lastE;
125 G4double HE2=HighEnergyJ2(x);
126 return dlg1*HighEnergyJ1(x)-lgoe*(HE2+HE2-HighEnergyJ3(x)/lastE);
127}
128
129inline G4double G4QElectronNuclearCrossSection::HighEnergyJ1(G4double lEn)
130{
131 static const G4double le=std::log(50000.); // log(E0)
132 static const G4double le2=le*le; // log(E0)^2
133 static const G4double a=.0375; // a
134 static const G4double ha=a*.5; // a/2
135 static const G4double ab=a*16.5; // a*b
136 static const G4double d=0.11; // d
137 static const G4double cd=1.0734/d; // c/d
138 static const G4double ele=std::exp(-d*le); // E0^(-d)
139 return ha*(lEn*lEn-le2)-ab*(lEn-le)-cd*(std::exp(-d*lEn)-ele);
140}
141
142inline G4double G4QElectronNuclearCrossSection::HighEnergyJ2(G4double lEn)
143{
144 static const G4double e=50000.; // E0
145 static const G4double le=std::log(e); // log(E0)
146 static const G4double le1=(le-1.)*e; // (log(E0)-1)*E0
147 static const G4double a=.0375; // a
148 static const G4double ab=a*16.5; // a*b
149 static const G4double d=1.-0.11; // 1-d
150 static const G4double cd=1.0734/d; // c/(1-d)
151 static const G4double ele=std::exp(d*le); // E0^(1-d)
152 G4double En=std::exp(lEn);
153 return a*((lEn-1.)*En-le1)-ab*(En-e)+cd*(std::exp(d*lEn)-ele);
154}
155
156inline G4double G4QElectronNuclearCrossSection::HighEnergyJ3(G4double lEn)
157{
158 static const G4double e=50000.; // E0
159 static const G4double le=std::log(e); // log(E0)
160 static const G4double e2=e*e; // E0^2
161 static const G4double leh=(le-.5)*e2; // (log(E0)-.5)*E0^2
162 static const G4double ha=.0375*.5; // a/2
163 static const G4double hab=ha*16.5; // a*b/2
164 static const G4double d=2.-.11; // 2-d
165 static const G4double cd=1.0734/d; // c/(2-d)
166 static const G4double ele=std::exp(d*le); // E0^(2-d)
167 G4double lastE2=std::exp(lEn+lEn);
168 return ha*((lEn-.5)*lastE2-leh)-hab*(lastE2-e2)+cd*(std::exp(d*lEn)-ele);
169}
170
171#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4double CalculateCrossSection(G4bool CS, G4int F, G4int I, G4int PDG, G4int Z, G4int N, G4double Momentum)
virtual G4double GetCrossSection(G4bool fCS, G4double pMom, G4int tgZ, G4int tgN, G4int pPDG=0)
G4double GetVirtualFactor(G4double nu, G4double Q2)
G4double ThresholdEnergy(G4int Z, G4int N, G4int PDG=11)