Geant4 Cross Reference |
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 // ------------------------------------------------------------------- 28 // 29 // GEANT4 Class file 30 // 31 // 32 // File name: G4BetheHeitlerModel 33 // 34 // Author: Vladimir Ivanchenko on base of Michel Maire code 35 // 36 // Creation date: 15.03.2005 37 // 38 // Modifications by Vladimir Ivanchenko, Michel Maire, Mihaly Novak 39 // 40 // Class Description: 41 // 42 // ------------------------------------------------------------------- 43 // 44 45 #include "G4BetheHeitlerModel.hh" 46 #include "G4PhysicalConstants.hh" 47 #include "G4SystemOfUnits.hh" 48 #include "G4Electron.hh" 49 #include "G4Positron.hh" 50 #include "G4Gamma.hh" 51 #include "Randomize.hh" 52 #include "G4ParticleChangeForGamma.hh" 53 #include "G4Pow.hh" 54 #include "G4Exp.hh" 55 #include "G4ModifiedTsai.hh" 56 #include "G4EmParameters.hh" 57 #include "G4EmElementXS.hh" 58 #include "G4AutoLock.hh" 59 60 const G4int G4BetheHeitlerModel::gMaxZet = 120; 61 std::vector<G4BetheHeitlerModel::ElementData*> G4BetheHeitlerModel::gElementData; 62 63 namespace 64 { 65 G4Mutex theBetheHMutex = G4MUTEX_INITIALIZER; 66 } 67 68 G4BetheHeitlerModel::G4BetheHeitlerModel(const G4ParticleDefinition*, 69 const G4String& nam) 70 : G4VEmModel(nam), 71 fG4Calc(G4Pow::GetInstance()), fTheGamma(G4Gamma::Gamma()), 72 fTheElectron(G4Electron::Electron()), fThePositron(G4Positron::Positron()), 73 fParticleChange(nullptr) 74 { 75 SetAngularDistribution(new G4ModifiedTsai()); 76 } 77 78 G4BetheHeitlerModel::~G4BetheHeitlerModel() 79 { 80 if (isFirstInstance) { 81 for (auto const & ptr : gElementData) { delete ptr; } 82 gElementData.clear(); 83 } 84 delete fXSection; 85 } 86 87 void G4BetheHeitlerModel::Initialise(const G4ParticleDefinition* p, 88 const G4DataVector& cuts) 89 { 90 if (!fParticleChange) { fParticleChange = GetParticleChangeForGamma(); } 91 92 if (isFirstInstance || gElementData.empty()) { 93 G4AutoLock l(&theBetheHMutex); 94 if (gElementData.empty()) { 95 isFirstInstance = true; 96 gElementData.resize(gMaxZet+1, nullptr); 97 98 // EPICS2017 flag should be checked only once 99 useEPICS2017 = G4EmParameters::Instance()->UseEPICS2017XS(); 100 if (useEPICS2017) { 101 fXSection = new G4EmElementXS(1, 100, "convEPICS2017", "/epics2017/pair/pp-cs-"); 102 } 103 } 104 // static data should be initialised only in the one instance 105 InitialiseElementData(); 106 l.unlock(); 107 } 108 // element selectors should be initialised in the master thread 109 if(IsMaster()) { 110 InitialiseElementSelectors(p, cuts); 111 } 112 } 113 114 void G4BetheHeitlerModel::InitialiseLocal(const G4ParticleDefinition*, 115 G4VEmModel* masterModel) 116 { 117 SetElementSelectors(masterModel->GetElementSelectors()); 118 } 119 120 // Calculates the microscopic cross section in GEANT4 internal units. 121 // A parametrized formula from L. Urban is used to estimate 122 // the total cross section. 123 // It gives a good description of the data from 1.5 MeV to 100 GeV. 124 // below 1.5 MeV: sigma=sigma(1.5MeV)*(GammaEnergy-2electronmass) 125 // *(GammaEnergy-2electronmass) 126 G4double 127 G4BetheHeitlerModel::ComputeCrossSectionPerAtom(const G4ParticleDefinition*, 128 G4double gammaEnergy, G4double Z, 129 G4double, G4double, G4double) 130 { 131 G4double xSection = 0.0 ; 132 // short versions 133 static const G4double kMC2 = CLHEP::electron_mass_c2; 134 // zero cross section below the kinematical limit: Eg<2mc^2 135 if (Z < 0.9 || gammaEnergy <= 2.0*kMC2) { return xSection; } 136 137 G4int iZ = G4lrint(Z); 138 if (useEPICS2017 && iZ < 101) { 139 return fXSection->GetXS(iZ, gammaEnergy); 140 } 141 142 // 143 static const G4double gammaEnergyLimit = 1.5*CLHEP::MeV; 144 // set coefficients a, b c 145 static const G4double a0 = 8.7842e+2*CLHEP::microbarn; 146 static const G4double a1 = -1.9625e+3*CLHEP::microbarn; 147 static const G4double a2 = 1.2949e+3*CLHEP::microbarn; 148 static const G4double a3 = -2.0028e+2*CLHEP::microbarn; 149 static const G4double a4 = 1.2575e+1*CLHEP::microbarn; 150 static const G4double a5 = -2.8333e-1*CLHEP::microbarn; 151 152 static const G4double b0 = -1.0342e+1*CLHEP::microbarn; 153 static const G4double b1 = 1.7692e+1*CLHEP::microbarn; 154 static const G4double b2 = -8.2381 *CLHEP::microbarn; 155 static const G4double b3 = 1.3063 *CLHEP::microbarn; 156 static const G4double b4 = -9.0815e-2*CLHEP::microbarn; 157 static const G4double b5 = 2.3586e-3*CLHEP::microbarn; 158 159 static const G4double c0 = -4.5263e+2*CLHEP::microbarn; 160 static const G4double c1 = 1.1161e+3*CLHEP::microbarn; 161 static const G4double c2 = -8.6749e+2*CLHEP::microbarn; 162 static const G4double c3 = 2.1773e+2*CLHEP::microbarn; 163 static const G4double c4 = -2.0467e+1*CLHEP::microbarn; 164 static const G4double c5 = 6.5372e-1*CLHEP::microbarn; 165 // check low energy limit of the approximation (1.5 MeV) 166 G4double gammaEnergyOrg = gammaEnergy; 167 if (gammaEnergy < gammaEnergyLimit) { gammaEnergy = gammaEnergyLimit; } 168 // compute gamma energy variables 169 const G4double x = G4Log(gammaEnergy/kMC2); 170 const G4double x2 = x *x; 171 const G4double x3 = x2*x; 172 const G4double x4 = x3*x; 173 const G4double x5 = x4*x; 174 // 175 const G4double F1 = a0 + a1*x + a2*x2 + a3*x3 + a4*x4 + a5*x5; 176 const G4double F2 = b0 + b1*x + b2*x2 + b3*x3 + b4*x4 + b5*x5; 177 const G4double F3 = c0 + c1*x + c2*x2 + c3*x3 + c4*x4 + c5*x5; 178 // compute the approximated cross section 179 xSection = (Z + 1.)*(F1*Z + F2*Z*Z + F3); 180 // check if we are below the limit of the approximation and apply correction 181 if (gammaEnergyOrg < gammaEnergyLimit) { 182 const G4double dum = (gammaEnergyOrg-2.*kMC2)/(gammaEnergyLimit-2.*kMC2); 183 xSection *= dum*dum; 184 } 185 // make sure that the cross section is never negative 186 xSection = std::max(xSection, 0.); 187 return xSection; 188 } 189 190 // The secondaries e+e- energies are sampled using the Bethe - Heitler 191 // cross sections with Coulomb correction. 192 // A modified version of the random number techniques of Butcher & Messel 193 // is used (Nuc Phys 20(1960),15). 194 // 195 // GEANT4 internal units. 196 // 197 // Note 1 : Effects due to the breakdown of the Born approximation at 198 // low energy are ignored. 199 // Note 2 : The differential cross section implicitly takes account of 200 // pair creation in both nuclear and atomic electron fields. 201 // However triplet prodution is not generated. 202 void G4BetheHeitlerModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect, 203 const G4MaterialCutsCouple* couple, 204 const G4DynamicParticle* aDynamicGamma, 205 G4double, G4double) 206 { 207 // set some constant values 208 const G4double gammaEnergy = aDynamicGamma->GetKineticEnergy(); 209 const G4double eps0 = CLHEP::electron_mass_c2/gammaEnergy; 210 // 211 // check kinematical limit: gamma energy(Eg) must be at least 2 e- rest mass 212 if (eps0 > 0.5) { return; } 213 // 214 // select target element of the material (probs. are based on partial x-secs) 215 const G4Element* anElement = SelectTargetAtom(couple, fTheGamma, gammaEnergy, 216 aDynamicGamma->GetLogKineticEnergy()); 217 218 // 219 // get the random engine 220 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine(); 221 // 222 // 'eps' is the total energy transferred to one of the e-/e+ pair in initial 223 // gamma energy units Eg. Since the corresponding DCS is symmetric on eps=0.5, 224 // the kinematical limits for eps0=mc^2/Eg <= eps <= 0.5 225 // 1. 'eps' is sampled uniformly on the [eps0, 0.5] inteval if Eg<Egsmall 226 // 2. otherwise, on the [eps_min, 0.5] interval according to the DCS (case 2.) 227 G4double eps; 228 // case 1. 229 static const G4double Egsmall = 2.*CLHEP::MeV; 230 if (gammaEnergy < Egsmall) { 231 eps = eps0 + (0.5-eps0)*rndmEngine->flat(); 232 } else { 233 // case 2. 234 // get the Coulomb factor for the target element (Z) and gamma energy (Eg) 235 // F(Z) = 8*ln(Z)/3 if Eg <= 50 [MeV] => no Coulomb correction 236 // F(Z) = 8*ln(Z)/3 + 8*fc(Z) if Eg > 50 [MeV] => fc(Z) is the Coulomb cor. 237 // 238 // The screening variable 'delta(eps)' = 136*Z^{-1/3}*eps0/[eps(1-eps)] 239 // Due to the Coulomb correction, the DCS can go below zero even at 240 // kinematicaly allowed eps > eps0 values. In order to exclude this eps 241 // range with negative DCS, the minimum eps value will be set to eps_min = 242 // max[eps0, epsp] with epsp is the solution of SF(delta(epsp)) - F(Z)/2 = 0 243 // with SF being the screening function (SF1=SF2 at high value of delta). 244 // The solution is epsp = 0.5 - 0.5*sqrt[ 1 - 4*136*Z^{-1/3}eps0/deltap] 245 // with deltap = Exp[(42.038-F(Z))/8.29]-0.958. So the limits are: 246 // - when eps=eps_max = 0.5 => delta_min = 136*Z^{-1/3}*eps0/4 247 // - epsp = 0.5 - 0.5*sqrt[ 1 - delta_min/deltap] 248 // - and eps_min = max[eps0, epsp] 249 static const G4double midEnergy = 50.*CLHEP::MeV; 250 const G4int iZet = std::min(gMaxZet, anElement->GetZasInt()); 251 const G4double deltaFactor = 136.*eps0/anElement->GetIonisation()->GetZ3(); 252 G4double deltaMax = gElementData[iZet]->fDeltaMaxLow; 253 G4double FZ = 8.*anElement->GetIonisation()->GetlogZ3(); 254 if (gammaEnergy > midEnergy) { 255 FZ += 8.*(anElement->GetfCoulomb()); 256 deltaMax = gElementData[iZet]->fDeltaMaxHigh; 257 } 258 const G4double deltaMin = 4.*deltaFactor; 259 // 260 // compute the limits of eps 261 const G4double epsp = 0.5 - 0.5*std::sqrt(1. - deltaMin/deltaMax) ; 262 const G4double epsMin = std::max(eps0,epsp); 263 const G4double epsRange = 0.5 - epsMin; 264 // 265 // sample the energy rate (eps) of the created electron (or positron) 266 G4double F10, F20; 267 ScreenFunction12(deltaMin, F10, F20); 268 F10 -= FZ; 269 F20 -= FZ; 270 const G4double NormF1 = std::max(F10 * epsRange * epsRange, 0.); 271 const G4double NormF2 = std::max(1.5 * F20 , 0.); 272 const G4double NormCond = NormF1/(NormF1 + NormF2); 273 // we will need 3 uniform random number for each trial of sampling 274 G4double rndmv[3]; 275 G4double greject = 0.; 276 do { 277 rndmEngine->flatArray(3, rndmv); 278 if (NormCond > rndmv[0]) { 279 eps = 0.5 - epsRange * fG4Calc->A13(rndmv[1]); 280 const G4double delta = deltaFactor/(eps*(1.-eps)); 281 greject = (ScreenFunction1(delta)-FZ)/F10; 282 } else { 283 eps = epsMin + epsRange*rndmv[1]; 284 const G4double delta = deltaFactor/(eps*(1.-eps)); 285 greject = (ScreenFunction2(delta)-FZ)/F20; 286 } 287 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko 288 } while (greject < rndmv[2]); 289 } // end of eps sampling 290 // 291 // select charges randomly 292 G4double eTotEnergy, pTotEnergy; 293 if (rndmEngine->flat() > 0.5) { 294 eTotEnergy = (1.-eps)*gammaEnergy; 295 pTotEnergy = eps*gammaEnergy; 296 } else { 297 pTotEnergy = (1.-eps)*gammaEnergy; 298 eTotEnergy = eps*gammaEnergy; 299 } 300 // 301 // sample pair kinematics 302 const G4double eKinEnergy = std::max(0.,eTotEnergy - CLHEP::electron_mass_c2); 303 const G4double pKinEnergy = std::max(0.,pTotEnergy - CLHEP::electron_mass_c2); 304 // 305 G4ThreeVector eDirection, pDirection; 306 // 307 GetAngularDistribution()->SamplePairDirections(aDynamicGamma, 308 eKinEnergy, pKinEnergy, 309 eDirection, pDirection); 310 // create G4DynamicParticle object for the particle1 311 auto aParticle1= new G4DynamicParticle(fTheElectron,eDirection,eKinEnergy); 312 // create G4DynamicParticle object for the particle2 313 auto aParticle2= new G4DynamicParticle(fThePositron,pDirection,pKinEnergy); 314 // Fill output vector 315 fvect->push_back(aParticle1); 316 fvect->push_back(aParticle2); 317 // kill incident photon 318 fParticleChange->SetProposedKineticEnergy(0.); 319 fParticleChange->ProposeTrackStatus(fStopAndKill); 320 } 321 322 // should be called only by the master and at initialisation 323 void G4BetheHeitlerModel::InitialiseElementData() 324 { 325 // create for all elements that are in the detector 326 auto elemTable = G4Element::GetElementTable(); 327 for (auto const & elem : *elemTable) { 328 const G4int Z = elem->GetZasInt(); 329 const G4int iz = std::min(gMaxZet, Z); 330 if (nullptr == gElementData[iz]) { // create it if doesn't exist yet 331 G4double FZLow = 8.*elem->GetIonisation()->GetlogZ3(); 332 G4double FZHigh = FZLow + 8.*elem->GetfCoulomb(); 333 auto elD = new ElementData(); 334 elD->fDeltaMaxLow = G4Exp((42.038 - FZLow )/8.29) - 0.958; 335 elD->fDeltaMaxHigh = G4Exp((42.038 - FZHigh)/8.29) - 0.958; 336 gElementData[iz] = elD; 337 } 338 if (useEPICS2017 && Z < 101) { 339 fXSection->Retrieve(Z); 340 } 341 } 342 343 } 344 345