[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: G4PenelopeComptonModel.cc,v 1.8 2009/10/23 09:29:24 pandola Exp $ |
---|
[1228] | 27 | // GEANT4 tag $Name: geant4-09-03 $ |
---|
[968] | 28 | // |
---|
| 29 | // Author: Luciano Pandola |
---|
| 30 | // |
---|
| 31 | // History: |
---|
| 32 | // -------- |
---|
| 33 | // 02 Oct 2008 L Pandola Migration from process to model |
---|
| 34 | // 28 Oct 2008 L Pandola Treat the database data from Penelope according to the |
---|
| 35 | // original model, namely merging levels below 15 eV in |
---|
| 36 | // a single one. Still, it is not fully compliant with the |
---|
| 37 | // original Penelope model, because plasma excitation is not |
---|
| 38 | // considered. |
---|
| 39 | // 22 Nov 2008 L Pandola Make unit of measurements explicit for binding energies |
---|
| 40 | // that are read from the external files. |
---|
| 41 | // 24 Nov 2008 L Pandola Find a cleaner way to delete vectors. |
---|
[1055] | 42 | // 16 Apr 2009 V Ivanchenko Cleanup initialisation and generation of secondaries: |
---|
| 43 | // - apply internal high-energy limit only in constructor |
---|
| 44 | // - do not apply low-energy limit (default is 0) |
---|
| 45 | // - do not apply production threshold on level of the model |
---|
[1192] | 46 | // 21 Oct 2009 L Pandola Remove un-necessary fUseAtomicDeexcitation flag - now managed by |
---|
| 47 | // G4VEmModel::DeexcitationFlag() |
---|
| 48 | // Add ActivateAuger() method |
---|
[968] | 49 | // |
---|
| 50 | |
---|
| 51 | #include "G4PenelopeComptonModel.hh" |
---|
| 52 | #include "G4ParticleDefinition.hh" |
---|
| 53 | #include "G4MaterialCutsCouple.hh" |
---|
| 54 | #include "G4ProductionCutsTable.hh" |
---|
| 55 | #include "G4DynamicParticle.hh" |
---|
| 56 | #include "G4VEMDataSet.hh" |
---|
| 57 | #include "G4PhysicsTable.hh" |
---|
| 58 | #include "G4ElementTable.hh" |
---|
| 59 | #include "G4Element.hh" |
---|
| 60 | #include "G4PhysicsLogVector.hh" |
---|
| 61 | #include "G4PenelopeIntegrator.hh" |
---|
| 62 | #include "G4AtomicTransitionManager.hh" |
---|
| 63 | #include "G4AtomicShell.hh" |
---|
| 64 | #include "G4Gamma.hh" |
---|
| 65 | #include "G4Electron.hh" |
---|
| 66 | |
---|
| 67 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 68 | |
---|
| 69 | |
---|
| 70 | G4PenelopeComptonModel::G4PenelopeComptonModel(const G4ParticleDefinition*, |
---|
| 71 | const G4String& nam) |
---|
| 72 | :G4VEmModel(nam),ionizationEnergy(0),hartreeFunction(0), |
---|
| 73 | occupationNumber(0),isInitialised(false) |
---|
| 74 | { |
---|
| 75 | fIntrinsicLowEnergyLimit = 100.0*eV; |
---|
| 76 | fIntrinsicHighEnergyLimit = 100.0*GeV; |
---|
[1055] | 77 | // SetLowEnergyLimit(fIntrinsicLowEnergyLimit); |
---|
[968] | 78 | SetHighEnergyLimit(fIntrinsicHighEnergyLimit); |
---|
| 79 | // |
---|
| 80 | energyForIntegration = 0.0; |
---|
| 81 | ZForIntegration = 1; |
---|
| 82 | |
---|
[1192] | 83 | //by default, the model will use atomic deexcitation |
---|
| 84 | SetDeexcitationFlag(true); |
---|
| 85 | ActivateAuger(false); |
---|
| 86 | |
---|
[968] | 87 | verboseLevel= 0; |
---|
| 88 | // Verbosity scale: |
---|
| 89 | // 0 = nothing |
---|
| 90 | // 1 = warning for energy non-conservation |
---|
| 91 | // 2 = details of energy budget |
---|
| 92 | // 3 = calculation of cross sections, file openings, sampling of atoms |
---|
| 93 | // 4 = entering in methods |
---|
| 94 | |
---|
| 95 | //These vectors do not change when materials or cut change. |
---|
| 96 | //Therefore I can read it at the constructor |
---|
| 97 | ionizationEnergy = new std::map<G4int,G4DataVector*>; |
---|
| 98 | hartreeFunction = new std::map<G4int,G4DataVector*>; |
---|
| 99 | occupationNumber = new std::map<G4int,G4DataVector*>; |
---|
| 100 | |
---|
| 101 | ReadData(); //Read data from file |
---|
| 102 | |
---|
| 103 | } |
---|
| 104 | |
---|
| 105 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 106 | |
---|
| 107 | G4PenelopeComptonModel::~G4PenelopeComptonModel() |
---|
| 108 | { |
---|
| 109 | std::map <G4int,G4DataVector*>::iterator i; |
---|
| 110 | for (i=ionizationEnergy->begin();i != ionizationEnergy->end();i++) |
---|
| 111 | if (i->second) delete i->second; |
---|
| 112 | for (i=hartreeFunction->begin();i != hartreeFunction->end();i++) |
---|
| 113 | if (i->second) delete i->second; |
---|
| 114 | for (i=occupationNumber->begin();i != occupationNumber->end();i++) |
---|
| 115 | if (i->second) delete i->second; |
---|
| 116 | |
---|
| 117 | |
---|
| 118 | if (ionizationEnergy) |
---|
| 119 | delete ionizationEnergy; |
---|
| 120 | if (hartreeFunction) |
---|
| 121 | delete hartreeFunction; |
---|
| 122 | if (occupationNumber) |
---|
| 123 | delete occupationNumber; |
---|
| 124 | } |
---|
| 125 | |
---|
| 126 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 127 | |
---|
| 128 | void G4PenelopeComptonModel::Initialise(const G4ParticleDefinition* particle, |
---|
[1055] | 129 | const G4DataVector& cuts) |
---|
[968] | 130 | { |
---|
| 131 | if (verboseLevel > 3) |
---|
| 132 | G4cout << "Calling G4PenelopeComptonModel::Initialise()" << G4endl; |
---|
| 133 | |
---|
| 134 | InitialiseElementSelectors(particle,cuts); |
---|
| 135 | |
---|
[1055] | 136 | if (verboseLevel > 0) { |
---|
| 137 | G4cout << "Penelope Compton model is initialized " << G4endl |
---|
| 138 | << "Energy range: " |
---|
| 139 | << LowEnergyLimit() / keV << " keV - " |
---|
| 140 | << HighEnergyLimit() / GeV << " GeV" |
---|
| 141 | << G4endl; |
---|
| 142 | } |
---|
[968] | 143 | |
---|
| 144 | if(isInitialised) return; |
---|
[1055] | 145 | fParticleChange = GetParticleChangeForGamma(); |
---|
[968] | 146 | isInitialised = true; |
---|
| 147 | } |
---|
| 148 | |
---|
| 149 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 150 | |
---|
| 151 | G4double G4PenelopeComptonModel::ComputeCrossSectionPerAtom( |
---|
| 152 | const G4ParticleDefinition*, |
---|
| 153 | G4double energy, |
---|
| 154 | G4double Z, G4double, |
---|
| 155 | G4double, G4double) |
---|
| 156 | { |
---|
| 157 | // Penelope model to calculate the Compton scattering cross section: |
---|
| 158 | // D. Brusa et al., Nucl. Instrum. Meth. A 379 (1996) 167 |
---|
| 159 | // |
---|
| 160 | // The cross section for Compton scattering is calculated according to the Klein-Nishina |
---|
| 161 | // formula for energy > 5 MeV. |
---|
| 162 | // For E < 5 MeV it is used a parametrization for the differential cross-section dSigma/dOmega, |
---|
| 163 | // which is integrated numerically in cos(theta), G4PenelopeComptonModel::DifferentialCrossSection(). |
---|
| 164 | // The parametrization includes the J(p) |
---|
| 165 | // distribution profiles for the atomic shells, that are tabulated from Hartree-Fock calculations |
---|
| 166 | // from F. Biggs et al., At. Data Nucl. Data Tables 16 (1975) 201 |
---|
| 167 | // |
---|
| 168 | if (verboseLevel > 3) |
---|
| 169 | G4cout << "Calling ComputeCrossSectionPerAtom() of G4PenelopeComptonModel" << G4endl; |
---|
| 170 | |
---|
| 171 | G4int iZ = (G4int) Z; |
---|
| 172 | G4double cs=0.0; |
---|
| 173 | energyForIntegration=energy; |
---|
| 174 | ZForIntegration = iZ; |
---|
| 175 | if (energy< 5*MeV) |
---|
| 176 | { |
---|
| 177 | // numerically integrate differential cross section dSigma/dOmega |
---|
| 178 | G4PenelopeIntegrator<G4PenelopeComptonModel,G4double (G4PenelopeComptonModel::*)(G4double)> |
---|
| 179 | theIntegrator; |
---|
| 180 | cs = theIntegrator.Calculate(this,&G4PenelopeComptonModel::DifferentialCrossSection,-1.0,1.0,1e-05); |
---|
| 181 | } |
---|
| 182 | else |
---|
| 183 | { |
---|
| 184 | // use Klein-Nishina formula |
---|
| 185 | G4double ki=energy/electron_mass_c2; |
---|
| 186 | G4double ki3=ki*ki; |
---|
| 187 | G4double ki2=1.0+2*ki; |
---|
| 188 | G4double ki1=ki3-ki2-1.0; |
---|
| 189 | G4double t0=1.0/(ki2); |
---|
| 190 | G4double csl = 0.5*ki3*t0*t0+ki2*t0+ki1*std::log(t0)-(1.0/t0); |
---|
| 191 | G4int nosc = occupationNumber->find(iZ)->second->size(); |
---|
| 192 | for (G4int i=0;i<nosc;i++) |
---|
| 193 | { |
---|
| 194 | G4double ionEnergy = (*(ionizationEnergy->find(iZ)->second))[i]; |
---|
| 195 | G4double tau=(energy-ionEnergy)/energy; |
---|
| 196 | if (tau > t0) |
---|
| 197 | { |
---|
| 198 | G4double csu = 0.5*ki3*tau*tau+ki2*tau+ki1*std::log(tau)-(1.0/tau); |
---|
| 199 | G4int f = (G4int) (*(occupationNumber->find(iZ)->second))[i]; |
---|
| 200 | cs = cs + f*(csu-csl); |
---|
| 201 | } |
---|
| 202 | } |
---|
| 203 | cs=pi*classic_electr_radius*classic_electr_radius*cs/(ki*ki3); |
---|
| 204 | } |
---|
| 205 | |
---|
| 206 | if (verboseLevel > 2) |
---|
| 207 | G4cout << "Compton cross Section at " << energy/keV << " keV for Z=" << Z << |
---|
| 208 | " = " << cs/barn << " barn" << G4endl; |
---|
| 209 | return cs; |
---|
| 210 | } |
---|
| 211 | |
---|
| 212 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 213 | |
---|
| 214 | void G4PenelopeComptonModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect, |
---|
| 215 | const G4MaterialCutsCouple* couple, |
---|
| 216 | const G4DynamicParticle* aDynamicGamma, |
---|
| 217 | G4double, |
---|
| 218 | G4double) |
---|
| 219 | { |
---|
| 220 | |
---|
| 221 | // Penelope model to sample the Compton scattering final state. |
---|
| 222 | // D. Brusa et al., Nucl. Instrum. Meth. A 379 (1996) 167 |
---|
| 223 | // The model determines also the original shell from which the electron is expelled, |
---|
| 224 | // in order to produce fluorescence de-excitation (from G4DeexcitationManager) |
---|
| 225 | // |
---|
| 226 | // The final state for Compton scattering is calculated according to the Klein-Nishina |
---|
| 227 | // formula for energy > 5 MeV. In this case, the Doppler broadening is negligible and |
---|
| 228 | // one can assume that the target electron is at rest. |
---|
| 229 | // For E < 5 MeV it is used the parametrization for the differential cross-section dSigma/dOmega, |
---|
| 230 | // to sample the scattering angle and the energy of the emerging electron, which is |
---|
| 231 | // G4PenelopeComptonModel::DifferentialCrossSection(). The rejection method is |
---|
| 232 | // used to sample cos(theta). The efficiency increases monotonically with photon energy and is |
---|
| 233 | // nearly independent on the Z; typical values are 35%, 80% and 95% for 1 keV, 1 MeV and 10 MeV, |
---|
| 234 | // respectively. |
---|
| 235 | // The parametrization includes the J(p) distribution profiles for the atomic shells, that are |
---|
| 236 | // tabulated |
---|
| 237 | // from Hartree-Fock calculations from F. Biggs et al., At. Data Nucl. Data Tables 16 (1975) 201. |
---|
| 238 | // Doppler broadening is included. |
---|
| 239 | // |
---|
| 240 | |
---|
| 241 | if (verboseLevel > 3) |
---|
| 242 | G4cout << "Calling SampleSecondaries() of G4PenelopeComptonModel" << G4endl; |
---|
| 243 | |
---|
| 244 | G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy(); |
---|
| 245 | |
---|
[1055] | 246 | if (photonEnergy0 <= fIntrinsicLowEnergyLimit) |
---|
| 247 | { |
---|
[968] | 248 | fParticleChange->ProposeTrackStatus(fStopAndKill); |
---|
| 249 | fParticleChange->SetProposedKineticEnergy(0.); |
---|
| 250 | fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0); |
---|
| 251 | return ; |
---|
[1055] | 252 | } |
---|
[968] | 253 | |
---|
| 254 | G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection(); |
---|
| 255 | |
---|
| 256 | // Select randomly one element in the current material |
---|
| 257 | if (verboseLevel > 2) |
---|
| 258 | G4cout << "Going to select element in " << couple->GetMaterial()->GetName() << G4endl; |
---|
| 259 | // atom can be selected effitiantly if element selectors are initialised |
---|
[1055] | 260 | const G4Element* anElement = |
---|
| 261 | SelectRandomAtom(couple,G4Gamma::GammaDefinition(),photonEnergy0); |
---|
[968] | 262 | G4int Z = (G4int) anElement->GetZ(); |
---|
| 263 | if (verboseLevel > 2) |
---|
| 264 | G4cout << "Selected " << anElement->GetName() << G4endl; |
---|
| 265 | |
---|
| 266 | const G4int nmax = 64; |
---|
| 267 | G4double rn[nmax],pac[nmax]; |
---|
| 268 | |
---|
| 269 | G4double ki,ki1,ki2,ki3,taumin,a1,a2; |
---|
| 270 | G4double tau,TST; |
---|
| 271 | G4double S=0.0; |
---|
| 272 | G4double epsilon,cosTheta; |
---|
| 273 | G4double harFunc = 0.0; |
---|
| 274 | G4int occupNb= 0; |
---|
| 275 | G4double ionEnergy=0.0; |
---|
| 276 | G4int nosc = occupationNumber->find(Z)->second->size(); |
---|
| 277 | G4int iosc = nosc; |
---|
| 278 | ki = photonEnergy0/electron_mass_c2; |
---|
| 279 | ki2 = 2*ki+1.0; |
---|
| 280 | ki3 = ki*ki; |
---|
| 281 | ki1 = ki3-ki2-1.0; |
---|
| 282 | taumin = 1.0/ki2; |
---|
| 283 | a1 = std::log(ki2); |
---|
| 284 | a2 = a1+2.0*ki*(1.0+ki)/(ki2*ki2); |
---|
| 285 | //If the incoming photon is above 5 MeV, the quicker approach based on the |
---|
| 286 | //pure Klein-Nishina formula is used |
---|
| 287 | if (photonEnergy0 > 5*MeV) |
---|
| 288 | { |
---|
| 289 | do{ |
---|
| 290 | do{ |
---|
| 291 | if ((a2*G4UniformRand()) < a1) |
---|
| 292 | { |
---|
| 293 | tau = std::pow(taumin,G4UniformRand()); |
---|
| 294 | } |
---|
| 295 | else |
---|
| 296 | { |
---|
| 297 | tau = std::sqrt(1.0+G4UniformRand()*(taumin*taumin-1.0)); |
---|
| 298 | } |
---|
| 299 | //rejection function |
---|
| 300 | TST = (1+tau*(ki1+tau*(ki2+tau*ki3)))/(ki3*tau*(1.0+tau*tau)); |
---|
| 301 | }while (G4UniformRand()> TST); |
---|
| 302 | epsilon=tau; |
---|
| 303 | cosTheta = 1.0 - (1.0-tau)/(ki*tau); |
---|
| 304 | //Target shell electrons |
---|
| 305 | TST = Z*G4UniformRand(); |
---|
| 306 | iosc = nosc; |
---|
| 307 | S=0.0; |
---|
| 308 | for (G4int j=0;j<nosc;j++) |
---|
| 309 | { |
---|
| 310 | occupNb = (G4int) (*(occupationNumber->find(Z)->second))[j]; |
---|
| 311 | S = S + occupNb; |
---|
| 312 | if (S > TST) iosc = j; |
---|
| 313 | if (S > TST) break; |
---|
| 314 | } |
---|
| 315 | ionEnergy = (*(ionizationEnergy->find(Z)->second))[iosc]; |
---|
| 316 | }while((epsilon*photonEnergy0-photonEnergy0+ionEnergy) >0); |
---|
| 317 | } |
---|
[1055] | 318 | else //photonEnergy0 < 5 MeV |
---|
[968] | 319 | { |
---|
| 320 | //Incoherent scattering function for theta=PI |
---|
| 321 | G4double s0=0.0; |
---|
| 322 | G4double pzomc=0.0,rni=0.0; |
---|
| 323 | G4double aux=0.0; |
---|
| 324 | for (G4int i=0;i<nosc;i++){ |
---|
| 325 | ionEnergy = (*(ionizationEnergy->find(Z)->second))[i]; |
---|
| 326 | if (photonEnergy0 > ionEnergy) |
---|
| 327 | { |
---|
| 328 | G4double aux = photonEnergy0*(photonEnergy0-ionEnergy)*2.0; |
---|
| 329 | harFunc = (*(hartreeFunction->find(Z)->second))[i]/fine_structure_const; |
---|
| 330 | occupNb = (G4int) (*(occupationNumber->find(Z)->second))[i]; |
---|
| 331 | pzomc = harFunc*(aux-electron_mass_c2*ionEnergy)/ |
---|
| 332 | (electron_mass_c2*std::sqrt(2.0*aux+ionEnergy*ionEnergy)); |
---|
| 333 | if (pzomc > 0) |
---|
| 334 | { |
---|
| 335 | rni = 1.0-0.5*std::exp(0.5-(std::sqrt(0.5)+std::sqrt(2.0)*pzomc)* |
---|
| 336 | (std::sqrt(0.5)+std::sqrt(2.0)*pzomc)); |
---|
| 337 | } |
---|
| 338 | else |
---|
| 339 | { |
---|
| 340 | rni = 0.5*std::exp(0.5-(std::sqrt(0.5)-std::sqrt(2.0)*pzomc)* |
---|
| 341 | (std::sqrt(0.5)-std::sqrt(2.0)*pzomc)); |
---|
| 342 | } |
---|
| 343 | s0 = s0 + occupNb*rni; |
---|
| 344 | } |
---|
| 345 | } |
---|
| 346 | |
---|
| 347 | //Sampling tau |
---|
| 348 | G4double cdt1; |
---|
| 349 | do |
---|
| 350 | { |
---|
| 351 | if ((G4UniformRand()*a2) < a1) |
---|
| 352 | { |
---|
| 353 | tau = std::pow(taumin,G4UniformRand()); |
---|
| 354 | } |
---|
| 355 | else |
---|
| 356 | { |
---|
| 357 | tau = std::sqrt(1.0+G4UniformRand()*(taumin*taumin-1.0)); |
---|
| 358 | } |
---|
| 359 | cdt1 = (1.0-tau)/(ki*tau); |
---|
| 360 | S=0.0; |
---|
| 361 | //Incoherent scattering function |
---|
| 362 | for (G4int i=0;i<nosc;i++){ |
---|
| 363 | ionEnergy = (*(ionizationEnergy->find(Z)->second))[i]; |
---|
| 364 | if (photonEnergy0 > ionEnergy) //sum only on excitable levels |
---|
| 365 | { |
---|
| 366 | aux = photonEnergy0*(photonEnergy0-ionEnergy)*cdt1; |
---|
| 367 | harFunc = (*(hartreeFunction->find(Z)->second))[i]/fine_structure_const; |
---|
| 368 | occupNb = (G4int) (*(occupationNumber->find(Z)->second))[i]; |
---|
| 369 | pzomc = harFunc*(aux-electron_mass_c2*ionEnergy)/ |
---|
| 370 | (electron_mass_c2*std::sqrt(2.0*aux+ionEnergy*ionEnergy)); |
---|
| 371 | if (pzomc > 0) |
---|
| 372 | { |
---|
| 373 | rn[i] = 1.0-0.5*std::exp(0.5-(std::sqrt(0.5)+std::sqrt(2.0)*pzomc)* |
---|
| 374 | (std::sqrt(0.5)+std::sqrt(2.0)*pzomc)); |
---|
| 375 | } |
---|
| 376 | else |
---|
| 377 | { |
---|
| 378 | rn[i] = 0.5*std::exp(0.5-(std::sqrt(0.5)-std::sqrt(2.0)*pzomc)* |
---|
| 379 | (std::sqrt(0.5)-std::sqrt(2.0)*pzomc)); |
---|
| 380 | } |
---|
| 381 | S = S + occupNb*rn[i]; |
---|
| 382 | pac[i] = S; |
---|
| 383 | } |
---|
| 384 | else |
---|
| 385 | { |
---|
| 386 | pac[i] = S-(1e-06); |
---|
| 387 | } |
---|
| 388 | } |
---|
| 389 | //Rejection function |
---|
| 390 | TST = S*(1.0+tau*(ki1+tau*(ki2+tau*ki3)))/(ki3*tau*(1.0+tau*tau)); |
---|
| 391 | }while ((G4UniformRand()*s0) > TST); |
---|
| 392 | //Target electron shell |
---|
| 393 | cosTheta = 1.0 - cdt1; |
---|
| 394 | G4double fpzmax=0.0,fpz=0.0; |
---|
| 395 | G4double A=0.0; |
---|
| 396 | do |
---|
| 397 | { |
---|
| 398 | do |
---|
| 399 | { |
---|
| 400 | TST =S*G4UniformRand(); |
---|
| 401 | iosc=nosc; |
---|
| 402 | for (G4int i=0;i<nosc;i++){ |
---|
| 403 | if (pac[i]>TST) iosc = i; |
---|
| 404 | if (pac[i]>TST) break; |
---|
| 405 | } |
---|
| 406 | A = G4UniformRand()*rn[iosc]; |
---|
| 407 | harFunc = (*(hartreeFunction->find(Z)->second))[iosc]/fine_structure_const; |
---|
| 408 | occupNb = (G4int) (*(occupationNumber->find(Z)->second))[iosc]; |
---|
| 409 | if (A < 0.5) { |
---|
| 410 | pzomc = (std::sqrt(0.5)-std::sqrt(0.5-std::log(2.0*A)))/ |
---|
| 411 | (std::sqrt(2.0)*harFunc); |
---|
| 412 | } |
---|
| 413 | else |
---|
| 414 | { |
---|
| 415 | pzomc = (std::sqrt(0.5-std::log(2.0-2.0*A))-std::sqrt(0.5))/ |
---|
| 416 | (std::sqrt(2.0)*harFunc); |
---|
| 417 | } |
---|
| 418 | } while (pzomc < -1); |
---|
| 419 | // F(EP) rejection |
---|
| 420 | G4double XQC = 1.0+tau*(tau-2.0*cosTheta); |
---|
| 421 | G4double AF = std::sqrt(XQC)*(1.0+tau*(tau-cosTheta)/XQC); |
---|
| 422 | if (AF > 0) { |
---|
| 423 | fpzmax = 1.0+AF*0.2; |
---|
| 424 | } |
---|
| 425 | else |
---|
| 426 | { |
---|
| 427 | fpzmax = 1.0-AF*0.2; |
---|
| 428 | } |
---|
| 429 | fpz = 1.0+AF*std::max(std::min(pzomc,0.2),-0.2); |
---|
| 430 | }while ((fpzmax*G4UniformRand())>fpz); |
---|
| 431 | |
---|
| 432 | //Energy of the scattered photon |
---|
| 433 | G4double T = pzomc*pzomc; |
---|
| 434 | G4double b1 = 1.0-T*tau*tau; |
---|
| 435 | G4double b2 = 1.0-T*tau*cosTheta; |
---|
| 436 | if (pzomc > 0.0) |
---|
| 437 | { |
---|
| 438 | epsilon = (tau/b1)*(b2+std::sqrt(std::abs(b2*b2-b1*(1.0-T)))); |
---|
| 439 | } |
---|
| 440 | else |
---|
| 441 | { |
---|
| 442 | epsilon = (tau/b1)*(b2-std::sqrt(std::abs(b2*b2-b1*(1.0-T)))); |
---|
| 443 | } |
---|
| 444 | } |
---|
| 445 | |
---|
| 446 | |
---|
| 447 | //Ok, the kinematics has been calculated. |
---|
| 448 | G4double sinTheta = std::sqrt(1-cosTheta*cosTheta); |
---|
| 449 | G4double phi = twopi * G4UniformRand() ; |
---|
| 450 | G4double dirx = sinTheta * std::cos(phi); |
---|
| 451 | G4double diry = sinTheta * std::sin(phi); |
---|
| 452 | G4double dirz = cosTheta ; |
---|
| 453 | |
---|
| 454 | // Update G4VParticleChange for the scattered photon |
---|
| 455 | G4ThreeVector photonDirection1(dirx,diry,dirz); |
---|
| 456 | photonDirection1.rotateUz(photonDirection0); |
---|
| 457 | fParticleChange->ProposeMomentumDirection(photonDirection1) ; |
---|
| 458 | G4double photonEnergy1 = epsilon * photonEnergy0; |
---|
| 459 | |
---|
| 460 | if (photonEnergy1 > 0.) |
---|
| 461 | { |
---|
| 462 | fParticleChange->SetProposedKineticEnergy(photonEnergy1) ; |
---|
| 463 | } |
---|
| 464 | else |
---|
| 465 | { |
---|
| 466 | fParticleChange->SetProposedKineticEnergy(0.) ; |
---|
| 467 | fParticleChange->ProposeTrackStatus(fStopAndKill); |
---|
| 468 | } |
---|
| 469 | |
---|
| 470 | //Create scattered electron |
---|
| 471 | G4double diffEnergy = photonEnergy0*(1-epsilon); |
---|
| 472 | ionEnergy = (*(ionizationEnergy->find(Z)->second))[iosc]; |
---|
[1055] | 473 | G4double Q2 = |
---|
| 474 | photonEnergy0*photonEnergy0+photonEnergy1*(photonEnergy1-2.0*photonEnergy0*cosTheta); |
---|
[968] | 475 | G4double cosThetaE; //scattering angle for the electron |
---|
| 476 | if (Q2 > 1.0e-12) |
---|
| 477 | { |
---|
| 478 | cosThetaE = (photonEnergy0-photonEnergy1*cosTheta)/std::sqrt(Q2); |
---|
| 479 | } |
---|
| 480 | else |
---|
| 481 | { |
---|
| 482 | cosThetaE = 1.0; |
---|
| 483 | } |
---|
| 484 | G4double sinThetaE = std::sqrt(1-cosThetaE*cosThetaE); |
---|
| 485 | |
---|
| 486 | //initialize here, then check photons created by Atomic-Deexcitation, and the final state e- |
---|
| 487 | std::vector<G4DynamicParticle*>* photonVector=0; |
---|
| 488 | |
---|
| 489 | const G4AtomicTransitionManager* transitionManager = G4AtomicTransitionManager::Instance(); |
---|
| 490 | const G4AtomicShell* shell = transitionManager->Shell(Z,iosc); |
---|
| 491 | G4double bindingEnergy = shell->BindingEnergy(); |
---|
| 492 | G4int shellId = shell->ShellId(); |
---|
| 493 | G4double ionEnergyInPenelopeDatabase = ionEnergy; |
---|
[1055] | 494 | //protection against energy non-conservation |
---|
| 495 | ionEnergy = std::max(bindingEnergy,ionEnergyInPenelopeDatabase); |
---|
[968] | 496 | |
---|
[1055] | 497 | //subtract the excitation energy. If not emitted by fluorescence |
---|
[968] | 498 | //the ionization energy is deposited as local energy deposition |
---|
[1055] | 499 | G4double eKineticEnergy = diffEnergy - ionEnergy; |
---|
[968] | 500 | G4double localEnergyDeposit = ionEnergy; |
---|
| 501 | G4double energyInFluorescence = 0.; //testing purposes only |
---|
| 502 | |
---|
| 503 | if (eKineticEnergy < 0) |
---|
| 504 | { |
---|
[1055] | 505 | //It means that there was some problem/mismatch between the two databases. |
---|
| 506 | //Try to make it work |
---|
[968] | 507 | //In this case available Energy (diffEnergy) < ionEnergy |
---|
| 508 | //Full residual energy is deposited locally |
---|
| 509 | localEnergyDeposit = diffEnergy; |
---|
| 510 | eKineticEnergy = 0.0; |
---|
| 511 | } |
---|
| 512 | |
---|
| 513 | //the local energy deposit is what remains: part of this may be spent for fluorescence. |
---|
[1055] | 514 | if(DeexcitationFlag() && Z > 5) { |
---|
[968] | 515 | |
---|
[1055] | 516 | const G4ProductionCutsTable* theCoupleTable= |
---|
| 517 | G4ProductionCutsTable::GetProductionCutsTable(); |
---|
[968] | 518 | |
---|
[1055] | 519 | size_t index = couple->GetIndex(); |
---|
| 520 | G4double cutg = (*(theCoupleTable->GetEnergyCutsVector(0)))[index]; |
---|
| 521 | G4double cute = (*(theCoupleTable->GetEnergyCutsVector(1)))[index]; |
---|
| 522 | |
---|
| 523 | // Generation of fluorescence |
---|
| 524 | // Data in EADL are available only for Z > 5 |
---|
| 525 | // Protection to avoid generating photons in the unphysical case of |
---|
| 526 | // shell binding energy > photon energy |
---|
| 527 | if (localEnergyDeposit > cutg || localEnergyDeposit > cute) |
---|
| 528 | { |
---|
| 529 | G4DynamicParticle* aPhoton; |
---|
| 530 | deexcitationManager.SetCutForSecondaryPhotons(cutg); |
---|
| 531 | deexcitationManager.SetCutForAugerElectrons(cute); |
---|
[968] | 532 | |
---|
[1055] | 533 | photonVector = deexcitationManager.GenerateParticles(Z,shellId); |
---|
| 534 | if(photonVector) |
---|
| 535 | { |
---|
| 536 | size_t nPhotons = photonVector->size(); |
---|
| 537 | for (size_t k=0; k<nPhotons; k++) |
---|
[968] | 538 | { |
---|
[1055] | 539 | aPhoton = (*photonVector)[k]; |
---|
| 540 | if (aPhoton) |
---|
[968] | 541 | { |
---|
[1055] | 542 | G4double itsEnergy = aPhoton->GetKineticEnergy(); |
---|
| 543 | if (itsEnergy <= localEnergyDeposit) |
---|
| 544 | { |
---|
| 545 | localEnergyDeposit -= itsEnergy; |
---|
| 546 | if (aPhoton->GetDefinition() == G4Gamma::Gamma()) |
---|
| 547 | energyInFluorescence += itsEnergy;; |
---|
| 548 | fvect->push_back(aPhoton); |
---|
| 549 | } |
---|
| 550 | else |
---|
| 551 | { |
---|
| 552 | delete aPhoton; |
---|
| 553 | (*photonVector)[k]=0; |
---|
| 554 | } |
---|
[968] | 555 | } |
---|
| 556 | } |
---|
[1055] | 557 | delete photonVector; |
---|
[968] | 558 | } |
---|
[1055] | 559 | } |
---|
| 560 | } |
---|
[968] | 561 | |
---|
| 562 | //Produce explicitely the electron only if its proposed kinetic energy is |
---|
| 563 | //above the cut, otherwise add local energy deposition |
---|
| 564 | G4DynamicParticle* electron = 0; |
---|
[1055] | 565 | // if (eKineticEnergy > cutE) // VI: may be consistency problem if cut is applied here |
---|
| 566 | if (eKineticEnergy > 0.0) |
---|
[968] | 567 | { |
---|
| 568 | G4double xEl = sinThetaE * std::cos(phi+pi); |
---|
| 569 | G4double yEl = sinThetaE * std::sin(phi+pi); |
---|
| 570 | G4double zEl = cosThetaE; |
---|
| 571 | G4ThreeVector eDirection(xEl,yEl,zEl); //electron direction |
---|
| 572 | eDirection.rotateUz(photonDirection0); |
---|
| 573 | electron = new G4DynamicParticle (G4Electron::Electron(), |
---|
| 574 | eDirection,eKineticEnergy) ; |
---|
| 575 | fvect->push_back(electron); |
---|
| 576 | } |
---|
| 577 | else |
---|
| 578 | { |
---|
| 579 | localEnergyDeposit += eKineticEnergy; |
---|
| 580 | } |
---|
| 581 | |
---|
| 582 | if (localEnergyDeposit < 0) |
---|
| 583 | { |
---|
| 584 | G4cout << "WARNING-" |
---|
| 585 | << "G4PenelopeComptonModel::SampleSecondaries - Negative energy deposit" |
---|
| 586 | << G4endl; |
---|
| 587 | localEnergyDeposit=0.; |
---|
| 588 | } |
---|
| 589 | fParticleChange->ProposeLocalEnergyDeposit(localEnergyDeposit); |
---|
| 590 | |
---|
| 591 | G4double electronEnergy = 0.; |
---|
| 592 | if (verboseLevel > 1) |
---|
| 593 | { |
---|
| 594 | G4cout << "-----------------------------------------------------------" << G4endl; |
---|
| 595 | G4cout << "Energy balance from G4PenelopeCompton" << G4endl; |
---|
| 596 | G4cout << "Incoming photon energy: " << photonEnergy0/keV << " keV" << G4endl; |
---|
| 597 | G4cout << "-----------------------------------------------------------" << G4endl; |
---|
| 598 | G4cout << "Scattered photon: " << photonEnergy1/keV << " keV" << G4endl; |
---|
| 599 | if (electron) |
---|
| 600 | electronEnergy = eKineticEnergy; |
---|
| 601 | G4cout << "Scattered electron " << electronEnergy/keV << " keV" << G4endl; |
---|
| 602 | G4cout << "Fluorescence: " << energyInFluorescence/keV << " keV" << G4endl; |
---|
| 603 | G4cout << "Local energy deposit " << localEnergyDeposit/keV << " keV" << G4endl; |
---|
| 604 | G4cout << "Total final state: " << (photonEnergy1+electronEnergy+energyInFluorescence+ |
---|
| 605 | localEnergyDeposit)/keV << |
---|
| 606 | " keV" << G4endl; |
---|
| 607 | G4cout << "-----------------------------------------------------------" << G4endl; |
---|
| 608 | } |
---|
| 609 | if (verboseLevel > 0) |
---|
| 610 | { |
---|
| 611 | G4double energyDiff = std::fabs(photonEnergy1+ |
---|
| 612 | electronEnergy+energyInFluorescence+ |
---|
| 613 | localEnergyDeposit-photonEnergy0); |
---|
| 614 | if (energyDiff > 0.05*keV) |
---|
| 615 | G4cout << "Warning from G4PenelopeCompton: problem with energy conservation: " << |
---|
| 616 | (photonEnergy1+electronEnergy+energyInFluorescence+localEnergyDeposit)/keV << |
---|
| 617 | " keV (final) vs. " << |
---|
| 618 | photonEnergy0/keV << " keV (initial)" << G4endl; |
---|
| 619 | } |
---|
| 620 | } |
---|
| 621 | |
---|
| 622 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 623 | |
---|
| 624 | void G4PenelopeComptonModel::ReadData() |
---|
| 625 | { |
---|
| 626 | char* path = getenv("G4LEDATA"); |
---|
| 627 | if (!path) |
---|
| 628 | { |
---|
| 629 | G4String excep = "G4PenelopeComptonModel - G4LEDATA environment variable not set!"; |
---|
| 630 | G4Exception(excep); |
---|
| 631 | } |
---|
| 632 | G4String pathString(path); |
---|
| 633 | G4String pathFile = pathString + "/penelope/compton-pen.dat"; |
---|
| 634 | std::ifstream file(pathFile); |
---|
| 635 | |
---|
| 636 | if (!file.is_open()) |
---|
| 637 | { |
---|
| 638 | G4String excep = "G4PenelopeComptonModel - data file " + pathFile + " not found!"; |
---|
| 639 | G4Exception(excep); |
---|
| 640 | } |
---|
| 641 | |
---|
| 642 | G4int k1,test,test1; |
---|
| 643 | G4double a1,a2; |
---|
| 644 | G4int Z=1,nLevels=0; |
---|
| 645 | |
---|
| 646 | if (!ionizationEnergy || !hartreeFunction || !occupationNumber) |
---|
| 647 | { |
---|
| 648 | G4String excep = "G4PenelopeComptonModel: problem with reading data from file"; |
---|
| 649 | G4Exception(excep); |
---|
| 650 | } |
---|
| 651 | |
---|
| 652 | do{ |
---|
| 653 | G4double harOfElectronsBelowThreshold = 0; |
---|
| 654 | G4int nbOfElectronsBelowThreshold = 0; |
---|
| 655 | G4DataVector* occVector = new G4DataVector; |
---|
| 656 | G4DataVector* harVector = new G4DataVector; |
---|
| 657 | G4DataVector* bindingEVector = new G4DataVector; |
---|
| 658 | file >> Z >> nLevels; |
---|
| 659 | for (G4int h=0;h<nLevels;h++) |
---|
| 660 | { |
---|
| 661 | file >> k1 >> a1 >> a2; |
---|
| 662 | //Make explicit unit of measurements for ionisation energy, which is MeV |
---|
| 663 | a1 *= MeV; |
---|
| 664 | if (a1 > 15*eV) |
---|
| 665 | { |
---|
| 666 | occVector->push_back((G4double) k1); |
---|
| 667 | bindingEVector->push_back(a1); |
---|
| 668 | harVector->push_back(a2); |
---|
| 669 | } |
---|
| 670 | else |
---|
| 671 | { |
---|
| 672 | nbOfElectronsBelowThreshold += k1; |
---|
| 673 | harOfElectronsBelowThreshold += k1*a2; |
---|
| 674 | } |
---|
| 675 | } |
---|
| 676 | //Add the "final" level |
---|
| 677 | if (nbOfElectronsBelowThreshold) |
---|
| 678 | { |
---|
| 679 | occVector->push_back(nbOfElectronsBelowThreshold); |
---|
| 680 | bindingEVector->push_back(0*eV); |
---|
| 681 | G4double averageHartree = |
---|
| 682 | harOfElectronsBelowThreshold/((G4double) nbOfElectronsBelowThreshold); |
---|
| 683 | harVector->push_back(averageHartree); |
---|
| 684 | } |
---|
| 685 | //Ok, done for element Z |
---|
| 686 | occupationNumber->insert(std::make_pair(Z,occVector)); |
---|
| 687 | ionizationEnergy->insert(std::make_pair(Z,bindingEVector)); |
---|
| 688 | hartreeFunction->insert(std::make_pair(Z,harVector)); |
---|
| 689 | file >> test >> test1; //-1 -1 close the data for each Z |
---|
| 690 | if (test > 0) { |
---|
| 691 | G4String excep = "G4PenelopeComptonModel - data file corrupted!"; |
---|
| 692 | G4Exception(excep); |
---|
| 693 | } |
---|
| 694 | }while (test != -2); //the very last Z is closed with -2 instead of -1 |
---|
| 695 | file.close(); |
---|
| 696 | if (verboseLevel > 2) |
---|
| 697 | { |
---|
| 698 | G4cout << "Data from G4PenelopeComptonModel read " << G4endl; |
---|
| 699 | } |
---|
| 700 | } |
---|
| 701 | |
---|
| 702 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 703 | |
---|
| 704 | G4double G4PenelopeComptonModel::DifferentialCrossSection(G4double cosTheta) |
---|
| 705 | { |
---|
| 706 | // |
---|
| 707 | // Penelope model for the Compton scattering differential cross section |
---|
| 708 | // dSigma/dOmega. |
---|
| 709 | // D. Brusa et al., Nucl. Instrum. Meth. A 379 (1996) 167 |
---|
| 710 | // The parametrization includes the J(p) distribution profiles for the atomic shells, |
---|
| 711 | // that are tabulated from Hartree-Fock calculations |
---|
| 712 | // from F. Biggs et al., At. Data Nucl. Data Tables 16 (1975) 201 |
---|
| 713 | // |
---|
| 714 | const G4double k2 = std::sqrt(2.0); |
---|
| 715 | const G4double k1 = std::sqrt(0.5); |
---|
| 716 | const G4double k12 = 0.5; |
---|
| 717 | G4double cdt1 = 1.0-cosTheta; |
---|
| 718 | G4double energy = energyForIntegration; |
---|
| 719 | G4int Z = ZForIntegration; |
---|
| 720 | //energy of Compton line; |
---|
| 721 | G4double EOEC = 1.0+(energy/electron_mass_c2)*cdt1; |
---|
| 722 | G4double ECOE = 1.0/EOEC; |
---|
| 723 | //Incoherent scattering function (analytical profile) |
---|
| 724 | G4double sia = 0.0; |
---|
| 725 | G4int nosc = occupationNumber->find(Z)->second->size(); |
---|
| 726 | for (G4int i=0;i<nosc;i++){ |
---|
| 727 | G4double ionEnergy = (*(ionizationEnergy->find(Z)->second))[i]; |
---|
| 728 | //Sum only of those shells for which E>Eion |
---|
| 729 | if (energy > ionEnergy) |
---|
| 730 | { |
---|
| 731 | G4double aux = energy * (energy-ionEnergy)*cdt1; |
---|
| 732 | G4double Pzimax = |
---|
| 733 | (aux - electron_mass_c2*ionEnergy)/(electron_mass_c2*std::sqrt(2*aux+ionEnergy*ionEnergy)); |
---|
| 734 | G4double harFunc = (*(hartreeFunction->find(Z)->second))[i]/fine_structure_const; |
---|
| 735 | G4int occupNb = (G4int) (*(occupationNumber->find(Z)->second))[i]; |
---|
| 736 | G4double x = harFunc*Pzimax; |
---|
| 737 | G4double siap = 0; |
---|
| 738 | if (x > 0) |
---|
| 739 | { |
---|
| 740 | siap = 1.0-0.5*std::exp(k12-(k1+k2*x)*(k1+k2*x)); |
---|
| 741 | } |
---|
| 742 | else |
---|
| 743 | { |
---|
| 744 | siap = 0.5*std::exp(k12-(k1-k2*x)*(k1-k2*x)); |
---|
| 745 | } |
---|
| 746 | sia = sia + occupNb*siap; //sum of all contributions; |
---|
| 747 | } |
---|
| 748 | } |
---|
| 749 | G4double XKN = EOEC+ECOE-1+cosTheta*cosTheta; |
---|
| 750 | G4double diffCS = pi*classic_electr_radius*classic_electr_radius*ECOE*ECOE*XKN*sia; |
---|
| 751 | return diffCS; |
---|
| 752 | } |
---|
| 753 | |
---|
[1192] | 754 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 755 | |
---|
| 756 | void G4PenelopeComptonModel::ActivateAuger(G4bool augerbool) |
---|
| 757 | { |
---|
| 758 | if (!DeexcitationFlag() && augerbool) |
---|
| 759 | { |
---|
| 760 | G4cout << "WARNING - G4PenelopeComptonModel" << G4endl; |
---|
| 761 | G4cout << "The use of the Atomic Deexcitation Manager is set to false " << G4endl; |
---|
| 762 | G4cout << "Therefore, Auger electrons will be not generated anyway" << G4endl; |
---|
| 763 | } |
---|
| 764 | deexcitationManager.ActivateAugerElectronProduction(augerbool); |
---|
| 765 | if (verboseLevel > 1) |
---|
| 766 | G4cout << "Auger production set to " << augerbool << G4endl; |
---|
| 767 | } |
---|
| 768 | |
---|
| 769 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|