| 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 | // $Id: G4LivermoreNuclearGammaConversionModel.cc,v 1.1 2010/11/10 17:09:16 flongo Exp $
|
|---|
| 27 | // GEANT4 tag $Name: geant4-09-04-ref-00 $
|
|---|
| 28 | //
|
|---|
| 29 | //
|
|---|
| 30 | // Author: Sebastien Inserti
|
|---|
| 31 | // 30 October 2008
|
|---|
| 32 | //
|
|---|
| 33 | // History:
|
|---|
| 34 | // --------
|
|---|
| 35 | // 12 Apr 2009 V Ivanchenko Cleanup initialisation and generation of secondaries:
|
|---|
| 36 | // - apply internal high-energy limit only in constructor
|
|---|
| 37 | // - do not apply low-energy limit (default is 0)
|
|---|
| 38 | // - use CLHEP electron mass for low-enegry limit
|
|---|
| 39 | // - remove MeanFreePath method and table
|
|---|
| 40 |
|
|---|
| 41 |
|
|---|
| 42 | #include "G4LivermoreNuclearGammaConversionModel.hh"
|
|---|
| 43 |
|
|---|
| 44 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 45 |
|
|---|
| 46 | using namespace std;
|
|---|
| 47 |
|
|---|
| 48 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 49 |
|
|---|
| 50 | G4LivermoreNuclearGammaConversionModel::G4LivermoreNuclearGammaConversionModel(const G4ParticleDefinition*,
|
|---|
| 51 | const G4String& nam)
|
|---|
| 52 | :G4VEmModel(nam),smallEnergy(2.*MeV),isInitialised(false),
|
|---|
| 53 | crossSectionHandler(0),meanFreePathTable(0)
|
|---|
| 54 | {
|
|---|
| 55 | lowEnergyLimit = 2.0*electron_mass_c2;
|
|---|
| 56 | highEnergyLimit = 100 * GeV;
|
|---|
| 57 | SetHighEnergyLimit(highEnergyLimit);
|
|---|
| 58 |
|
|---|
| 59 | verboseLevel= 0;
|
|---|
| 60 | // Verbosity scale:
|
|---|
| 61 | // 0 = nothing
|
|---|
| 62 | // 1 = warning for energy non-conservation
|
|---|
| 63 | // 2 = details of energy budget
|
|---|
| 64 | // 3 = calculation of cross sections, file openings, sampling of atoms
|
|---|
| 65 | // 4 = entering in methods
|
|---|
| 66 |
|
|---|
| 67 | if(verboseLevel > 0) {
|
|---|
| 68 | G4cout << "Livermore Nuclear Gamma conversion is constructed " << G4endl
|
|---|
| 69 | << "Energy range: "
|
|---|
| 70 | << lowEnergyLimit / MeV << " MeV - "
|
|---|
| 71 | << highEnergyLimit / GeV << " GeV"
|
|---|
| 72 | << G4endl;
|
|---|
| 73 | }
|
|---|
| 74 | }
|
|---|
| 75 |
|
|---|
| 76 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 77 |
|
|---|
| 78 | G4LivermoreNuclearGammaConversionModel::~G4LivermoreNuclearGammaConversionModel()
|
|---|
| 79 | {
|
|---|
| 80 | if (crossSectionHandler) delete crossSectionHandler;
|
|---|
| 81 | }
|
|---|
| 82 |
|
|---|
| 83 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 84 |
|
|---|
| 85 | void
|
|---|
| 86 | G4LivermoreNuclearGammaConversionModel::Initialise(const G4ParticleDefinition*,
|
|---|
| 87 | const G4DataVector&)
|
|---|
| 88 | {
|
|---|
| 89 | if (verboseLevel > 3)
|
|---|
| 90 | G4cout << "Calling G4LivermoreNuclearGammaConversionModel::Initialise()" << G4endl;
|
|---|
| 91 |
|
|---|
| 92 | if (crossSectionHandler)
|
|---|
| 93 | {
|
|---|
| 94 | crossSectionHandler->Clear();
|
|---|
| 95 | delete crossSectionHandler;
|
|---|
| 96 | }
|
|---|
| 97 |
|
|---|
| 98 | // Read data tables for all materials
|
|---|
| 99 |
|
|---|
| 100 | crossSectionHandler = new G4CrossSectionHandler();
|
|---|
| 101 | crossSectionHandler->Initialise(0,lowEnergyLimit,100.*GeV,400);
|
|---|
| 102 | G4String crossSectionFile = "pairdata/pp-pair-cs-"; // here only pair in nuclear field cs should be used
|
|---|
| 103 | crossSectionHandler->LoadData(crossSectionFile);
|
|---|
| 104 |
|
|---|
| 105 | //
|
|---|
| 106 |
|
|---|
| 107 | if (verboseLevel > 0) {
|
|---|
| 108 | G4cout << "Loaded cross section files for Livermore GammaConversion" << G4endl;
|
|---|
| 109 | G4cout << "To obtain the total cross section this should be used only " << G4endl
|
|---|
| 110 | << "in connection with G4ElectronGammaConversion " << G4endl;
|
|---|
| 111 | }
|
|---|
| 112 |
|
|---|
| 113 | if (verboseLevel > 0) {
|
|---|
| 114 | G4cout << "Livermore Nuclear Gamma Conversion model is initialized " << G4endl
|
|---|
| 115 | << "Energy range: "
|
|---|
| 116 | << LowEnergyLimit() / MeV << " MeV - "
|
|---|
| 117 | << HighEnergyLimit() / GeV << " GeV"
|
|---|
| 118 | << G4endl;
|
|---|
| 119 | }
|
|---|
| 120 |
|
|---|
| 121 | if(isInitialised) return;
|
|---|
| 122 | fParticleChange = GetParticleChangeForGamma();
|
|---|
| 123 | isInitialised = true;
|
|---|
| 124 | }
|
|---|
| 125 |
|
|---|
| 126 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 127 |
|
|---|
| 128 | G4double
|
|---|
| 129 | G4LivermoreNuclearGammaConversionModel::ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
|
|---|
| 130 | G4double GammaEnergy,
|
|---|
| 131 | G4double Z, G4double,
|
|---|
| 132 | G4double, G4double)
|
|---|
| 133 | {
|
|---|
| 134 | if (verboseLevel > 3) {
|
|---|
| 135 | G4cout << "Calling ComputeCrossSectionPerAtom() of G4LivermoreNuclearGammaConversionModel"
|
|---|
| 136 | << G4endl;
|
|---|
| 137 | }
|
|---|
| 138 | if (GammaEnergy < lowEnergyLimit || GammaEnergy > highEnergyLimit) return 0;
|
|---|
| 139 |
|
|---|
| 140 | G4double cs = crossSectionHandler->FindValue(G4int(Z), GammaEnergy);
|
|---|
| 141 | return cs;
|
|---|
| 142 | }
|
|---|
| 143 |
|
|---|
| 144 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 145 |
|
|---|
| 146 | void G4LivermoreNuclearGammaConversionModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
|
|---|
| 147 | const G4MaterialCutsCouple* couple,
|
|---|
| 148 | const G4DynamicParticle* aDynamicGamma,
|
|---|
| 149 | G4double,
|
|---|
| 150 | G4double)
|
|---|
| 151 | {
|
|---|
| 152 |
|
|---|
| 153 | // The energies of the e+ e- secondaries are sampled using the Bethe - Heitler
|
|---|
| 154 | // cross sections with Coulomb correction. A modified version of the random
|
|---|
| 155 | // number techniques of Butcher & Messel is used (Nuc Phys 20(1960),15).
|
|---|
| 156 |
|
|---|
| 157 | // Note 1 : Effects due to the breakdown of the Born approximation at low
|
|---|
| 158 | // energy are ignored.
|
|---|
| 159 | // Note 2 : The differential cross section implicitly takes account of
|
|---|
| 160 | // pair creation in both nuclear and atomic electron fields. However triplet
|
|---|
| 161 | // prodution is not generated.
|
|---|
| 162 |
|
|---|
| 163 | if (verboseLevel > 3)
|
|---|
| 164 | G4cout << "Calling SampleSecondaries() of G4LivermoreNuclearGammaConversionModel" << G4endl;
|
|---|
| 165 |
|
|---|
| 166 | G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
|
|---|
| 167 | G4ParticleMomentum photonDirection = aDynamicGamma->GetMomentumDirection();
|
|---|
| 168 |
|
|---|
| 169 | G4double epsilon ;
|
|---|
| 170 | G4double epsilon0 = electron_mass_c2 / photonEnergy ;
|
|---|
| 171 |
|
|---|
| 172 | // Do it fast if photon energy < 2. MeV
|
|---|
| 173 | if (photonEnergy < smallEnergy )
|
|---|
| 174 | {
|
|---|
| 175 | epsilon = epsilon0 + (0.5 - epsilon0) * G4UniformRand();
|
|---|
| 176 | }
|
|---|
| 177 | else
|
|---|
| 178 | {
|
|---|
| 179 | // Select randomly one element in the current material
|
|---|
| 180 | //const G4Element* element = crossSectionHandler->SelectRandomElement(couple,photonEnergy);
|
|---|
| 181 | const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
|
|---|
| 182 | const G4Element* element = SelectRandomAtom(couple,particle,photonEnergy);
|
|---|
| 183 |
|
|---|
| 184 | if (element == 0)
|
|---|
| 185 | {
|
|---|
| 186 | G4cout << "G4LivermoreNuclearGammaConversionModel::SampleSecondaries - element = 0"
|
|---|
| 187 | << G4endl;
|
|---|
| 188 | return;
|
|---|
| 189 | }
|
|---|
| 190 | G4IonisParamElm* ionisation = element->GetIonisation();
|
|---|
| 191 | if (ionisation == 0)
|
|---|
| 192 | {
|
|---|
| 193 | G4cout << "G4LivermoreNuclearGammaConversionModel::SampleSecondaries - ionisation = 0"
|
|---|
| 194 | << G4endl;
|
|---|
| 195 | return;
|
|---|
| 196 | }
|
|---|
| 197 |
|
|---|
| 198 | // Extract Coulomb factor for this Element
|
|---|
| 199 | G4double fZ = 8. * (ionisation->GetlogZ3());
|
|---|
| 200 | if (photonEnergy > 50. * MeV) fZ += 8. * (element->GetfCoulomb());
|
|---|
| 201 |
|
|---|
| 202 | // Limits of the screening variable
|
|---|
| 203 | G4double screenFactor = 136. * epsilon0 / (element->GetIonisation()->GetZ3()) ;
|
|---|
| 204 | G4double screenMax = std::exp ((42.24 - fZ)/8.368) - 0.952 ;
|
|---|
| 205 | G4double screenMin = std::min(4.*screenFactor,screenMax) ;
|
|---|
| 206 |
|
|---|
| 207 | // Limits of the energy sampling
|
|---|
| 208 | G4double epsilon1 = 0.5 - 0.5 * std::sqrt(1. - screenMin / screenMax) ;
|
|---|
| 209 | G4double epsilonMin = std::max(epsilon0,epsilon1);
|
|---|
| 210 | G4double epsilonRange = 0.5 - epsilonMin ;
|
|---|
| 211 |
|
|---|
| 212 | // Sample the energy rate of the created electron (or positron)
|
|---|
| 213 | G4double screen;
|
|---|
| 214 | G4double gReject ;
|
|---|
| 215 |
|
|---|
| 216 | G4double f10 = ScreenFunction1(screenMin) - fZ;
|
|---|
| 217 | G4double f20 = ScreenFunction2(screenMin) - fZ;
|
|---|
| 218 | G4double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
|
|---|
| 219 | G4double normF2 = std::max(1.5 * f20,0.);
|
|---|
| 220 |
|
|---|
| 221 | do {
|
|---|
| 222 | if (normF1 / (normF1 + normF2) > G4UniformRand() )
|
|---|
| 223 | {
|
|---|
| 224 | epsilon = 0.5 - epsilonRange * std::pow(G4UniformRand(), 0.3333) ;
|
|---|
| 225 | screen = screenFactor / (epsilon * (1. - epsilon));
|
|---|
| 226 | gReject = (ScreenFunction1(screen) - fZ) / f10 ;
|
|---|
| 227 | }
|
|---|
| 228 | else
|
|---|
| 229 | {
|
|---|
| 230 | epsilon = epsilonMin + epsilonRange * G4UniformRand();
|
|---|
| 231 | screen = screenFactor / (epsilon * (1 - epsilon));
|
|---|
| 232 | gReject = (ScreenFunction2(screen) - fZ) / f20 ;
|
|---|
| 233 | }
|
|---|
| 234 | } while ( gReject < G4UniformRand() );
|
|---|
| 235 |
|
|---|
| 236 | } // End of epsilon sampling
|
|---|
| 237 |
|
|---|
| 238 | // Fix charges randomly
|
|---|
| 239 |
|
|---|
| 240 | G4double electronTotEnergy;
|
|---|
| 241 | G4double positronTotEnergy;
|
|---|
| 242 |
|
|---|
| 243 | if (CLHEP::RandBit::shootBit())
|
|---|
| 244 | {
|
|---|
| 245 | electronTotEnergy = (1. - epsilon) * photonEnergy;
|
|---|
| 246 | positronTotEnergy = epsilon * photonEnergy;
|
|---|
| 247 | }
|
|---|
| 248 | else
|
|---|
| 249 | {
|
|---|
| 250 | positronTotEnergy = (1. - epsilon) * photonEnergy;
|
|---|
| 251 | electronTotEnergy = epsilon * photonEnergy;
|
|---|
| 252 | }
|
|---|
| 253 |
|
|---|
| 254 | // Scattered electron (positron) angles. ( Z - axis along the parent photon)
|
|---|
| 255 | // Universal distribution suggested by L. Urban (Geant3 manual (1993) Phys211),
|
|---|
| 256 | // derived from Tsai distribution (Rev. Mod. Phys. 49, 421 (1977)
|
|---|
| 257 |
|
|---|
| 258 | G4double u;
|
|---|
| 259 | const G4double a1 = 0.625;
|
|---|
| 260 | G4double a2 = 3. * a1;
|
|---|
| 261 | // G4double d = 27. ;
|
|---|
| 262 |
|
|---|
| 263 | // if (9. / (9. + d) > G4UniformRand())
|
|---|
| 264 | if (0.25 > G4UniformRand())
|
|---|
| 265 | {
|
|---|
| 266 | u = - std::log(G4UniformRand() * G4UniformRand()) / a1 ;
|
|---|
| 267 | }
|
|---|
| 268 | else
|
|---|
| 269 | {
|
|---|
| 270 | u = - std::log(G4UniformRand() * G4UniformRand()) / a2 ;
|
|---|
| 271 | }
|
|---|
| 272 |
|
|---|
| 273 | G4double thetaEle = u*electron_mass_c2/electronTotEnergy;
|
|---|
| 274 | G4double thetaPos = u*electron_mass_c2/positronTotEnergy;
|
|---|
| 275 | G4double phi = twopi * G4UniformRand();
|
|---|
| 276 |
|
|---|
| 277 | G4double dxEle= std::sin(thetaEle)*std::cos(phi),dyEle= std::sin(thetaEle)*std::sin(phi),dzEle=std::cos(thetaEle);
|
|---|
| 278 | G4double dxPos=-std::sin(thetaPos)*std::cos(phi),dyPos=-std::sin(thetaPos)*std::sin(phi),dzPos=std::cos(thetaPos);
|
|---|
| 279 |
|
|---|
| 280 |
|
|---|
| 281 | // Kinematics of the created pair:
|
|---|
| 282 | // the electron and positron are assumed to have a symetric angular
|
|---|
| 283 | // distribution with respect to the Z axis along the parent photon
|
|---|
| 284 |
|
|---|
| 285 | G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
|
|---|
| 286 |
|
|---|
| 287 | // SI - The range test has been removed wrt original G4LowEnergyGammaconversion class
|
|---|
| 288 |
|
|---|
| 289 | G4ThreeVector electronDirection (dxEle, dyEle, dzEle);
|
|---|
| 290 | electronDirection.rotateUz(photonDirection);
|
|---|
| 291 |
|
|---|
| 292 | G4DynamicParticle* particle1 = new G4DynamicParticle (G4Electron::Electron(),
|
|---|
| 293 | electronDirection,
|
|---|
| 294 | electronKineEnergy);
|
|---|
| 295 |
|
|---|
| 296 | // The e+ is always created (even with kinetic energy = 0) for further annihilation
|
|---|
| 297 | G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
|
|---|
| 298 |
|
|---|
| 299 | // SI - The range test has been removed wrt original G4LowEnergyGammaconversion class
|
|---|
| 300 |
|
|---|
| 301 | G4ThreeVector positronDirection (dxPos, dyPos, dzPos);
|
|---|
| 302 | positronDirection.rotateUz(photonDirection);
|
|---|
| 303 |
|
|---|
| 304 | // Create G4DynamicParticle object for the particle2
|
|---|
| 305 | G4DynamicParticle* particle2 = new G4DynamicParticle(G4Positron::Positron(),
|
|---|
| 306 | positronDirection, positronKineEnergy);
|
|---|
| 307 | // Fill output vector
|
|---|
| 308 |
|
|---|
| 309 | fvect->push_back(particle1);
|
|---|
| 310 | fvect->push_back(particle2);
|
|---|
| 311 |
|
|---|
| 312 | // kill incident photon
|
|---|
| 313 | fParticleChange->SetProposedKineticEnergy(0.);
|
|---|
| 314 | fParticleChange->ProposeTrackStatus(fStopAndKill);
|
|---|
| 315 |
|
|---|
| 316 | }
|
|---|
| 317 |
|
|---|
| 318 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 319 |
|
|---|
| 320 | G4double G4LivermoreNuclearGammaConversionModel::ScreenFunction1(G4double screenVariable)
|
|---|
| 321 | {
|
|---|
| 322 | // Compute the value of the screening function 3*phi1 - phi2
|
|---|
| 323 |
|
|---|
| 324 | G4double value;
|
|---|
| 325 |
|
|---|
| 326 | if (screenVariable > 1.)
|
|---|
| 327 | value = 42.24 - 8.368 * std::log(screenVariable + 0.952);
|
|---|
| 328 | else
|
|---|
| 329 | value = 42.392 - screenVariable * (7.796 - 1.961 * screenVariable);
|
|---|
| 330 |
|
|---|
| 331 | return value;
|
|---|
| 332 | }
|
|---|
| 333 |
|
|---|
| 334 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 335 |
|
|---|
| 336 | G4double G4LivermoreNuclearGammaConversionModel::ScreenFunction2(G4double screenVariable)
|
|---|
| 337 | {
|
|---|
| 338 | // Compute the value of the screening function 1.5*phi1 - 0.5*phi2
|
|---|
| 339 |
|
|---|
| 340 | G4double value;
|
|---|
| 341 |
|
|---|
| 342 | if (screenVariable > 1.)
|
|---|
| 343 | value = 42.24 - 8.368 * std::log(screenVariable + 0.952);
|
|---|
| 344 | else
|
|---|
| 345 | value = 41.405 - screenVariable * (5.828 - 0.8945 * screenVariable);
|
|---|
| 346 |
|
|---|
| 347 | return value;
|
|---|
| 348 | }
|
|---|
| 349 |
|
|---|