[1058] | 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 | // |
---|
[1196] | 26 | // $Id: G4PairProductionRelModel.hh,v 1.3 2009/06/04 13:45:53 gunter Exp $ |
---|
[1228] | 27 | // GEANT4 tag $Name: geant4-09-03 $ |
---|
[1058] | 28 | // |
---|
| 29 | // ------------------------------------------------------------------- |
---|
| 30 | // |
---|
| 31 | // GEANT4 Class header file |
---|
| 32 | // |
---|
| 33 | // |
---|
| 34 | // File name: G4PairProductionRelModel |
---|
| 35 | // |
---|
| 36 | // Author: Andreas Schaelicke |
---|
| 37 | // |
---|
| 38 | // Creation date: 02.04.2009 |
---|
| 39 | // |
---|
| 40 | // Modifications: |
---|
| 41 | // |
---|
| 42 | // Class Description: |
---|
| 43 | // |
---|
| 44 | // Implementation of gamma convertion to e+e- in the field of a nucleus |
---|
| 45 | // relativistic approximation |
---|
| 46 | // |
---|
| 47 | |
---|
| 48 | // ------------------------------------------------------------------- |
---|
| 49 | // |
---|
| 50 | |
---|
| 51 | #ifndef G4PairProductionRelModel_h |
---|
| 52 | #define G4PairProductionRelModel_h 1 |
---|
| 53 | |
---|
| 54 | #include "G4VEmModel.hh" |
---|
| 55 | #include "G4PhysicsTable.hh" |
---|
| 56 | #include "G4NistManager.hh" |
---|
| 57 | |
---|
| 58 | class G4ParticleChangeForGamma; |
---|
| 59 | |
---|
| 60 | class G4PairProductionRelModel : public G4VEmModel |
---|
| 61 | { |
---|
| 62 | |
---|
| 63 | public: |
---|
| 64 | |
---|
| 65 | G4PairProductionRelModel(const G4ParticleDefinition* p = 0, |
---|
| 66 | const G4String& nam = "Bethe-Heitler"); |
---|
| 67 | |
---|
| 68 | virtual ~G4PairProductionRelModel(); |
---|
| 69 | |
---|
| 70 | virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&); |
---|
| 71 | |
---|
| 72 | virtual G4double ComputeCrossSectionPerAtom( |
---|
| 73 | const G4ParticleDefinition*, |
---|
| 74 | G4double kinEnergy, |
---|
| 75 | G4double Z, |
---|
| 76 | G4double A=0., |
---|
| 77 | G4double cut=0., |
---|
| 78 | G4double emax=DBL_MAX); |
---|
| 79 | |
---|
| 80 | virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*, |
---|
| 81 | const G4MaterialCutsCouple*, |
---|
| 82 | const G4DynamicParticle*, |
---|
| 83 | G4double tmin, |
---|
| 84 | G4double maxEnergy); |
---|
| 85 | |
---|
| 86 | virtual void SetupForMaterial(const G4ParticleDefinition*, |
---|
| 87 | const G4Material*,G4double); |
---|
| 88 | |
---|
| 89 | // * fast inline functions * |
---|
| 90 | inline void SetCurrentElement(const G4double); |
---|
| 91 | |
---|
| 92 | // set / get methods |
---|
| 93 | inline void SetLPMconstant(G4double val); |
---|
| 94 | inline G4double LPMconstant() const; |
---|
| 95 | |
---|
| 96 | inline void SetLPMflag(G4bool); |
---|
| 97 | inline G4bool LPMflag() const; |
---|
| 98 | |
---|
| 99 | protected: |
---|
| 100 | // screening functions |
---|
| 101 | inline G4double Phi1(G4double delta) const; |
---|
| 102 | inline G4double Phi2(G4double delta) const; |
---|
| 103 | inline G4double DeltaMax() const; |
---|
| 104 | inline G4double DeltaMin(G4double) const; |
---|
| 105 | // lpm functions |
---|
| 106 | void CalcLPMFunctions(G4double k, G4double eplus); |
---|
| 107 | |
---|
| 108 | // obsolete |
---|
| 109 | G4double ScreenFunction1(G4double ScreenVariable); |
---|
| 110 | G4double ScreenFunction2(G4double ScreenVariable); |
---|
| 111 | |
---|
| 112 | |
---|
| 113 | |
---|
| 114 | |
---|
| 115 | G4double ComputeXSectionPerAtom(G4double totalEnergy, G4double Z); |
---|
| 116 | |
---|
| 117 | G4double ComputeDXSectionPerAtom(G4double eplusEnergy, G4double totalEnergy, G4double Z); |
---|
| 118 | G4double ComputeRelDXSectionPerAtom(G4double eplusEnergy, G4double totalEnergy, G4double Z); |
---|
| 119 | |
---|
| 120 | // hide assignment operator |
---|
| 121 | G4PairProductionRelModel & operator=(const G4PairProductionRelModel &right); |
---|
| 122 | G4PairProductionRelModel(const G4PairProductionRelModel&); |
---|
| 123 | |
---|
| 124 | G4NistManager* nist; |
---|
| 125 | |
---|
| 126 | G4ParticleDefinition* theGamma; |
---|
| 127 | G4ParticleDefinition* theElectron; |
---|
| 128 | G4ParticleDefinition* thePositron; |
---|
| 129 | G4ParticleChangeForGamma* fParticleChange; |
---|
| 130 | G4PhysicsTable* theCrossSectionTable; |
---|
| 131 | |
---|
| 132 | G4double lowGammaEnergy; |
---|
| 133 | G4double highGammaEnergy; |
---|
| 134 | |
---|
| 135 | G4int nbins; |
---|
| 136 | size_t indexZ[120]; |
---|
| 137 | |
---|
| 138 | G4double fLPMconstant; |
---|
| 139 | G4bool fLPMflag; |
---|
| 140 | |
---|
| 141 | // cash |
---|
| 142 | G4double z13, z23, lnZ; |
---|
| 143 | G4double Fel, Finel, fCoulomb; |
---|
| 144 | G4double currentZ; |
---|
| 145 | |
---|
| 146 | // LPM effect |
---|
| 147 | G4double lpmEnergy; |
---|
| 148 | G4double xiLPM, phiLPM, gLPM; |
---|
| 149 | |
---|
| 150 | // consts |
---|
| 151 | G4bool use_completescreening; |
---|
| 152 | |
---|
| 153 | static const G4double xgi[8], wgi[8]; |
---|
| 154 | static const G4double Fel_light[5]; |
---|
| 155 | static const G4double Finel_light[5]; |
---|
| 156 | static const G4double facFel; |
---|
| 157 | static const G4double facFinel; |
---|
| 158 | |
---|
| 159 | static const G4double preS1, logTwo; |
---|
| 160 | |
---|
| 161 | }; |
---|
| 162 | |
---|
| 163 | |
---|
| 164 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 165 | |
---|
| 166 | inline |
---|
| 167 | void G4PairProductionRelModel::SetLPMconstant(G4double val) |
---|
| 168 | { |
---|
| 169 | fLPMconstant = val; |
---|
| 170 | } |
---|
| 171 | |
---|
| 172 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 173 | |
---|
| 174 | inline |
---|
| 175 | G4double G4PairProductionRelModel::LPMconstant() const |
---|
| 176 | { |
---|
| 177 | return fLPMconstant; |
---|
| 178 | } |
---|
| 179 | |
---|
| 180 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 181 | |
---|
| 182 | inline |
---|
| 183 | void G4PairProductionRelModel::SetLPMflag(G4bool val) |
---|
| 184 | { |
---|
| 185 | fLPMflag = val; |
---|
| 186 | } |
---|
| 187 | |
---|
| 188 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 189 | |
---|
| 190 | inline |
---|
| 191 | G4bool G4PairProductionRelModel::LPMflag() const |
---|
| 192 | { |
---|
| 193 | return fLPMflag; |
---|
| 194 | } |
---|
| 195 | |
---|
| 196 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 197 | |
---|
| 198 | inline void G4PairProductionRelModel::SetCurrentElement(const G4double Z) |
---|
| 199 | { |
---|
| 200 | if(Z != currentZ) { |
---|
| 201 | currentZ = Z; |
---|
| 202 | |
---|
| 203 | G4int iz = G4int(Z); |
---|
| 204 | z13 = nist->GetZ13(iz); |
---|
| 205 | z23 = z13*z13; |
---|
| 206 | lnZ = nist->GetLOGZ(iz); |
---|
| 207 | |
---|
| 208 | if (iz <= 4) { |
---|
| 209 | Fel = Fel_light[iz]; |
---|
| 210 | Finel = Finel_light[iz] ; |
---|
| 211 | } |
---|
| 212 | else { |
---|
| 213 | Fel = facFel - lnZ/3. ; |
---|
| 214 | Finel = facFinel - 2.*lnZ/3. ; |
---|
| 215 | } |
---|
| 216 | |
---|
| 217 | fCoulomb=GetCurrentElement()->GetfCoulomb(); |
---|
| 218 | } |
---|
| 219 | } |
---|
| 220 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 221 | |
---|
| 222 | inline G4double G4PairProductionRelModel::Phi1(G4double delta) const |
---|
| 223 | { |
---|
| 224 | G4double screenVal; |
---|
| 225 | |
---|
| 226 | if (delta > 1.) |
---|
| 227 | screenVal = 21.12 - 4.184*std::log(delta+0.952); |
---|
| 228 | else |
---|
| 229 | screenVal = 20.868 - delta*(3.242 - 0.625*delta); |
---|
| 230 | |
---|
| 231 | return screenVal; |
---|
| 232 | } |
---|
| 233 | |
---|
| 234 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 235 | |
---|
| 236 | inline G4double G4PairProductionRelModel::Phi2(G4double delta) const |
---|
| 237 | { |
---|
| 238 | G4double screenVal; |
---|
| 239 | |
---|
| 240 | if (delta > 1.) |
---|
| 241 | screenVal = 21.12 - 4.184*std::log(delta+0.952); |
---|
| 242 | else |
---|
| 243 | screenVal = 20.209 - delta*(1.930 + 0.086*delta); |
---|
| 244 | |
---|
| 245 | return screenVal; |
---|
| 246 | } |
---|
| 247 | |
---|
| 248 | inline G4double G4PairProductionRelModel::ScreenFunction1(G4double ScreenVariable) |
---|
| 249 | |
---|
| 250 | // compute the value of the screening function 3*PHI1 - PHI2 |
---|
| 251 | |
---|
| 252 | { |
---|
| 253 | G4double screenVal; |
---|
| 254 | |
---|
| 255 | if (ScreenVariable > 1.) |
---|
| 256 | screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952); |
---|
| 257 | else |
---|
| 258 | screenVal = 42.392 - ScreenVariable*(7.796 - 1.961*ScreenVariable); |
---|
| 259 | |
---|
| 260 | return screenVal; |
---|
| 261 | } |
---|
| 262 | |
---|
| 263 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 264 | |
---|
| 265 | inline G4double G4PairProductionRelModel::ScreenFunction2(G4double ScreenVariable) |
---|
| 266 | |
---|
| 267 | // compute the value of the screening function 1.5*PHI1 + 0.5*PHI2 |
---|
| 268 | |
---|
| 269 | { |
---|
| 270 | G4double screenVal; |
---|
| 271 | |
---|
| 272 | if (ScreenVariable > 1.) |
---|
| 273 | screenVal = 42.24 - 8.368*std::log(ScreenVariable+0.952); |
---|
| 274 | else |
---|
| 275 | screenVal = 41.405 - ScreenVariable*(5.828 - 0.8945*ScreenVariable); |
---|
| 276 | |
---|
| 277 | return screenVal; |
---|
| 278 | } |
---|
| 279 | |
---|
| 280 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 281 | |
---|
| 282 | |
---|
| 283 | inline G4double G4PairProductionRelModel::DeltaMax() const |
---|
| 284 | { |
---|
| 285 | // k > 50 MeV |
---|
| 286 | G4double FZ = 8.*(lnZ/3. + fCoulomb); |
---|
[1196] | 287 | return std::exp( (42.24-FZ)/8.368 ) + 0.952; |
---|
[1058] | 288 | } |
---|
| 289 | |
---|
| 290 | inline G4double G4PairProductionRelModel::DeltaMin(G4double k) const |
---|
| 291 | { |
---|
| 292 | return 4.*136./z13*(electron_mass_c2/k); |
---|
| 293 | } |
---|
| 294 | |
---|
| 295 | |
---|
| 296 | |
---|
| 297 | #endif |
---|