[819] | 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 | // |
---|
[1007] | 27 | // GEANT4 tag $Name: geant4-09-02 $ |
---|
[819] | 28 | // |
---|
| 29 | // |
---|
| 30 | // GEANT4 physics class: G4ElectroNuclearCrossSection -- 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 | #ifndef G4ElectroNuclearCrossSection_h |
---|
| 36 | #define G4ElectroNuclearCrossSection_h 1 |
---|
| 37 | |
---|
| 38 | #include "G4VCrossSectionDataSet.hh" |
---|
| 39 | #include "G4DynamicParticle.hh" |
---|
| 40 | #include "G4Element.hh" |
---|
| 41 | #include "G4ParticleTable.hh" |
---|
| 42 | #include "G4NucleiProperties.hh" |
---|
| 43 | #include <vector> |
---|
| 44 | #include "Randomize.hh" |
---|
| 45 | #include "G4Electron.hh" |
---|
| 46 | #include "G4Positron.hh" |
---|
| 47 | |
---|
| 48 | class G4ElectroNuclearCrossSection : public G4VCrossSectionDataSet |
---|
| 49 | { |
---|
| 50 | public: |
---|
| 51 | |
---|
| 52 | G4ElectroNuclearCrossSection(); |
---|
| 53 | virtual ~G4ElectroNuclearCrossSection(); |
---|
| 54 | |
---|
| 55 | G4bool IsApplicable(const G4DynamicParticle* aParticle, const G4Element* ) |
---|
| 56 | { |
---|
| 57 | return IsZAApplicable(aParticle, 0., 0.); |
---|
| 58 | } |
---|
| 59 | |
---|
| 60 | G4bool IsZAApplicable(const G4DynamicParticle* aParticle, G4double /*ZZ*/, |
---|
| 61 | G4double /*AA*/) |
---|
| 62 | { |
---|
| 63 | G4bool result = false; |
---|
| 64 | if( aParticle->GetDefinition()==G4Electron::ElectronDefinition()) result = true; |
---|
| 65 | if( aParticle->GetDefinition()==G4Positron::PositronDefinition()) result = true; |
---|
| 66 | return result; |
---|
| 67 | } |
---|
| 68 | |
---|
| 69 | |
---|
| 70 | G4double GetCrossSection(const G4DynamicParticle* aParticle, |
---|
| 71 | const G4Element* anElement, G4double T=0.); |
---|
| 72 | |
---|
| 73 | G4double GetIsoZACrossSection(const G4DynamicParticle* aParticle, |
---|
| 74 | G4double ZZ, G4double AA, G4double T=0.); |
---|
| 75 | |
---|
| 76 | |
---|
| 77 | G4double GetEquivalentPhotonEnergy(); |
---|
| 78 | |
---|
| 79 | G4double GetVirtualFactor(G4double nu, G4double Q2); |
---|
| 80 | |
---|
| 81 | G4double GetEquivalentPhotonQ2(G4double nu); |
---|
| 82 | |
---|
| 83 | void BuildPhysicsTable(const G4ParticleDefinition&) {} |
---|
| 84 | |
---|
| 85 | void DumpPhysicsTable(const G4ParticleDefinition&) {} |
---|
| 86 | |
---|
| 87 | private: |
---|
| 88 | G4int GetFunctions(G4double a, G4double* x, G4double* y, G4double* z); |
---|
| 89 | //G4double LinearFit(G4double X, G4int N, const G4double* XN, const G4double* YN); |
---|
| 90 | G4double ThresholdEnergy(G4int Z, G4int N); |
---|
| 91 | G4double HighEnergyJ1(G4double lE); |
---|
| 92 | G4double HighEnergyJ2(G4double lE); |
---|
| 93 | G4double HighEnergyJ3(G4double lE); |
---|
| 94 | G4double SolveTheEquation(G4double f); |
---|
| 95 | G4double Fun(G4double x); |
---|
| 96 | G4double DFun(G4double x); |
---|
| 97 | |
---|
| 98 | // Body |
---|
| 99 | private: |
---|
| 100 | static G4int lastN; // The last N of calculated nucleus |
---|
| 101 | static G4int lastZ; // The last Z of calculated nucleus |
---|
| 102 | static G4int lastF; // Last used in the cross section TheFirstBin |
---|
| 103 | static G4double* lastJ1; // Pointer to the last array of the J1 function |
---|
| 104 | static G4double* lastJ2; // Pointer to the last array of the J2 function |
---|
| 105 | static G4double* lastJ3; // Pointer to the last array of the J3 function |
---|
| 106 | static G4int lastL; // Last used in the cross section TheLastBin |
---|
| 107 | static G4double lastE; // Last used in the cross section Energy |
---|
| 108 | static G4double lastTH; // Last value of the Energy Threshold |
---|
| 109 | static G4double lastSig; // Last value of the Cross Section |
---|
| 110 | static G4double lastG; // Last value of gamma=lnE-ln(me) |
---|
| 111 | static G4double lastH; // Last value of the High energy A-dependence |
---|
| 112 | |
---|
| 113 | static std::vector <G4double*> J1; // Vector of pointers to the J1 tabulated functions |
---|
| 114 | static std::vector <G4double*> J2; // Vector of pointers to the J2 tabulated functions |
---|
| 115 | static std::vector <G4double*> J3; // Vector of pointers to the J3 tabulated functions |
---|
| 116 | }; |
---|
| 117 | |
---|
| 118 | inline G4double G4ElectroNuclearCrossSection::DFun(G4double x)// Parametrization of the PhotoNucCS |
---|
| 119 | { |
---|
| 120 | static const G4double shd=1.0734; // HE PomShadowing(D) |
---|
| 121 | static const G4double poc=0.0375; // HE Pomeron coefficient |
---|
| 122 | static const G4double pos=16.5; // HE Pomeron shift |
---|
| 123 | static const G4double reg=.11; // HE Reggeon slope |
---|
| 124 | static const G4double mel=0.5109989; // Mass of an electron in MeV |
---|
| 125 | static const G4double lmel=std::log(mel); // Log of an electron mass |
---|
| 126 | G4double y=std::exp(x-lastG-lmel); // y for the x |
---|
| 127 | G4double flux=lastG*(2.-y*(2.-y))-1.; // flux factor |
---|
| 128 | return (poc*(x-pos)+shd*std::exp(-reg*x))*flux; |
---|
| 129 | } |
---|
| 130 | |
---|
| 131 | inline G4double G4ElectroNuclearCrossSection::Fun(G4double x) // Integrated PhoNuc cross section |
---|
| 132 | { |
---|
| 133 | G4double dlg1=lastG+lastG-1.; |
---|
| 134 | G4double lgoe=lastG/lastE; |
---|
| 135 | G4double HE2=HighEnergyJ2(x); |
---|
| 136 | return dlg1*HighEnergyJ1(x)-lgoe*(HE2+HE2-HighEnergyJ3(x)/lastE); |
---|
| 137 | } |
---|
| 138 | |
---|
| 139 | inline G4double G4ElectroNuclearCrossSection::HighEnergyJ1(G4double lEn) |
---|
| 140 | { |
---|
| 141 | static const G4double le=std::log(50000.); // std::log(E0) |
---|
| 142 | static const G4double le2=le*le; // std::log(E0)^2 |
---|
| 143 | static const G4double a=.0375; // a |
---|
| 144 | static const G4double ha=a*.5; // a/2 |
---|
| 145 | static const G4double ab=a*16.5; // a*b |
---|
| 146 | static const G4double d=0.11; // d |
---|
| 147 | static const G4double cd=1.0734/d; // c/d |
---|
| 148 | static const G4double ele=std::exp(-d*le); // E0^(-d) |
---|
| 149 | return ha*(lEn*lEn-le2)-ab*(lEn-le)-cd*(std::exp(-d*lEn)-ele); |
---|
| 150 | } |
---|
| 151 | |
---|
| 152 | inline G4double G4ElectroNuclearCrossSection::HighEnergyJ2(G4double lEn) |
---|
| 153 | { |
---|
| 154 | static const G4double e=50000.; // E0 |
---|
| 155 | static const G4double le=std::log(e); // std::log(E0) |
---|
| 156 | static const G4double le1=(le-1.)*e; // (std::log(E0)-1)*E0 |
---|
| 157 | static const G4double a=.0375; // a |
---|
| 158 | static const G4double ab=a*16.5; // a*b |
---|
| 159 | static const G4double d=1.-0.11; // 1-d |
---|
| 160 | static const G4double cd=1.0734/d; // c/(1-d) |
---|
| 161 | static const G4double ele=std::exp(d*le); // E0^(1-d) |
---|
| 162 | G4double En=std::exp(lEn); |
---|
| 163 | return a*((lEn-1.)*En-le1)-ab*(En-e)+cd*(std::exp(d*lEn)-ele); |
---|
| 164 | } |
---|
| 165 | |
---|
| 166 | inline G4double G4ElectroNuclearCrossSection::HighEnergyJ3(G4double lEn) |
---|
| 167 | { |
---|
| 168 | static const G4double e=50000.; // E0 |
---|
| 169 | static const G4double le=std::log(e); // std::log(E0) |
---|
| 170 | static const G4double e2=e*e; // E0^2 |
---|
| 171 | static const G4double leh=(le-.5)*e2; // (std::log(E0)-.5)*E0^2 |
---|
| 172 | static const G4double ha=.0375*.5; // a/2 |
---|
| 173 | static const G4double hab=ha*16.5; // a*b/2 |
---|
| 174 | static const G4double d=2.-.11; // 2-d |
---|
| 175 | static const G4double cd=1.0734/d; // c/(2-d) |
---|
| 176 | static const G4double ele=std::exp(d*le); // E0^(2-d) |
---|
| 177 | G4double lastE2=std::exp(lEn+lEn); |
---|
| 178 | return ha*((lEn-.5)*lastE2-leh)-hab*(lastE2-e2)+cd*(std::exp(d*lEn)-ele); |
---|
| 179 | } |
---|
| 180 | |
---|
| 181 | #endif |
---|