[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 | // |
---|
| 26 | // $Id: G4LivermoreComptonModel.cc,v 1.2 2009/01/21 10:58:13 sincerti Exp $ |
---|
| 27 | // GEANT4 tag $Name: geant4-09-02-ref-02 $ |
---|
| 28 | // |
---|
| 29 | |
---|
| 30 | #include "G4LivermoreComptonModel.hh" |
---|
| 31 | |
---|
| 32 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 33 | |
---|
| 34 | using namespace std; |
---|
| 35 | |
---|
| 36 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 37 | |
---|
| 38 | G4LivermoreComptonModel::G4LivermoreComptonModel(const G4ParticleDefinition*, |
---|
| 39 | const G4String& nam) |
---|
| 40 | :G4VEmModel(nam),isInitialised(false),meanFreePathTable(0),scatterFunctionData(0),crossSectionHandler(0) |
---|
| 41 | { |
---|
| 42 | lowEnergyLimit = 250 * eV; // SI - Could be 10 eV ? |
---|
| 43 | highEnergyLimit = 100 * GeV; |
---|
| 44 | SetLowEnergyLimit(lowEnergyLimit); |
---|
| 45 | SetHighEnergyLimit(highEnergyLimit); |
---|
| 46 | |
---|
| 47 | verboseLevel=0 ; |
---|
| 48 | // Verbosity scale: |
---|
| 49 | // 0 = nothing |
---|
| 50 | // 1 = warning for energy non-conservation |
---|
| 51 | // 2 = details of energy budget |
---|
| 52 | // 3 = calculation of cross sections, file openings, sampling of atoms |
---|
| 53 | // 4 = entering in methods |
---|
| 54 | |
---|
| 55 | G4cout << "Livermore Compton model is constructed " << G4endl |
---|
| 56 | << "Energy range: " |
---|
| 57 | << lowEnergyLimit / keV << " keV - " |
---|
| 58 | << highEnergyLimit / GeV << " GeV" |
---|
| 59 | << G4endl; |
---|
| 60 | |
---|
| 61 | } |
---|
| 62 | |
---|
| 63 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 64 | |
---|
| 65 | G4LivermoreComptonModel::~G4LivermoreComptonModel() |
---|
| 66 | { |
---|
| 67 | |
---|
| 68 | if (meanFreePathTable) delete meanFreePathTable; |
---|
| 69 | if (crossSectionHandler) delete crossSectionHandler; |
---|
| 70 | if (scatterFunctionData) delete scatterFunctionData; |
---|
| 71 | |
---|
| 72 | } |
---|
| 73 | |
---|
| 74 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 75 | |
---|
| 76 | void G4LivermoreComptonModel::Initialise(const G4ParticleDefinition* particle, |
---|
| 77 | const G4DataVector& cuts) |
---|
| 78 | { |
---|
| 79 | |
---|
| 80 | if (verboseLevel > 3) |
---|
| 81 | G4cout << "Calling G4LivermoreComptonModel::Initialise()" << G4endl; |
---|
| 82 | |
---|
| 83 | if (crossSectionHandler) |
---|
| 84 | { |
---|
| 85 | crossSectionHandler->Clear(); |
---|
| 86 | delete crossSectionHandler; |
---|
| 87 | } |
---|
| 88 | |
---|
| 89 | // Energy limits |
---|
| 90 | |
---|
| 91 | if (LowEnergyLimit() < lowEnergyLimit) |
---|
| 92 | { |
---|
| 93 | G4cout << "G4LivermoreComptonModel: low energy limit increased from " << |
---|
| 94 | LowEnergyLimit()/eV << " eV to " << lowEnergyLimit << " eV" << G4endl; |
---|
| 95 | SetLowEnergyLimit(lowEnergyLimit); |
---|
| 96 | } |
---|
| 97 | |
---|
| 98 | if (HighEnergyLimit() > highEnergyLimit) |
---|
| 99 | { |
---|
| 100 | G4cout << "G4LivermoreComptonModel: high energy limit decreased from " << |
---|
| 101 | HighEnergyLimit()/GeV << " GeV to " << highEnergyLimit << " GeV" << G4endl; |
---|
| 102 | SetHighEnergyLimit(highEnergyLimit); |
---|
| 103 | } |
---|
| 104 | |
---|
| 105 | // Reading of data files - all materials are read |
---|
| 106 | |
---|
| 107 | crossSectionHandler = new G4CrossSectionHandler; |
---|
| 108 | crossSectionHandler->Clear(); |
---|
| 109 | G4String crossSectionFile = "comp/ce-cs-"; |
---|
| 110 | crossSectionHandler->LoadData(crossSectionFile); |
---|
| 111 | |
---|
| 112 | G4VDataSetAlgorithm* scatterInterpolation = new G4LogLogInterpolation; |
---|
| 113 | G4String scatterFile = "comp/ce-sf-"; |
---|
| 114 | scatterFunctionData = new G4CompositeEMDataSet(scatterInterpolation, 1., 1.); |
---|
| 115 | scatterFunctionData->LoadData(scatterFile); |
---|
| 116 | |
---|
| 117 | // For Doppler broadening |
---|
| 118 | shellData.SetOccupancyData(); |
---|
| 119 | G4String file = "/doppler/shell-doppler"; |
---|
| 120 | shellData.LoadData(file); |
---|
| 121 | |
---|
| 122 | meanFreePathTable = 0; |
---|
| 123 | meanFreePathTable = crossSectionHandler->BuildMeanFreePathForMaterials(); |
---|
| 124 | |
---|
| 125 | if (verboseLevel > 2) |
---|
| 126 | G4cout << "Loaded cross section files for Livermore Compton model" << G4endl; |
---|
| 127 | |
---|
| 128 | InitialiseElementSelectors(particle,cuts); |
---|
| 129 | |
---|
| 130 | G4cout << "Livermore Compton model is initialized " << G4endl |
---|
| 131 | << "Energy range: " |
---|
| 132 | << LowEnergyLimit() / keV << " keV - " |
---|
| 133 | << HighEnergyLimit() / GeV << " GeV" |
---|
| 134 | << G4endl; |
---|
| 135 | |
---|
| 136 | // |
---|
| 137 | |
---|
| 138 | if(isInitialised) return; |
---|
| 139 | |
---|
| 140 | if(pParticleChange) |
---|
| 141 | fParticleChange = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange); |
---|
| 142 | else |
---|
| 143 | fParticleChange = new G4ParticleChangeForGamma(); |
---|
| 144 | |
---|
| 145 | isInitialised = true; |
---|
| 146 | |
---|
| 147 | } |
---|
| 148 | |
---|
| 149 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 150 | |
---|
| 151 | G4double G4LivermoreComptonModel::ComputeCrossSectionPerAtom( |
---|
| 152 | const G4ParticleDefinition*, |
---|
| 153 | G4double GammaEnergy, |
---|
| 154 | G4double Z, G4double, |
---|
| 155 | G4double, G4double) |
---|
| 156 | { |
---|
| 157 | |
---|
| 158 | if (verboseLevel > 3) |
---|
| 159 | G4cout << "Calling ComputeCrossSectionPerAtom() of G4LivermoreComptonModel" << G4endl; |
---|
| 160 | |
---|
| 161 | G4double cs = crossSectionHandler->FindValue(G4int(Z), GammaEnergy); |
---|
| 162 | |
---|
| 163 | return cs; |
---|
| 164 | |
---|
| 165 | } |
---|
| 166 | |
---|
| 167 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 168 | |
---|
| 169 | void G4LivermoreComptonModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect, |
---|
| 170 | const G4MaterialCutsCouple* couple, |
---|
| 171 | const G4DynamicParticle* aDynamicGamma, |
---|
| 172 | G4double, |
---|
| 173 | G4double) |
---|
| 174 | { |
---|
| 175 | |
---|
| 176 | // The scattered gamma energy is sampled according to Klein - Nishina formula. |
---|
| 177 | // then accepted or rejected depending on the Scattering Function multiplied |
---|
| 178 | // by factor from Klein - Nishina formula. |
---|
| 179 | // Expression of the angular distribution as Klein Nishina |
---|
| 180 | // angular and energy distribution and Scattering fuctions is taken from |
---|
| 181 | // D. E. Cullen "A simple model of photon transport" Nucl. Instr. Meth. |
---|
| 182 | // Phys. Res. B 101 (1995). Method of sampling with form factors is different |
---|
| 183 | // data are interpolated while in the article they are fitted. |
---|
| 184 | // Reference to the article is from J. Stepanek New Photon, Positron |
---|
| 185 | // and Electron Interaction Data for GEANT in Energy Range from 1 eV to 10 |
---|
| 186 | // TeV (draft). |
---|
| 187 | // The random number techniques of Butcher & Messel are used |
---|
| 188 | // (Nucl Phys 20(1960),15). |
---|
| 189 | |
---|
| 190 | if (verboseLevel > 3) |
---|
| 191 | G4cout << "Calling SampleSecondaries() of G4LivermoreComptonModel" << G4endl; |
---|
| 192 | |
---|
| 193 | G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy(); |
---|
| 194 | |
---|
| 195 | if (photonEnergy0 <= lowEnergyLimit) |
---|
| 196 | { |
---|
| 197 | fParticleChange->ProposeTrackStatus(fStopAndKill); |
---|
| 198 | fParticleChange->SetProposedKineticEnergy(0.); |
---|
| 199 | fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0); |
---|
| 200 | return ; |
---|
| 201 | } |
---|
| 202 | |
---|
| 203 | G4double e0m = photonEnergy0 / electron_mass_c2 ; |
---|
| 204 | G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection(); |
---|
| 205 | |
---|
| 206 | // Select randomly one element in the current material |
---|
| 207 | G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy0); |
---|
| 208 | |
---|
| 209 | G4double epsilon0 = 1. / (1. + 2. * e0m); |
---|
| 210 | G4double epsilon0Sq = epsilon0 * epsilon0; |
---|
| 211 | G4double alpha1 = -std::log(epsilon0); |
---|
| 212 | G4double alpha2 = 0.5 * (1. - epsilon0Sq); |
---|
| 213 | |
---|
| 214 | G4double wlPhoton = h_Planck*c_light/photonEnergy0; |
---|
| 215 | |
---|
| 216 | // Sample the energy of the scattered photon |
---|
| 217 | G4double epsilon; |
---|
| 218 | G4double epsilonSq; |
---|
| 219 | G4double oneCosT; |
---|
| 220 | G4double sinT2; |
---|
| 221 | G4double gReject; |
---|
| 222 | |
---|
| 223 | do |
---|
| 224 | { |
---|
| 225 | if ( alpha1/(alpha1+alpha2) > G4UniformRand()) |
---|
| 226 | { |
---|
| 227 | epsilon = std::exp(-alpha1 * G4UniformRand()); // std::pow(epsilon0,G4UniformRand()) |
---|
| 228 | epsilonSq = epsilon * epsilon; |
---|
| 229 | } |
---|
| 230 | else |
---|
| 231 | { |
---|
| 232 | epsilonSq = epsilon0Sq + (1. - epsilon0Sq) * G4UniformRand(); |
---|
| 233 | epsilon = std::sqrt(epsilonSq); |
---|
| 234 | } |
---|
| 235 | |
---|
| 236 | oneCosT = (1. - epsilon) / ( epsilon * e0m); |
---|
| 237 | sinT2 = oneCosT * (2. - oneCosT); |
---|
| 238 | G4double x = std::sqrt(oneCosT/2.) / (wlPhoton/cm); |
---|
| 239 | G4double scatteringFunction = scatterFunctionData->FindValue(x,Z-1); |
---|
| 240 | gReject = (1. - epsilon * sinT2 / (1. + epsilonSq)) * scatteringFunction; |
---|
| 241 | |
---|
| 242 | } while(gReject < G4UniformRand()*Z); |
---|
| 243 | |
---|
| 244 | G4double cosTheta = 1. - oneCosT; |
---|
| 245 | G4double sinTheta = std::sqrt (sinT2); |
---|
| 246 | G4double phi = twopi * G4UniformRand() ; |
---|
| 247 | G4double dirx = sinTheta * std::cos(phi); |
---|
| 248 | G4double diry = sinTheta * std::sin(phi); |
---|
| 249 | G4double dirz = cosTheta ; |
---|
| 250 | |
---|
| 251 | // Doppler broadening - Method based on: |
---|
| 252 | // Y. Namito, S. Ban and H. Hirayama, |
---|
| 253 | // "Implementation of the Doppler Broadening of a Compton-Scattered Photon Into the EGS4 Code" |
---|
| 254 | // NIM A 349, pp. 489-494, 1994 |
---|
| 255 | |
---|
| 256 | // Maximum number of sampling iterations |
---|
| 257 | G4int maxDopplerIterations = 1000; |
---|
| 258 | G4double bindingE = 0.; |
---|
| 259 | G4double photonEoriginal = epsilon * photonEnergy0; |
---|
| 260 | G4double photonE = -1.; |
---|
| 261 | G4int iteration = 0; |
---|
| 262 | G4double eMax = photonEnergy0; |
---|
| 263 | do |
---|
| 264 | { |
---|
| 265 | iteration++; |
---|
| 266 | // Select shell based on shell occupancy |
---|
| 267 | G4int shell = shellData.SelectRandomShell(Z); |
---|
| 268 | bindingE = shellData.BindingEnergy(Z,shell); |
---|
| 269 | |
---|
| 270 | eMax = photonEnergy0 - bindingE; |
---|
| 271 | |
---|
| 272 | // Randomly sample bound electron momentum (memento: the data set is in Atomic Units) |
---|
| 273 | G4double pSample = profileData.RandomSelectMomentum(Z,shell); |
---|
| 274 | // Rescale from atomic units |
---|
| 275 | G4double pDoppler = pSample * fine_structure_const; |
---|
| 276 | G4double pDoppler2 = pDoppler * pDoppler; |
---|
| 277 | G4double var2 = 1. + oneCosT * e0m; |
---|
| 278 | G4double var3 = var2*var2 - pDoppler2; |
---|
| 279 | G4double var4 = var2 - pDoppler2 * cosTheta; |
---|
| 280 | G4double var = var4*var4 - var3 + pDoppler2 * var3; |
---|
| 281 | if (var > 0.) |
---|
| 282 | { |
---|
| 283 | G4double varSqrt = std::sqrt(var); |
---|
| 284 | G4double scale = photonEnergy0 / var3; |
---|
| 285 | // Random select either root |
---|
| 286 | if (G4UniformRand() < 0.5) photonE = (var4 - varSqrt) * scale; |
---|
| 287 | else photonE = (var4 + varSqrt) * scale; |
---|
| 288 | } |
---|
| 289 | else |
---|
| 290 | { |
---|
| 291 | photonE = -1.; |
---|
| 292 | } |
---|
| 293 | } while ( iteration <= maxDopplerIterations && |
---|
| 294 | (photonE < 0. || photonE > eMax || photonE < eMax*G4UniformRand()) ); |
---|
| 295 | |
---|
| 296 | // End of recalculation of photon energy with Doppler broadening |
---|
| 297 | // Revert to original if maximum number of iterations threshold has been reached |
---|
| 298 | |
---|
| 299 | if (iteration >= maxDopplerIterations) |
---|
| 300 | { |
---|
| 301 | photonE = photonEoriginal; |
---|
| 302 | bindingE = 0.; |
---|
| 303 | } |
---|
| 304 | |
---|
| 305 | // Update G4VParticleChange for the scattered photon |
---|
| 306 | |
---|
| 307 | G4ThreeVector photonDirection1(dirx,diry,dirz); |
---|
| 308 | photonDirection1.rotateUz(photonDirection0); |
---|
| 309 | fParticleChange->ProposeMomentumDirection(photonDirection1) ; |
---|
| 310 | |
---|
| 311 | G4double photonEnergy1 = photonE; |
---|
| 312 | |
---|
| 313 | if (photonEnergy1 > 0.) |
---|
| 314 | { |
---|
| 315 | fParticleChange->SetProposedKineticEnergy(photonEnergy1) ; |
---|
| 316 | } |
---|
| 317 | else |
---|
| 318 | { |
---|
| 319 | fParticleChange->SetProposedKineticEnergy(0.) ; |
---|
| 320 | fParticleChange->ProposeTrackStatus(fStopAndKill); |
---|
| 321 | } |
---|
| 322 | |
---|
| 323 | // Kinematics of the scattered electron |
---|
| 324 | G4double eKineticEnergy = photonEnergy0 - photonEnergy1 - bindingE; |
---|
| 325 | G4double eTotalEnergy = eKineticEnergy + electron_mass_c2; |
---|
| 326 | |
---|
| 327 | G4double electronE = photonEnergy0 * (1. - epsilon) + electron_mass_c2; |
---|
| 328 | G4double electronP2 = electronE*electronE - electron_mass_c2*electron_mass_c2; |
---|
| 329 | G4double sinThetaE = -1.; |
---|
| 330 | G4double cosThetaE = 0.; |
---|
| 331 | if (electronP2 > 0.) |
---|
| 332 | { |
---|
| 333 | cosThetaE = (eTotalEnergy + photonEnergy1 )* (1. - epsilon) / std::sqrt(electronP2); |
---|
| 334 | sinThetaE = -1. * sqrt(1. - cosThetaE * cosThetaE); |
---|
| 335 | } |
---|
| 336 | |
---|
| 337 | G4double eDirX = sinThetaE * std::cos(phi); |
---|
| 338 | G4double eDirY = sinThetaE * std::sin(phi); |
---|
| 339 | G4double eDirZ = cosThetaE; |
---|
| 340 | |
---|
| 341 | G4ThreeVector eDirection(eDirX,eDirY,eDirZ); |
---|
| 342 | eDirection.rotateUz(photonDirection0); |
---|
| 343 | |
---|
| 344 | // SI - The range test has been removed wrt original G4LowEnergyCompton class |
---|
| 345 | |
---|
| 346 | fParticleChange->ProposeLocalEnergyDeposit(bindingE); |
---|
| 347 | |
---|
| 348 | G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),eDirection,eKineticEnergy) ; |
---|
| 349 | fvect->push_back(dp); |
---|
| 350 | |
---|
| 351 | } |
---|
| 352 | |
---|
| 353 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 354 | |
---|
| 355 | G4double G4LivermoreComptonModel::GetMeanFreePath(const G4Track& track, |
---|
| 356 | G4double, // previousStepSize |
---|
| 357 | G4ForceCondition*) |
---|
| 358 | { |
---|
| 359 | const G4DynamicParticle* photon = track.GetDynamicParticle(); |
---|
| 360 | G4double energy = photon->GetKineticEnergy(); |
---|
| 361 | const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple(); |
---|
| 362 | size_t materialIndex = couple->GetIndex(); |
---|
| 363 | |
---|
| 364 | G4double meanFreePath; |
---|
| 365 | if (energy > highEnergyLimit) meanFreePath = meanFreePathTable->FindValue(highEnergyLimit,materialIndex); |
---|
| 366 | else if (energy < lowEnergyLimit) meanFreePath = DBL_MAX; |
---|
| 367 | else meanFreePath = meanFreePathTable->FindValue(energy,materialIndex); |
---|
| 368 | return meanFreePath; |
---|
| 369 | } |
---|
| 370 | |
---|