[968] | 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 | // |
---|
[1192] | 26 | // $Id: G4PenelopePhotoElectricModel.cc,v 1.10 2009/10/23 09:29:24 pandola Exp $ |
---|
| 27 | // GEANT4 tag $Name: emlowen-V09-02-64 $ |
---|
[968] | 28 | // |
---|
| 29 | // Author: Luciano Pandola |
---|
| 30 | // |
---|
| 31 | // History: |
---|
| 32 | // -------- |
---|
[1055] | 33 | // 08 Oct 2008 L Pandola Migration from process to model |
---|
| 34 | // 08 Jan 2009 L Pandola Check shell index to avoid mismatch between |
---|
| 35 | // the Penelope cross section database and the |
---|
| 36 | // G4AtomicTransitionManager database. It suppresses |
---|
| 37 | // a warning from G4AtomicTransitionManager only. |
---|
| 38 | // Results are unchanged. |
---|
| 39 | // 25 Mar 2009 L Pandola Small fix to avoid wrong energy-violation warnings |
---|
| 40 | // 17 Apr 2009 V Ivanchenko Cleanup initialisation and generation of secondaries: |
---|
| 41 | // - apply internal high-energy limit only in constructor |
---|
| 42 | // - do not apply low-energy limit (default is 0) |
---|
| 43 | // - do not apply production threshold on secondaries |
---|
| 44 | // 19 May 2009 L Pandola Explicitely set to zero pointers deleted in |
---|
| 45 | // Initialise(), since they might be checked later on |
---|
[1192] | 46 | // 21 Oct 2009 L Pandola Remove un-necessary fUseAtomicDeexcitation flag - now managed by |
---|
| 47 | // G4VEmModel::DeexcitationFlag() |
---|
[968] | 48 | // |
---|
| 49 | |
---|
| 50 | #include "G4PenelopePhotoElectricModel.hh" |
---|
| 51 | #include "G4ParticleDefinition.hh" |
---|
| 52 | #include "G4MaterialCutsCouple.hh" |
---|
| 53 | #include "G4ProductionCutsTable.hh" |
---|
| 54 | #include "G4DynamicParticle.hh" |
---|
| 55 | #include "G4PhysicsTable.hh" |
---|
| 56 | #include "G4ElementTable.hh" |
---|
| 57 | #include "G4Element.hh" |
---|
| 58 | #include "G4CrossSectionHandler.hh" |
---|
| 59 | #include "G4AtomicTransitionManager.hh" |
---|
| 60 | #include "G4AtomicShell.hh" |
---|
| 61 | #include "G4Gamma.hh" |
---|
| 62 | #include "G4Electron.hh" |
---|
| 63 | |
---|
| 64 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 65 | |
---|
| 66 | |
---|
| 67 | G4PenelopePhotoElectricModel::G4PenelopePhotoElectricModel(const G4ParticleDefinition*, |
---|
| 68 | const G4String& nam) |
---|
| 69 | :G4VEmModel(nam),isInitialised(false),crossSectionHandler(0), |
---|
| 70 | shellCrossSectionHandler(0) |
---|
| 71 | { |
---|
| 72 | fIntrinsicLowEnergyLimit = 100.0*eV; |
---|
| 73 | fIntrinsicHighEnergyLimit = 100.0*GeV; |
---|
[1055] | 74 | // SetLowEnergyLimit(fIntrinsicLowEnergyLimit); |
---|
[968] | 75 | SetHighEnergyLimit(fIntrinsicHighEnergyLimit); |
---|
| 76 | // |
---|
[1192] | 77 | //by default the model will inkove the atomic deexcitation |
---|
| 78 | SetDeexcitationFlag(true); |
---|
| 79 | |
---|
[968] | 80 | verboseLevel= 0; |
---|
| 81 | // Verbosity scale: |
---|
| 82 | // 0 = nothing |
---|
| 83 | // 1 = warning for energy non-conservation |
---|
| 84 | // 2 = details of energy budget |
---|
| 85 | // 3 = calculation of cross sections, file openings, sampling of atoms |
---|
| 86 | // 4 = entering in methods |
---|
| 87 | } |
---|
| 88 | |
---|
| 89 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 90 | |
---|
| 91 | G4PenelopePhotoElectricModel::~G4PenelopePhotoElectricModel() |
---|
| 92 | { |
---|
| 93 | if (crossSectionHandler) delete crossSectionHandler; |
---|
| 94 | if (shellCrossSectionHandler) delete shellCrossSectionHandler; |
---|
| 95 | } |
---|
| 96 | |
---|
| 97 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 98 | |
---|
| 99 | void G4PenelopePhotoElectricModel::Initialise(const G4ParticleDefinition*, |
---|
| 100 | const G4DataVector& ) |
---|
| 101 | { |
---|
| 102 | if (verboseLevel > 3) |
---|
| 103 | G4cout << "Calling G4PenelopePhotoElectricModel::Initialise()" << G4endl; |
---|
| 104 | if (crossSectionHandler) |
---|
| 105 | { |
---|
| 106 | crossSectionHandler->Clear(); |
---|
| 107 | delete crossSectionHandler; |
---|
[1055] | 108 | crossSectionHandler = 0; |
---|
[968] | 109 | } |
---|
| 110 | if (shellCrossSectionHandler) |
---|
| 111 | { |
---|
| 112 | shellCrossSectionHandler->Clear(); |
---|
| 113 | delete shellCrossSectionHandler; |
---|
[1055] | 114 | shellCrossSectionHandler =0; |
---|
[968] | 115 | } |
---|
| 116 | |
---|
| 117 | //Re-initialize cross section handlers |
---|
| 118 | crossSectionHandler = new G4CrossSectionHandler(); |
---|
| 119 | crossSectionHandler->Clear(); |
---|
| 120 | G4String crossSectionFile = "penelope/ph-cs-pen-"; |
---|
| 121 | crossSectionHandler->LoadData(crossSectionFile); |
---|
| 122 | shellCrossSectionHandler = new G4CrossSectionHandler(); |
---|
| 123 | shellCrossSectionHandler->Clear(); |
---|
| 124 | crossSectionFile = "penelope/ph-ss-cs-pen-"; |
---|
| 125 | shellCrossSectionHandler->LoadShellData(crossSectionFile); |
---|
| 126 | //This is used to retrieve cross section values later on |
---|
| 127 | crossSectionHandler->BuildMeanFreePathForMaterials(); |
---|
| 128 | |
---|
| 129 | if (verboseLevel > 2) |
---|
| 130 | G4cout << "Loaded cross section files for PenelopePhotoElectric" << G4endl; |
---|
| 131 | |
---|
[1055] | 132 | if (verboseLevel > 0) { |
---|
| 133 | G4cout << "Penelope Photo-Electric model is initialized " << G4endl |
---|
| 134 | << "Energy range: " |
---|
| 135 | << LowEnergyLimit() / MeV << " MeV - " |
---|
| 136 | << HighEnergyLimit() / GeV << " GeV" |
---|
| 137 | << G4endl; |
---|
| 138 | } |
---|
[968] | 139 | |
---|
| 140 | if(isInitialised) return; |
---|
[1055] | 141 | fParticleChange = GetParticleChangeForGamma(); |
---|
[968] | 142 | isInitialised = true; |
---|
| 143 | } |
---|
| 144 | |
---|
| 145 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 146 | |
---|
| 147 | G4double G4PenelopePhotoElectricModel::ComputeCrossSectionPerAtom( |
---|
| 148 | const G4ParticleDefinition*, |
---|
| 149 | G4double energy, |
---|
| 150 | G4double Z, G4double, |
---|
| 151 | G4double, G4double) |
---|
| 152 | { |
---|
| 153 | // |
---|
| 154 | // Penelope model. Use data-driven approach for cross section estimate (and |
---|
| 155 | // also shell sampling from a given atom). Data are from the Livermore database |
---|
| 156 | // D.E. Cullen et al., Report UCRL-50400 (1989) |
---|
| 157 | // |
---|
| 158 | |
---|
| 159 | if (verboseLevel > 3) |
---|
| 160 | G4cout << "Calling ComputeCrossSectionPerAtom() of G4PenelopePhotoElectricModel" << G4endl; |
---|
| 161 | |
---|
| 162 | G4int iZ = (G4int) Z; |
---|
[1055] | 163 | // if (!crossSectionHandler) // VI: should not be |
---|
| 164 | // { |
---|
| 165 | // G4cout << "G4PenelopePhotoElectricModel::ComputeCrossSectionPerAtom" << G4endl; |
---|
| 166 | // G4cout << "The cross section handler is not correctly initialized" << G4endl; |
---|
| 167 | // G4Exception(); |
---|
| 168 | // } |
---|
[968] | 169 | G4double cs = crossSectionHandler->FindValue(iZ,energy); |
---|
| 170 | |
---|
| 171 | if (verboseLevel > 2) |
---|
| 172 | G4cout << "Photoelectric cross section at " << energy/MeV << " MeV for Z=" << Z << |
---|
| 173 | " = " << cs/barn << " barn" << G4endl; |
---|
| 174 | return cs; |
---|
| 175 | } |
---|
| 176 | |
---|
| 177 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 178 | |
---|
| 179 | void G4PenelopePhotoElectricModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect, |
---|
| 180 | const G4MaterialCutsCouple* couple, |
---|
| 181 | const G4DynamicParticle* aDynamicGamma, |
---|
| 182 | G4double, |
---|
| 183 | G4double) |
---|
| 184 | { |
---|
| 185 | // |
---|
| 186 | // Photoelectric effect, Penelope model |
---|
| 187 | // |
---|
| 188 | // The target atom and the target shell are sampled according to the Livermore |
---|
| 189 | // database |
---|
| 190 | // D.E. Cullen et al., Report UCRL-50400 (1989) |
---|
| 191 | // The angular distribution of the electron in the final state is sampled |
---|
| 192 | // according to the Sauter distribution from |
---|
| 193 | // F. Sauter, Ann. Phys. 11 (1931) 454 |
---|
| 194 | // The energy of the final electron is given by the initial photon energy minus |
---|
| 195 | // the binding energy. Fluorescence de-excitation is subsequently produced |
---|
| 196 | // (to fill the vacancy) according to the general Geant4 G4DeexcitationManager: |
---|
| 197 | // J. Stepanek, Comp. Phys. Comm. 1206 pp 1-1-9 (1997) |
---|
| 198 | |
---|
| 199 | if (verboseLevel > 3) |
---|
| 200 | G4cout << "Calling SamplingSecondaries() of G4PenelopePhotoElectricModel" << G4endl; |
---|
| 201 | |
---|
| 202 | G4double photonEnergy = aDynamicGamma->GetKineticEnergy(); |
---|
| 203 | |
---|
[1055] | 204 | // always kill primary |
---|
| 205 | fParticleChange->ProposeTrackStatus(fStopAndKill); |
---|
| 206 | fParticleChange->SetProposedKineticEnergy(0.); |
---|
| 207 | |
---|
| 208 | if (photonEnergy <= fIntrinsicLowEnergyLimit) |
---|
| 209 | { |
---|
[968] | 210 | fParticleChange->ProposeLocalEnergyDeposit(photonEnergy); |
---|
| 211 | return ; |
---|
[1055] | 212 | } |
---|
[968] | 213 | |
---|
| 214 | G4ParticleMomentum photonDirection = aDynamicGamma->GetMomentumDirection(); |
---|
| 215 | |
---|
| 216 | // Select randomly one element in the current material |
---|
| 217 | if (verboseLevel > 2) |
---|
| 218 | G4cout << "Going to select element in " << couple->GetMaterial()->GetName() << G4endl; |
---|
| 219 | //use crossSectionHandler instead of G4EmElementSelector because in this case |
---|
[1055] | 220 | //the dimension of the table is equal to the dimension of the database |
---|
| 221 | //(less interpolation errors) |
---|
[968] | 222 | G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy); |
---|
| 223 | if (verboseLevel > 2) |
---|
| 224 | G4cout << "Selected Z = " << Z << G4endl; |
---|
| 225 | |
---|
| 226 | // Select the ionised shell in the current atom according to shell cross sections |
---|
| 227 | size_t shellIndex = shellCrossSectionHandler->SelectRandomShell(Z,photonEnergy); |
---|
| 228 | |
---|
| 229 | // Retrieve the corresponding identifier and binding energy of the selected shell |
---|
| 230 | const G4AtomicTransitionManager* transitionManager = G4AtomicTransitionManager::Instance(); |
---|
[1055] | 231 | |
---|
| 232 | //The number of shell cross section possibly reported in the Penelope database |
---|
| 233 | //might be different from the number of shells in the G4AtomicTransitionManager |
---|
| 234 | //(namely, Penelope may contain more shell, especially for very light elements). |
---|
| 235 | //In order to avoid a warning message from the G4AtomicTransitionManager, I |
---|
| 236 | //add this protection. Results are anyway changed, because when G4AtomicTransitionManager |
---|
| 237 | //has a shellID>maxID, it sets the shellID to the last valid shell. |
---|
| 238 | size_t numberOfShells = (size_t) transitionManager->NumberOfShells(Z); |
---|
| 239 | if (shellIndex >= numberOfShells) |
---|
| 240 | shellIndex = numberOfShells-1; |
---|
| 241 | |
---|
[968] | 242 | const G4AtomicShell* shell = transitionManager->Shell(Z,shellIndex); |
---|
| 243 | G4double bindingEnergy = shell->BindingEnergy(); |
---|
| 244 | G4int shellId = shell->ShellId(); |
---|
| 245 | |
---|
| 246 | G4double localEnergyDeposit = 0.0; |
---|
| 247 | |
---|
| 248 | // Primary outcoming electron |
---|
| 249 | G4double eKineticEnergy = photonEnergy - bindingEnergy; |
---|
| 250 | |
---|
| 251 | // There may be cases where the binding energy of the selected shell is > photon energy |
---|
| 252 | // In such cases do not generate secondaries |
---|
| 253 | if (eKineticEnergy > 0.) |
---|
| 254 | { |
---|
| 255 | //Now check if the electron is above cuts: if so, it is created explicitely |
---|
[1055] | 256 | //VI: checking cut here provides inconsistency in testing |
---|
| 257 | // if (eKineticEnergy > cutE) |
---|
| 258 | // The electron is created |
---|
| 259 | // Direction sampled from the Sauter distribution |
---|
| 260 | G4double cosTheta = SampleElectronDirection(eKineticEnergy); |
---|
| 261 | G4double sinTheta = std::sqrt(1-cosTheta*cosTheta); |
---|
| 262 | G4double phi = twopi * G4UniformRand() ; |
---|
| 263 | G4double dirx = sinTheta * std::cos(phi); |
---|
| 264 | G4double diry = sinTheta * std::sin(phi); |
---|
| 265 | G4double dirz = cosTheta ; |
---|
| 266 | G4ThreeVector electronDirection(dirx,diry,dirz); //electron direction |
---|
| 267 | electronDirection.rotateUz(photonDirection); |
---|
| 268 | G4DynamicParticle* electron = new G4DynamicParticle (G4Electron::Electron(), |
---|
| 269 | electronDirection, |
---|
| 270 | eKineticEnergy); |
---|
| 271 | fvect->push_back(electron); |
---|
| 272 | } |
---|
| 273 | // else |
---|
| 274 | // { |
---|
| 275 | // localEnergyDeposit += eKineticEnergy; |
---|
| 276 | // eKineticEnergy = 0; |
---|
| 277 | // } |
---|
| 278 | // } |
---|
[968] | 279 | else |
---|
[1055] | 280 | { |
---|
[968] | 281 | bindingEnergy = photonEnergy; |
---|
[1055] | 282 | } |
---|
[968] | 283 | G4double energyInFluorescence = 0; //testing purposes |
---|
| 284 | |
---|
| 285 | //Now, take care of fluorescence, if required |
---|
[1055] | 286 | if(DeexcitationFlag() && Z > 5) |
---|
[968] | 287 | { |
---|
[1055] | 288 | const G4ProductionCutsTable* theCoupleTable= |
---|
| 289 | G4ProductionCutsTable::GetProductionCutsTable(); |
---|
| 290 | size_t indx = couple->GetIndex(); |
---|
[968] | 291 | G4double cutG = (*(theCoupleTable->GetEnergyCutsVector(0)))[indx]; |
---|
[1055] | 292 | G4double cutE = (*(theCoupleTable->GetEnergyCutsVector(1)))[indx]; |
---|
[968] | 293 | |
---|
| 294 | // Protection to avoid generating photons in the unphysical case of |
---|
| 295 | // shell binding energy > photon energy |
---|
[1055] | 296 | if (bindingEnergy > cutG || bindingEnergy > cutE) |
---|
[968] | 297 | { |
---|
[1055] | 298 | deexcitationManager.SetCutForSecondaryPhotons(cutG); |
---|
| 299 | deexcitationManager.SetCutForAugerElectrons(cutE); |
---|
| 300 | std::vector<G4DynamicParticle*>* photonVector = |
---|
| 301 | deexcitationManager.GenerateParticles(Z,shellId); |
---|
| 302 | //Check for secondaries |
---|
| 303 | if(photonVector) |
---|
[968] | 304 | { |
---|
[1055] | 305 | for (size_t k=0; k< photonVector->size(); k++) |
---|
[968] | 306 | { |
---|
[1055] | 307 | G4DynamicParticle* aPhoton = (*photonVector)[k]; |
---|
| 308 | if (aPhoton) |
---|
[968] | 309 | { |
---|
[1055] | 310 | G4double itsEnergy = aPhoton->GetKineticEnergy(); |
---|
| 311 | if (itsEnergy <= bindingEnergy) |
---|
| 312 | { |
---|
| 313 | if(aPhoton->GetDefinition() == G4Gamma::Gamma()) |
---|
| 314 | energyInFluorescence += itsEnergy; |
---|
| 315 | bindingEnergy -= itsEnergy; |
---|
| 316 | fvect->push_back(aPhoton); |
---|
| 317 | } |
---|
| 318 | else |
---|
| 319 | { |
---|
| 320 | delete aPhoton; |
---|
| 321 | (*photonVector)[k] = 0; |
---|
| 322 | } |
---|
[968] | 323 | } |
---|
| 324 | } |
---|
[1055] | 325 | delete photonVector; |
---|
[968] | 326 | } |
---|
| 327 | } |
---|
| 328 | } |
---|
| 329 | //Residual energy is deposited locally |
---|
| 330 | localEnergyDeposit += bindingEnergy; |
---|
| 331 | |
---|
| 332 | if (localEnergyDeposit < 0) |
---|
| 333 | { |
---|
| 334 | G4cout << "WARNING - " |
---|
| 335 | << "G4PenelopePhotoElectric::PostStepDoIt - Negative energy deposit" |
---|
| 336 | << G4endl; |
---|
| 337 | localEnergyDeposit = 0; |
---|
| 338 | } |
---|
| 339 | |
---|
| 340 | fParticleChange->ProposeLocalEnergyDeposit(localEnergyDeposit); |
---|
| 341 | |
---|
| 342 | if (verboseLevel > 1) |
---|
| 343 | { |
---|
| 344 | G4cout << "-----------------------------------------------------------" << G4endl; |
---|
| 345 | G4cout << "Energy balance from G4PenelopePhotoElectric" << G4endl; |
---|
| 346 | G4cout << "Incoming photon energy: " << photonEnergy/keV << " keV" << G4endl; |
---|
| 347 | G4cout << "-----------------------------------------------------------" << G4endl; |
---|
| 348 | if (eKineticEnergy) |
---|
| 349 | G4cout << "Outgoing electron " << eKineticEnergy/keV << " keV" << G4endl; |
---|
| 350 | G4cout << "Fluorescence: " << energyInFluorescence/keV << " keV" << G4endl; |
---|
| 351 | G4cout << "Local energy deposit " << localEnergyDeposit/keV << " keV" << G4endl; |
---|
| 352 | G4cout << "Total final state: " << (eKineticEnergy+energyInFluorescence+localEnergyDeposit)/keV << |
---|
| 353 | " keV" << G4endl; |
---|
| 354 | G4cout << "-----------------------------------------------------------" << G4endl; |
---|
| 355 | } |
---|
| 356 | if (verboseLevel > 0) |
---|
| 357 | { |
---|
[1055] | 358 | G4double energyDiff = |
---|
| 359 | std::fabs(eKineticEnergy+energyInFluorescence+localEnergyDeposit-photonEnergy); |
---|
[968] | 360 | if (energyDiff > 0.05*keV) |
---|
| 361 | G4cout << "Warning from G4PenelopePhotoElectric: problem with energy conservation: " << |
---|
[1055] | 362 | (eKineticEnergy+energyInFluorescence+localEnergyDeposit)/keV |
---|
| 363 | << " keV (final) vs. " << |
---|
[968] | 364 | photonEnergy/keV << " keV (initial)" << G4endl; |
---|
| 365 | } |
---|
| 366 | } |
---|
| 367 | |
---|
| 368 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 369 | |
---|
| 370 | void G4PenelopePhotoElectricModel::ActivateAuger(G4bool augerbool) |
---|
| 371 | { |
---|
[1192] | 372 | if (!DeexcitationFlag() && augerbool) |
---|
[968] | 373 | { |
---|
| 374 | G4cout << "WARNING - G4PenelopePhotoElectricModel" << G4endl; |
---|
| 375 | G4cout << "The use of the Atomic Deexcitation Manager is set to false " << G4endl; |
---|
| 376 | G4cout << "Therefore, Auger electrons will be not generated anyway" << G4endl; |
---|
| 377 | } |
---|
| 378 | deexcitationManager.ActivateAugerElectronProduction(augerbool); |
---|
| 379 | if (verboseLevel > 1) |
---|
| 380 | G4cout << "Auger production set to " << augerbool << G4endl; |
---|
| 381 | } |
---|
| 382 | |
---|
| 383 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 384 | |
---|
| 385 | G4double G4PenelopePhotoElectricModel::SampleElectronDirection(G4double energy) |
---|
| 386 | { |
---|
| 387 | G4double costheta = 1.0; |
---|
| 388 | if (energy>1*GeV) return costheta; |
---|
| 389 | |
---|
| 390 | //1) initialize energy-dependent variables |
---|
| 391 | // Variable naming according to Eq. (2.24) of Penelope Manual |
---|
| 392 | // (pag. 44) |
---|
| 393 | G4double gamma = 1.0 + energy/electron_mass_c2; |
---|
| 394 | G4double gamma2 = gamma*gamma; |
---|
| 395 | G4double beta = std::sqrt((gamma2-1.0)/gamma2); |
---|
| 396 | |
---|
| 397 | // ac corresponds to "A" of Eq. (2.31) |
---|
| 398 | // |
---|
| 399 | G4double ac = (1.0/beta) - 1.0; |
---|
| 400 | G4double a1 = 0.5*beta*gamma*(gamma-1.0)*(gamma-2.0); |
---|
| 401 | G4double a2 = ac + 2.0; |
---|
| 402 | G4double gtmax = 2.0*(a1 + 1.0/ac); |
---|
| 403 | |
---|
| 404 | G4double tsam = 0; |
---|
| 405 | G4double gtr = 0; |
---|
| 406 | |
---|
| 407 | //2) sampling. Eq. (2.31) of Penelope Manual |
---|
| 408 | // tsam = 1-std::cos(theta) |
---|
| 409 | // gtr = rejection function according to Eq. (2.28) |
---|
| 410 | do{ |
---|
| 411 | G4double rand = G4UniformRand(); |
---|
| 412 | tsam = 2.0*ac * (2.0*rand + a2*std::sqrt(rand)) / (a2*a2 - 4.0*rand); |
---|
| 413 | gtr = (2.0 - tsam) * (a1 + 1.0/(ac+tsam)); |
---|
| 414 | }while(G4UniformRand()*gtmax > gtr); |
---|
| 415 | costheta = 1.0-tsam; |
---|
| 416 | return costheta; |
---|
| 417 | } |
---|
| 418 | |
---|