[819] | 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 | // |
---|
[961] | 26 | // $Id: G4PenelopeCompton.cc,v 1.33 2008/06/03 15:44:25 pandola Exp $ |
---|
[1007] | 27 | // GEANT4 tag $Name: geant4-09-02 $ |
---|
[819] | 28 | // |
---|
| 29 | // Author: Luciano Pandola |
---|
| 30 | // |
---|
| 31 | // History: |
---|
| 32 | // -------- |
---|
| 33 | // 12 Feb 2003 MG Pia const argument in SelectRandomAtomForCompton |
---|
| 34 | // Migration to "cuts per region" |
---|
| 35 | // 14 Feb 2003 MG Pia Corrected compilation errors and warnings |
---|
| 36 | // from SUN |
---|
| 37 | // Modified some variables to lowercase initial |
---|
| 38 | // 10 Mar 2003 V.Ivanchenko Remove CutPerMaterial warning |
---|
| 39 | // 13 Mar 2003 L.Pandola Code "cleaned" |
---|
| 40 | // 20 Mar 2003 L.Pandola ReadData() changed (performance improved) |
---|
| 41 | // 26 Mar 2003 L.Pandola Added fluorescence |
---|
| 42 | // 24 May 2003 MGP Removed memory leak |
---|
| 43 | // 09 Mar 2004 L.Pandola Bug fixed in the generation of final state |
---|
| 44 | // (bug report # 585) |
---|
| 45 | // 17 Mar 2004 L.Pandola Removed unnecessary calls to std::pow(a,b) |
---|
| 46 | // 18 Mar 2004 L.Pandola Use of std::map (code review) |
---|
[961] | 47 | // 26 Mar 2008 L.Pandola Add boolean flag to control atomic de-excitation |
---|
| 48 | // 27 Mar 2008 L.Pandola Re-named some variables to improve readability, |
---|
| 49 | // and check for strict energy conservation |
---|
| 50 | // 03 Jun 2008 L.Pandola Added further protection against non-conservation |
---|
| 51 | // of energy: it may happen because ionization energy |
---|
| 52 | // from de-excitation manager and from Penelope internal |
---|
| 53 | // database do not match (difference is <10 eV, but may |
---|
| 54 | // give a e- with negative kinetic energy). |
---|
[819] | 55 | // |
---|
| 56 | // ------------------------------------------------------------------- |
---|
| 57 | |
---|
| 58 | #include "G4PenelopeCompton.hh" |
---|
| 59 | #include "Randomize.hh" |
---|
| 60 | #include "G4ParticleDefinition.hh" |
---|
| 61 | #include "G4Track.hh" |
---|
| 62 | #include "G4Step.hh" |
---|
| 63 | #include "G4ForceCondition.hh" |
---|
| 64 | #include "G4Gamma.hh" |
---|
| 65 | #include "G4Electron.hh" |
---|
| 66 | #include "G4DynamicParticle.hh" |
---|
| 67 | #include "G4VParticleChange.hh" |
---|
| 68 | #include "G4ThreeVector.hh" |
---|
| 69 | #include "G4EnergyLossTables.hh" |
---|
| 70 | #include "G4VCrossSectionHandler.hh" |
---|
| 71 | #include "G4CrossSectionHandler.hh" |
---|
| 72 | #include "G4VEMDataSet.hh" |
---|
| 73 | #include "G4EMDataSet.hh" |
---|
| 74 | #include "G4CompositeEMDataSet.hh" |
---|
| 75 | #include "G4VDataSetAlgorithm.hh" |
---|
| 76 | #include "G4LogLogInterpolation.hh" |
---|
| 77 | #include "G4VRangeTest.hh" |
---|
| 78 | #include "G4RangeTest.hh" |
---|
| 79 | #include "G4ProductionCutsTable.hh" |
---|
| 80 | #include "G4AtomicTransitionManager.hh" |
---|
| 81 | #include "G4AtomicShell.hh" |
---|
| 82 | #include "G4AtomicDeexcitation.hh" |
---|
| 83 | #include "G4PenelopeIntegrator.hh" |
---|
| 84 | #include "G4MaterialCutsCouple.hh" |
---|
| 85 | |
---|
| 86 | |
---|
| 87 | G4PenelopeCompton::G4PenelopeCompton(const G4String& processName) |
---|
| 88 | : G4VDiscreteProcess(processName), |
---|
| 89 | lowEnergyLimit(250*eV), |
---|
| 90 | highEnergyLimit(100*GeV), |
---|
| 91 | intrinsicLowEnergyLimit(10*eV), |
---|
| 92 | intrinsicHighEnergyLimit(100*GeV), |
---|
| 93 | energyForIntegration(0.0), |
---|
| 94 | ZForIntegration(1), |
---|
| 95 | nBins(200), |
---|
[961] | 96 | cutForLowEnergySecondaryPhotons(250.0*eV), |
---|
| 97 | fUseAtomicDeexcitation(true) |
---|
[819] | 98 | { |
---|
| 99 | if (lowEnergyLimit < intrinsicLowEnergyLimit || |
---|
| 100 | highEnergyLimit > intrinsicHighEnergyLimit) |
---|
| 101 | { |
---|
| 102 | G4Exception("G4PenelopeCompton::G4PenelopeCompton - energy outside intrinsic process validity range"); |
---|
| 103 | } |
---|
| 104 | |
---|
| 105 | meanFreePathTable = 0; |
---|
| 106 | ionizationEnergy = new std::map<G4int,G4DataVector*>; |
---|
| 107 | hartreeFunction = new std::map<G4int,G4DataVector*>; |
---|
| 108 | occupationNumber = new std::map<G4int,G4DataVector*>; |
---|
| 109 | |
---|
| 110 | rangeTest = new G4RangeTest; |
---|
| 111 | |
---|
| 112 | ReadData(); //Read data from file |
---|
| 113 | |
---|
| 114 | if (verboseLevel > 0) |
---|
| 115 | { |
---|
| 116 | G4cout << GetProcessName() << " is created " << G4endl |
---|
| 117 | << "Energy range: " |
---|
| 118 | << lowEnergyLimit / keV << " keV - " |
---|
| 119 | << highEnergyLimit / GeV << " GeV" |
---|
| 120 | << G4endl; |
---|
| 121 | } |
---|
| 122 | } |
---|
| 123 | |
---|
| 124 | G4PenelopeCompton::~G4PenelopeCompton() |
---|
| 125 | { |
---|
| 126 | delete meanFreePathTable; |
---|
| 127 | delete rangeTest; |
---|
| 128 | |
---|
[961] | 129 | for (size_t i=0;i<matCrossSections->size();i++) |
---|
[819] | 130 | { |
---|
[961] | 131 | delete (*matCrossSections)[i]; |
---|
[819] | 132 | } |
---|
| 133 | |
---|
| 134 | delete matCrossSections; |
---|
| 135 | |
---|
| 136 | for (G4int Z=1;Z<100;Z++) |
---|
| 137 | { |
---|
| 138 | if (ionizationEnergy->count(Z)) delete (ionizationEnergy->find(Z)->second); |
---|
| 139 | if (hartreeFunction->count(Z)) delete (hartreeFunction->find(Z)->second); |
---|
| 140 | if (occupationNumber->count(Z)) delete (occupationNumber->find(Z)->second); |
---|
| 141 | } |
---|
| 142 | delete ionizationEnergy; |
---|
| 143 | delete hartreeFunction; |
---|
| 144 | delete occupationNumber; |
---|
| 145 | } |
---|
| 146 | |
---|
| 147 | void G4PenelopeCompton::BuildPhysicsTable(const G4ParticleDefinition& ) |
---|
| 148 | { |
---|
| 149 | G4DataVector energyVector; |
---|
| 150 | G4double dBin = std::log10(highEnergyLimit/lowEnergyLimit)/nBins; |
---|
[961] | 151 | for (G4int i=0;i<nBins;i++) |
---|
[819] | 152 | { |
---|
| 153 | energyVector.push_back(std::pow(10.,std::log10(lowEnergyLimit)+i*dBin)); |
---|
| 154 | } |
---|
| 155 | |
---|
| 156 | const G4MaterialTable* materialTable = G4Material::GetMaterialTable(); |
---|
| 157 | G4int nMaterials = G4Material::GetNumberOfMaterials(); |
---|
| 158 | G4VDataSetAlgorithm* algo = new G4LogLogInterpolation(); |
---|
| 159 | |
---|
[961] | 160 | //size_t nOfBins = energyVector.size(); |
---|
| 161 | //size_t bin=0; |
---|
[819] | 162 | |
---|
| 163 | G4DataVector* energies; |
---|
| 164 | G4DataVector* data; |
---|
| 165 | |
---|
| 166 | matCrossSections = new std::vector<G4VEMDataSet*>; |
---|
| 167 | |
---|
[961] | 168 | for (G4int m=0; m<nMaterials; m++) |
---|
[819] | 169 | { |
---|
| 170 | const G4Material* material= (*materialTable)[m]; |
---|
| 171 | G4int nElements = material->GetNumberOfElements(); |
---|
| 172 | const G4ElementVector* elementVector = material->GetElementVector(); |
---|
| 173 | const G4double* nAtomsPerVolume = material->GetAtomicNumDensityVector(); |
---|
| 174 | |
---|
| 175 | G4VEMDataSet* setForMat = new G4CompositeEMDataSet(algo,1.,1.); |
---|
| 176 | |
---|
[961] | 177 | for (G4int i=0; i<nElements; i++) { |
---|
[819] | 178 | |
---|
| 179 | G4int Z = (G4int) (*elementVector)[i]->GetZ(); |
---|
| 180 | G4double density = nAtomsPerVolume[i]; |
---|
| 181 | G4double cross=0.0; |
---|
| 182 | energies = new G4DataVector; |
---|
| 183 | data = new G4DataVector; |
---|
| 184 | |
---|
| 185 | |
---|
[961] | 186 | for (size_t bin=0; bin<energyVector.size(); bin++) |
---|
[819] | 187 | { |
---|
| 188 | G4double e = energyVector[bin]; |
---|
| 189 | energies->push_back(e); |
---|
| 190 | cross = density * CrossSection(e,Z); |
---|
| 191 | data->push_back(cross); |
---|
| 192 | } |
---|
| 193 | |
---|
| 194 | G4VEMDataSet* elSet = new G4EMDataSet(i,energies,data,algo,1.,1.); |
---|
| 195 | setForMat->AddComponent(elSet); |
---|
| 196 | } |
---|
| 197 | |
---|
| 198 | matCrossSections->push_back(setForMat); |
---|
| 199 | } |
---|
| 200 | |
---|
| 201 | |
---|
| 202 | //Build the mean free path table! |
---|
| 203 | G4double matCS = 0.0; |
---|
| 204 | G4VEMDataSet* matCrossSet = new G4CompositeEMDataSet(algo,1.,1.); |
---|
| 205 | G4VEMDataSet* materialSet = new G4CompositeEMDataSet(algo,1.,1.); |
---|
| 206 | |
---|
| 207 | |
---|
[961] | 208 | for (G4int m=0; m<nMaterials; m++) |
---|
[819] | 209 | { |
---|
| 210 | energies = new G4DataVector; |
---|
| 211 | data = new G4DataVector; |
---|
| 212 | const G4Material* material= (*materialTable)[m]; |
---|
| 213 | material= (*materialTable)[m]; |
---|
[961] | 214 | for (size_t bin=0; bin<energyVector.size(); bin++) |
---|
[819] | 215 | { |
---|
| 216 | G4double energy = energyVector[bin]; |
---|
| 217 | energies->push_back(energy); |
---|
| 218 | matCrossSet = (*matCrossSections)[m]; |
---|
| 219 | matCS = 0.0; |
---|
| 220 | G4int nElm = matCrossSet->NumberOfComponents(); |
---|
| 221 | for(G4int j=0; j<nElm; j++) { |
---|
| 222 | matCS += matCrossSet->GetComponent(j)->FindValue(energy); |
---|
| 223 | } |
---|
| 224 | if (matCS > 0.) |
---|
| 225 | { |
---|
| 226 | data->push_back(1./matCS); |
---|
| 227 | } |
---|
| 228 | else |
---|
| 229 | { |
---|
| 230 | data->push_back(DBL_MAX); |
---|
| 231 | } |
---|
| 232 | } |
---|
| 233 | G4VEMDataSet* dataSet = new G4EMDataSet(m,energies,data,algo,1.,1.); |
---|
| 234 | materialSet->AddComponent(dataSet); |
---|
| 235 | } |
---|
| 236 | meanFreePathTable = materialSet; |
---|
| 237 | } |
---|
| 238 | |
---|
| 239 | G4VParticleChange* G4PenelopeCompton::PostStepDoIt(const G4Track& aTrack, |
---|
| 240 | const G4Step& aStep) |
---|
| 241 | { |
---|
| 242 | //Penelope model |
---|
| 243 | |
---|
| 244 | aParticleChange.Initialize(aTrack); |
---|
| 245 | |
---|
| 246 | // Dynamic particle quantities |
---|
| 247 | const G4DynamicParticle* incidentPhoton = aTrack.GetDynamicParticle(); |
---|
| 248 | G4double photonEnergy0 = incidentPhoton->GetKineticEnergy(); |
---|
| 249 | |
---|
| 250 | if (photonEnergy0 <= lowEnergyLimit) |
---|
| 251 | { |
---|
| 252 | aParticleChange.ProposeTrackStatus(fStopAndKill); |
---|
| 253 | aParticleChange.ProposeEnergy(0.); |
---|
| 254 | aParticleChange.ProposeLocalEnergyDeposit(photonEnergy0); |
---|
| 255 | return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep); |
---|
| 256 | } |
---|
| 257 | |
---|
| 258 | G4ParticleMomentum photonDirection0 = incidentPhoton->GetMomentumDirection(); |
---|
| 259 | |
---|
| 260 | const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple(); |
---|
| 261 | const G4Material* material = couple->GetMaterial(); |
---|
[961] | 262 | |
---|
[819] | 263 | G4int Z = SelectRandomAtomForCompton(material,photonEnergy0); |
---|
| 264 | const G4int nmax = 64; |
---|
| 265 | G4double rn[nmax],pac[nmax]; |
---|
| 266 | |
---|
| 267 | G4double ki,ki1,ki2,ki3,taumin,a1,a2; |
---|
| 268 | G4double tau,TST; |
---|
| 269 | G4double S=0.0; |
---|
| 270 | G4double epsilon,cosTheta; |
---|
| 271 | G4double harFunc = 0.0; |
---|
| 272 | G4int occupNb= 0; |
---|
| 273 | G4double ionEnergy=0.0; |
---|
| 274 | G4int nosc = occupationNumber->find(Z)->second->size(); |
---|
| 275 | G4int iosc = nosc; |
---|
| 276 | ki = photonEnergy0/electron_mass_c2; |
---|
| 277 | ki2 = 2*ki+1.0; |
---|
| 278 | ki3 = ki*ki; |
---|
| 279 | ki1 = ki3-ki2-1.0; |
---|
| 280 | taumin = 1.0/ki2; |
---|
| 281 | a1 = std::log(ki2); |
---|
| 282 | a2 = a1+2.0*ki*(1.0+ki)/(ki2*ki2); |
---|
[961] | 283 | //If the incoming photon is above 5 MeV, the quicker approach based on the |
---|
| 284 | //pure Klein-Nishina formula is used |
---|
[819] | 285 | if (photonEnergy0 > 5*MeV) |
---|
| 286 | { |
---|
| 287 | do{ |
---|
| 288 | do{ |
---|
| 289 | if ((a2*G4UniformRand()) < a1) |
---|
| 290 | { |
---|
| 291 | tau = std::pow(taumin,G4UniformRand()); |
---|
| 292 | } |
---|
| 293 | else |
---|
| 294 | { |
---|
| 295 | tau = std::sqrt(1.0+G4UniformRand()*(taumin*taumin-1.0)); |
---|
| 296 | } |
---|
| 297 | //rejection function |
---|
| 298 | TST = (1+tau*(ki1+tau*(ki2+tau*ki3)))/(ki3*tau*(1.0+tau*tau)); |
---|
| 299 | }while (G4UniformRand()> TST); |
---|
| 300 | epsilon=tau; |
---|
| 301 | cosTheta = 1.0 - (1.0-tau)/(ki*tau); |
---|
| 302 | //Target shell electrons |
---|
| 303 | TST = Z*G4UniformRand(); |
---|
| 304 | iosc = nosc; |
---|
| 305 | S=0.0; |
---|
| 306 | for (G4int j=0;j<nosc;j++) |
---|
| 307 | { |
---|
| 308 | occupNb = (G4int) (*(occupationNumber->find(Z)->second))[j]; |
---|
| 309 | S = S + occupNb; |
---|
| 310 | if (S > TST) iosc = j; |
---|
| 311 | if (S > TST) break; |
---|
| 312 | } |
---|
| 313 | ionEnergy = (*(ionizationEnergy->find(Z)->second))[iosc]; |
---|
| 314 | }while((epsilon*photonEnergy0-photonEnergy0+ionEnergy) >0); |
---|
| 315 | } |
---|
| 316 | else //photonEnergy0<5 MeV |
---|
| 317 | { |
---|
| 318 | //Incoherent scattering function for theta=PI |
---|
| 319 | G4double s0=0.0; |
---|
| 320 | G4double pzomc=0.0,rni=0.0; |
---|
| 321 | G4double aux=0.0; |
---|
| 322 | for (G4int i=0;i<nosc;i++){ |
---|
| 323 | ionEnergy = (*(ionizationEnergy->find(Z)->second))[i]; |
---|
| 324 | if (photonEnergy0 > ionEnergy) |
---|
| 325 | { |
---|
| 326 | G4double aux = photonEnergy0*(photonEnergy0-ionEnergy)*2.0; |
---|
| 327 | harFunc = (*(hartreeFunction->find(Z)->second))[i]/fine_structure_const; |
---|
| 328 | occupNb = (G4int) (*(occupationNumber->find(Z)->second))[i]; |
---|
| 329 | pzomc = harFunc*(aux-electron_mass_c2*ionEnergy)/ |
---|
| 330 | (electron_mass_c2*std::sqrt(2.0*aux+ionEnergy*ionEnergy)); |
---|
| 331 | if (pzomc > 0) |
---|
| 332 | { |
---|
| 333 | rni = 1.0-0.5*std::exp(0.5-(std::sqrt(0.5)+std::sqrt(2.0)*pzomc)*(std::sqrt(0.5)+std::sqrt(2.0)*pzomc)); |
---|
| 334 | } |
---|
| 335 | else |
---|
| 336 | { |
---|
| 337 | rni = 0.5*std::exp(0.5-(std::sqrt(0.5)-std::sqrt(2.0)*pzomc)*(std::sqrt(0.5)-std::sqrt(2.0)*pzomc)); |
---|
| 338 | } |
---|
| 339 | s0 = s0 + occupNb*rni; |
---|
| 340 | } |
---|
| 341 | } |
---|
| 342 | |
---|
| 343 | //Sampling tau |
---|
| 344 | G4double cdt1; |
---|
| 345 | do |
---|
| 346 | { |
---|
| 347 | if ((G4UniformRand()*a2) < a1) |
---|
| 348 | { |
---|
| 349 | tau = std::pow(taumin,G4UniformRand()); |
---|
| 350 | } |
---|
| 351 | else |
---|
| 352 | { |
---|
| 353 | tau = std::sqrt(1.0+G4UniformRand()*(taumin*taumin-1.0)); |
---|
| 354 | } |
---|
| 355 | cdt1 = (1.0-tau)/(ki*tau); |
---|
| 356 | S=0.0; |
---|
| 357 | //Incoherent scattering function |
---|
| 358 | for (G4int i=0;i<nosc;i++){ |
---|
| 359 | ionEnergy = (*(ionizationEnergy->find(Z)->second))[i]; |
---|
| 360 | if (photonEnergy0 > ionEnergy) //sum only on excitable levels |
---|
| 361 | { |
---|
| 362 | aux = photonEnergy0*(photonEnergy0-ionEnergy)*cdt1; |
---|
| 363 | harFunc = (*(hartreeFunction->find(Z)->second))[i]/fine_structure_const; |
---|
| 364 | occupNb = (G4int) (*(occupationNumber->find(Z)->second))[i]; |
---|
| 365 | pzomc = harFunc*(aux-electron_mass_c2*ionEnergy)/ |
---|
| 366 | (electron_mass_c2*std::sqrt(2.0*aux+ionEnergy*ionEnergy)); |
---|
| 367 | if (pzomc > 0) |
---|
| 368 | { |
---|
[961] | 369 | rn[i] = 1.0-0.5*std::exp(0.5-(std::sqrt(0.5)+std::sqrt(2.0)*pzomc)* |
---|
| 370 | (std::sqrt(0.5)+std::sqrt(2.0)*pzomc)); |
---|
[819] | 371 | } |
---|
| 372 | else |
---|
| 373 | { |
---|
[961] | 374 | rn[i] = 0.5*std::exp(0.5-(std::sqrt(0.5)-std::sqrt(2.0)*pzomc)* |
---|
| 375 | (std::sqrt(0.5)-std::sqrt(2.0)*pzomc)); |
---|
[819] | 376 | } |
---|
| 377 | S = S + occupNb*rn[i]; |
---|
| 378 | pac[i] = S; |
---|
| 379 | } |
---|
| 380 | else |
---|
| 381 | { |
---|
| 382 | pac[i] = S-(1e-06); |
---|
| 383 | } |
---|
| 384 | } |
---|
| 385 | //Rejection function |
---|
| 386 | TST = S*(1.0+tau*(ki1+tau*(ki2+tau*ki3)))/(ki3*tau*(1.0+tau*tau)); |
---|
| 387 | }while ((G4UniformRand()*s0) > TST); |
---|
| 388 | |
---|
| 389 | //Target electron shell |
---|
| 390 | cosTheta = 1.0 - cdt1; |
---|
| 391 | G4double fpzmax=0.0,fpz=0.0; |
---|
| 392 | G4double A=0.0; |
---|
| 393 | do |
---|
| 394 | { |
---|
| 395 | do |
---|
| 396 | { |
---|
| 397 | TST =S*G4UniformRand(); |
---|
| 398 | iosc=nosc; |
---|
| 399 | for (G4int i=0;i<nosc;i++){ |
---|
| 400 | if (pac[i]>TST) iosc = i; |
---|
| 401 | if (pac[i]>TST) break; |
---|
| 402 | } |
---|
| 403 | A = G4UniformRand()*rn[iosc]; |
---|
| 404 | harFunc = (*(hartreeFunction->find(Z)->second))[iosc]/fine_structure_const; |
---|
| 405 | occupNb = (G4int) (*(occupationNumber->find(Z)->second))[iosc]; |
---|
| 406 | if (A < 0.5) { |
---|
| 407 | pzomc = (std::sqrt(0.5)-std::sqrt(0.5-std::log(2.0*A)))/ |
---|
| 408 | (std::sqrt(2.0)*harFunc); |
---|
| 409 | } |
---|
| 410 | else |
---|
| 411 | { |
---|
| 412 | pzomc = (std::sqrt(0.5-std::log(2.0-2.0*A))-std::sqrt(0.5))/ |
---|
| 413 | (std::sqrt(2.0)*harFunc); |
---|
| 414 | } |
---|
| 415 | } while (pzomc < -1); |
---|
| 416 | // F(EP) rejection |
---|
| 417 | G4double XQC = 1.0+tau*(tau-2.0*cosTheta); |
---|
| 418 | G4double AF = std::sqrt(XQC)*(1.0+tau*(tau-cosTheta)/XQC); |
---|
| 419 | if (AF > 0) { |
---|
| 420 | fpzmax = 1.0+AF*0.2; |
---|
| 421 | } |
---|
| 422 | else |
---|
| 423 | { |
---|
| 424 | fpzmax = 1.0-AF*0.2; |
---|
| 425 | } |
---|
| 426 | fpz = 1.0+AF*std::max(std::min(pzomc,0.2),-0.2); |
---|
| 427 | }while ((fpzmax*G4UniformRand())>fpz); |
---|
| 428 | |
---|
| 429 | //Energy of the scattered photon |
---|
| 430 | G4double T = pzomc*pzomc; |
---|
| 431 | G4double b1 = 1.0-T*tau*tau; |
---|
| 432 | G4double b2 = 1.0-T*tau*cosTheta; |
---|
| 433 | if (pzomc > 0.0) |
---|
| 434 | { |
---|
| 435 | epsilon = (tau/b1)*(b2+std::sqrt(std::abs(b2*b2-b1*(1.0-T)))); |
---|
| 436 | } |
---|
| 437 | else |
---|
| 438 | { |
---|
| 439 | epsilon = (tau/b1)*(b2-std::sqrt(std::abs(b2*b2-b1*(1.0-T)))); |
---|
| 440 | } |
---|
| 441 | } |
---|
| 442 | |
---|
| 443 | G4double sinTheta = std::sqrt(1-cosTheta*cosTheta); |
---|
| 444 | G4double phi = twopi * G4UniformRand() ; |
---|
| 445 | G4double dirx = sinTheta * std::cos(phi); |
---|
| 446 | G4double diry = sinTheta * std::sin(phi); |
---|
| 447 | G4double dirz = cosTheta ; |
---|
| 448 | |
---|
| 449 | // Update G4VParticleChange for the scattered photon |
---|
| 450 | |
---|
| 451 | G4ThreeVector photonDirection1(dirx,diry,dirz); |
---|
| 452 | photonDirection1.rotateUz(photonDirection0); |
---|
| 453 | aParticleChange.ProposeMomentumDirection(photonDirection1) ; |
---|
| 454 | G4double photonEnergy1 = epsilon * photonEnergy0; |
---|
| 455 | |
---|
| 456 | if (photonEnergy1 > 0.) |
---|
| 457 | { |
---|
| 458 | aParticleChange.ProposeEnergy(photonEnergy1) ; |
---|
| 459 | } |
---|
| 460 | else |
---|
| 461 | { |
---|
| 462 | aParticleChange.ProposeEnergy(0.) ; |
---|
| 463 | aParticleChange.ProposeTrackStatus(fStopAndKill); |
---|
| 464 | } |
---|
| 465 | |
---|
| 466 | |
---|
[961] | 467 | // Kinematics of the scattered electron |
---|
[819] | 468 | G4double diffEnergy = photonEnergy0*(1-epsilon); |
---|
| 469 | ionEnergy = (*(ionizationEnergy->find(Z)->second))[iosc]; |
---|
| 470 | G4double Q2 = photonEnergy0*photonEnergy0+photonEnergy1*(photonEnergy1-2.0*photonEnergy0*cosTheta); |
---|
| 471 | G4double cosThetaE; //scattering angle for the electron |
---|
| 472 | if (Q2 > 1.0e-12) |
---|
| 473 | { |
---|
| 474 | cosThetaE = (photonEnergy0-photonEnergy1*cosTheta)/std::sqrt(Q2); |
---|
| 475 | } |
---|
| 476 | else |
---|
| 477 | { |
---|
| 478 | cosThetaE = 1.0; |
---|
| 479 | } |
---|
| 480 | G4double sinThetaE = std::sqrt(1-cosThetaE*cosThetaE); |
---|
[961] | 481 | //initialize here, then check photons created by Atomic-Deexcitation, and the final state e- |
---|
| 482 | G4int nbOfSecondaries = 0; |
---|
| 483 | |
---|
| 484 | std::vector<G4DynamicParticle*>* photonVector=0; |
---|
[819] | 485 | |
---|
| 486 | const G4AtomicTransitionManager* transitionManager = G4AtomicTransitionManager::Instance(); |
---|
| 487 | const G4AtomicShell* shell = transitionManager->Shell(Z,iosc); |
---|
| 488 | G4double bindingEnergy = shell->BindingEnergy(); |
---|
| 489 | G4int shellId = shell->ShellId(); |
---|
[961] | 490 | G4double ionEnergyInPenelopeDatabase = ionEnergy; |
---|
| 491 | ionEnergy = std::max(bindingEnergy,ionEnergyInPenelopeDatabase); //protection against energy non-conservation |
---|
[819] | 492 | |
---|
[961] | 493 | G4double eKineticEnergy = diffEnergy - ionEnergy; //subtract the excitation energy. If not emitted by fluorescence, |
---|
| 494 | //the ionization energy is deposited as local energy deposition |
---|
| 495 | G4double localEnergyDeposit = ionEnergy; |
---|
| 496 | G4double energyInFluorescence = 0.; //testing purposes only |
---|
[819] | 497 | |
---|
[961] | 498 | if (eKineticEnergy < 0) |
---|
| 499 | { |
---|
| 500 | //It means that there was some problem/mismatch between the two databases. Try to make it work |
---|
| 501 | //In this case available Energy (diffEnergy) < ionEnergy |
---|
| 502 | //Full residual energy is deposited locally |
---|
| 503 | localEnergyDeposit = diffEnergy; |
---|
| 504 | eKineticEnergy = 0.0; |
---|
| 505 | } |
---|
[819] | 506 | |
---|
[961] | 507 | //the local energy deposit is what remains: part of this may be spent for fluorescence. |
---|
| 508 | |
---|
| 509 | if (fUseAtomicDeexcitation) |
---|
| 510 | { |
---|
| 511 | G4int nPhotons=0; |
---|
| 512 | |
---|
| 513 | const G4ProductionCutsTable* theCoupleTable= |
---|
| 514 | G4ProductionCutsTable::GetProductionCutsTable(); |
---|
| 515 | size_t indx = couple->GetIndex(); |
---|
[819] | 516 | |
---|
[961] | 517 | G4double cutg = (*(theCoupleTable->GetEnergyCutsVector(0)))[indx]; |
---|
| 518 | cutg = std::max(cutForLowEnergySecondaryPhotons,cutg); |
---|
| 519 | |
---|
| 520 | G4double cute = (*(theCoupleTable->GetEnergyCutsVector(1)))[indx]; |
---|
| 521 | cute = std::max(cutForLowEnergySecondaryPhotons,cute); |
---|
| 522 | |
---|
| 523 | G4DynamicParticle* aPhoton; |
---|
| 524 | G4AtomicDeexcitation deexcitationManager; |
---|
| 525 | |
---|
| 526 | if (Z>5 && (localEnergyDeposit > cutg || localEnergyDeposit > cute)) |
---|
| 527 | { |
---|
| 528 | photonVector = deexcitationManager.GenerateParticles(Z,shellId); |
---|
| 529 | for (size_t k=0;k<photonVector->size();k++){ |
---|
| 530 | aPhoton = (*photonVector)[k]; |
---|
| 531 | if (aPhoton) |
---|
[819] | 532 | { |
---|
[961] | 533 | G4double itsCut = cutg; |
---|
| 534 | if (aPhoton->GetDefinition() == G4Electron::Electron()) itsCut = cute; |
---|
| 535 | G4double itsEnergy = aPhoton->GetKineticEnergy(); |
---|
| 536 | if (itsEnergy > itsCut && itsEnergy <= localEnergyDeposit) |
---|
| 537 | { |
---|
| 538 | nPhotons++; |
---|
| 539 | localEnergyDeposit -= itsEnergy; |
---|
| 540 | energyInFluorescence += itsEnergy; |
---|
| 541 | } |
---|
| 542 | else |
---|
| 543 | { |
---|
| 544 | delete aPhoton; |
---|
| 545 | (*photonVector)[k]=0; |
---|
| 546 | } |
---|
[819] | 547 | } |
---|
| 548 | } |
---|
[961] | 549 | } |
---|
| 550 | nbOfSecondaries=nPhotons; |
---|
[819] | 551 | } |
---|
| 552 | |
---|
[961] | 553 | |
---|
[819] | 554 | // Generate the electron only if with large enough range w.r.t. cuts and safety |
---|
| 555 | G4double safety = aStep.GetPostStepPoint()->GetSafety(); |
---|
| 556 | G4DynamicParticle* electron = 0; |
---|
[961] | 557 | if (rangeTest->Escape(G4Electron::Electron(),couple,eKineticEnergy,safety) && |
---|
| 558 | eKineticEnergy>cutForLowEnergySecondaryPhotons) |
---|
[819] | 559 | { |
---|
| 560 | G4double xEl = sinThetaE * std::cos(phi+pi); |
---|
| 561 | G4double yEl = sinThetaE * std::sin(phi+pi); |
---|
| 562 | G4double zEl = cosThetaE; |
---|
| 563 | G4ThreeVector eDirection(xEl,yEl,zEl); //electron direction |
---|
| 564 | eDirection.rotateUz(photonDirection0); |
---|
| 565 | electron = new G4DynamicParticle (G4Electron::Electron(), |
---|
| 566 | eDirection,eKineticEnergy) ; |
---|
| 567 | nbOfSecondaries++; |
---|
| 568 | } |
---|
| 569 | else |
---|
| 570 | { |
---|
[961] | 571 | localEnergyDeposit += eKineticEnergy; |
---|
[819] | 572 | } |
---|
| 573 | |
---|
| 574 | aParticleChange.SetNumberOfSecondaries(nbOfSecondaries); |
---|
| 575 | if (electron) aParticleChange.AddSecondary(electron); |
---|
[961] | 576 | |
---|
| 577 | //This block below is executed only if there is at least one secondary photon produced by |
---|
| 578 | //AtomicDeexcitation |
---|
| 579 | if (photonVector) |
---|
[819] | 580 | { |
---|
[961] | 581 | for (size_t ll=0;ll<photonVector->size();ll++) |
---|
| 582 | { |
---|
| 583 | if ((*photonVector)[ll]) aParticleChange.AddSecondary((*photonVector)[ll]); |
---|
| 584 | } |
---|
[819] | 585 | } |
---|
| 586 | delete photonVector; |
---|
[961] | 587 | if (localEnergyDeposit < 0) |
---|
[819] | 588 | { |
---|
| 589 | G4cout << "WARNING-" |
---|
| 590 | << "G4PenelopeCompton::PostStepDoIt - Negative energy deposit" |
---|
| 591 | << G4endl; |
---|
[961] | 592 | localEnergyDeposit=0.; |
---|
[819] | 593 | } |
---|
[961] | 594 | aParticleChange.ProposeLocalEnergyDeposit(localEnergyDeposit); |
---|
[819] | 595 | |
---|
| 596 | |
---|
[961] | 597 | if (verboseLevel > 1) |
---|
| 598 | { |
---|
| 599 | G4cout << "-----------------------------------------------------------" << G4endl; |
---|
| 600 | G4cout << "Energy balance from G4PenelopeCompton" << G4endl; |
---|
| 601 | G4cout << "Incoming photon energy: " << photonEnergy0/keV << " keV" << G4endl; |
---|
| 602 | G4cout << "-----------------------------------------------------------" << G4endl; |
---|
| 603 | G4cout << "Scattered photon: " << photonEnergy1/keV << " keV" << G4endl; |
---|
| 604 | G4double electronEnergy = 0.; |
---|
| 605 | if (electron) |
---|
| 606 | electronEnergy = eKineticEnergy; |
---|
| 607 | G4cout << "Scattered electron " << electronEnergy/keV << " keV" << G4endl; |
---|
| 608 | G4cout << "Fluorescence: " << energyInFluorescence/keV << " keV" << G4endl; |
---|
| 609 | G4cout << "Local energy deposit " << localEnergyDeposit/keV << " keV" << G4endl; |
---|
| 610 | G4cout << "Total final state: " << (photonEnergy1+electronEnergy+energyInFluorescence+localEnergyDeposit)/keV << |
---|
| 611 | " keV" << G4endl; |
---|
| 612 | G4cout << "-----------------------------------------------------------" << G4endl; |
---|
| 613 | } |
---|
| 614 | |
---|
[819] | 615 | return G4VDiscreteProcess::PostStepDoIt( aTrack, aStep); |
---|
| 616 | } |
---|
| 617 | |
---|
| 618 | G4bool G4PenelopeCompton::IsApplicable(const G4ParticleDefinition& particle) |
---|
| 619 | { |
---|
| 620 | return ( &particle == G4Gamma::Gamma() ); |
---|
| 621 | } |
---|
| 622 | |
---|
| 623 | G4double G4PenelopeCompton::GetMeanFreePath(const G4Track& track, |
---|
| 624 | G4double, // previousStepSize |
---|
| 625 | G4ForceCondition*) |
---|
| 626 | { |
---|
| 627 | const G4DynamicParticle* photon = track.GetDynamicParticle(); |
---|
| 628 | G4double energy = photon->GetKineticEnergy(); |
---|
| 629 | G4Material* material = track.GetMaterial(); |
---|
| 630 | size_t materialIndex = material->GetIndex(); |
---|
| 631 | |
---|
| 632 | G4double meanFreePath; |
---|
| 633 | if (energy > highEnergyLimit) meanFreePath = meanFreePathTable->FindValue(highEnergyLimit,materialIndex); |
---|
| 634 | else if (energy < lowEnergyLimit) meanFreePath = DBL_MAX; |
---|
| 635 | else meanFreePath = meanFreePathTable->FindValue(energy,materialIndex); |
---|
| 636 | return meanFreePath; |
---|
| 637 | } |
---|
| 638 | |
---|
| 639 | |
---|
| 640 | void G4PenelopeCompton::ReadData() |
---|
| 641 | { |
---|
| 642 | char* path = getenv("G4LEDATA"); |
---|
| 643 | if (!path) |
---|
| 644 | { |
---|
| 645 | G4String excep = "G4PenelopeCompton - G4LEDATA environment variable not set!"; |
---|
| 646 | G4Exception(excep); |
---|
| 647 | } |
---|
| 648 | G4String pathString(path); |
---|
| 649 | G4String pathFile = pathString + "/penelope/compton-pen.dat"; |
---|
| 650 | std::ifstream file(pathFile); |
---|
| 651 | std::filebuf* lsdp = file.rdbuf(); |
---|
| 652 | |
---|
| 653 | if (!(lsdp->is_open())) |
---|
| 654 | { |
---|
| 655 | G4String excep = "G4PenelopeCompton - data file " + pathFile + " not found!"; |
---|
| 656 | G4Exception(excep); |
---|
| 657 | } |
---|
| 658 | |
---|
| 659 | G4int k1,test,test1; |
---|
| 660 | G4double a1,a2; |
---|
| 661 | G4int Z=1,nLevels=0; |
---|
| 662 | G4DataVector* f; |
---|
| 663 | G4DataVector* u; |
---|
| 664 | G4DataVector* j; |
---|
| 665 | |
---|
| 666 | do{ |
---|
| 667 | f = new G4DataVector; |
---|
| 668 | u = new G4DataVector; |
---|
| 669 | j = new G4DataVector; |
---|
| 670 | file >> Z >> nLevels; |
---|
| 671 | for (G4int h=0;h<nLevels;h++){ |
---|
| 672 | file >> k1 >> a1 >> a2; |
---|
| 673 | f->push_back((G4double) k1); |
---|
| 674 | u->push_back(a1); |
---|
| 675 | j->push_back(a2); |
---|
| 676 | } |
---|
| 677 | ionizationEnergy->insert(std::make_pair(Z,u)); |
---|
| 678 | hartreeFunction->insert(std::make_pair(Z,j)); |
---|
| 679 | occupationNumber->insert(std::make_pair(Z,f)); |
---|
| 680 | file >> test >> test1; //-1 -1 close the data for each Z |
---|
| 681 | if (test > 0) { |
---|
| 682 | G4String excep = "G4PenelopeCompton - data file corrupted!"; |
---|
| 683 | G4Exception(excep); |
---|
| 684 | } |
---|
| 685 | }while (test != -2); //the very last Z is closed with -2 instead of -1 |
---|
| 686 | } |
---|
| 687 | |
---|
| 688 | G4double G4PenelopeCompton::CrossSection(G4double energy,G4int Z) |
---|
| 689 | { |
---|
| 690 | G4double cs=0.0; |
---|
| 691 | energyForIntegration=energy; |
---|
| 692 | ZForIntegration = Z; |
---|
| 693 | if (energy< 5*MeV) |
---|
| 694 | { |
---|
| 695 | G4PenelopeIntegrator<G4PenelopeCompton,G4double (G4PenelopeCompton::*)(G4double)> theIntegrator; |
---|
| 696 | cs = theIntegrator.Calculate(this,&G4PenelopeCompton::DifferentialCrossSection,-1.0,1.0,1e-05); |
---|
| 697 | } |
---|
| 698 | else |
---|
| 699 | { |
---|
| 700 | G4double ki=energy/electron_mass_c2; |
---|
| 701 | G4double ki3=ki*ki; |
---|
| 702 | G4double ki2=1.0+2*ki; |
---|
| 703 | G4double ki1=ki3-ki2-1.0; |
---|
| 704 | G4double t0=1.0/(ki2); |
---|
| 705 | G4double csl = 0.5*ki3*t0*t0+ki2*t0+ki1*std::log(t0)-(1.0/t0); |
---|
| 706 | G4int nosc = occupationNumber->find(Z)->second->size(); |
---|
| 707 | for (G4int i=0;i<nosc;i++) |
---|
| 708 | { |
---|
| 709 | G4double ionEnergy = (*(ionizationEnergy->find(Z)->second))[i]; |
---|
| 710 | G4double tau=(energy-ionEnergy)/energy; |
---|
| 711 | if (tau > t0) |
---|
| 712 | { |
---|
| 713 | G4double csu = 0.5*ki3*tau*tau+ki2*tau+ki1*std::log(tau)-(1.0/tau); |
---|
| 714 | G4int f = (G4int) (*(occupationNumber->find(Z)->second))[i]; |
---|
| 715 | cs = cs + f*(csu-csl); |
---|
| 716 | } |
---|
| 717 | } |
---|
| 718 | cs=pi*classic_electr_radius*classic_electr_radius*cs/(ki*ki3); |
---|
| 719 | } |
---|
| 720 | return cs; |
---|
| 721 | } |
---|
| 722 | |
---|
| 723 | |
---|
| 724 | G4double G4PenelopeCompton::DifferentialCrossSection(G4double cosTheta) |
---|
| 725 | { |
---|
| 726 | const G4double k2 = std::sqrt(2.0); |
---|
| 727 | const G4double k1 = std::sqrt(0.5); |
---|
| 728 | const G4double k12 = 0.5; |
---|
| 729 | G4double cdt1 = 1.0-cosTheta; |
---|
| 730 | G4double energy = energyForIntegration; |
---|
| 731 | G4int Z = ZForIntegration; |
---|
| 732 | G4double ionEnergy=0.0,Pzimax=0.0,XKN=0.0; |
---|
| 733 | G4double diffCS=0.0; |
---|
| 734 | G4double x=0.0,siap=0.0; |
---|
| 735 | G4double harFunc=0.0; |
---|
| 736 | G4int occupNb; |
---|
| 737 | //energy of Compton line; |
---|
| 738 | G4double EOEC = 1.0+(energy/electron_mass_c2)*cdt1; |
---|
| 739 | G4double ECOE = 1.0/EOEC; |
---|
| 740 | //Incoherent scattering function (analytical profile) |
---|
| 741 | G4double sia = 0.0; |
---|
| 742 | G4int nosc = occupationNumber->find(Z)->second->size(); |
---|
| 743 | for (G4int i=0;i<nosc;i++){ |
---|
| 744 | ionEnergy = (*(ionizationEnergy->find(Z)->second))[i]; |
---|
| 745 | //Sum only of those shells for which E>Eion |
---|
| 746 | if (energy > ionEnergy) |
---|
| 747 | { |
---|
| 748 | G4double aux = energy * (energy-ionEnergy)*cdt1; |
---|
| 749 | Pzimax = (aux - electron_mass_c2*ionEnergy)/(electron_mass_c2*std::sqrt(2*aux+ionEnergy*ionEnergy)); |
---|
| 750 | harFunc = (*(hartreeFunction->find(Z)->second))[i]/fine_structure_const; |
---|
| 751 | occupNb = (G4int) (*(occupationNumber->find(Z)->second))[i]; |
---|
| 752 | x = harFunc*Pzimax; |
---|
| 753 | if (x > 0) |
---|
| 754 | { |
---|
| 755 | siap = 1.0-0.5*std::exp(k12-(k1+k2*x)*(k1+k2*x)); |
---|
| 756 | } |
---|
| 757 | else |
---|
| 758 | { |
---|
| 759 | siap = 0.5*std::exp(k12-(k1-k2*x)*(k1-k2*x)); |
---|
| 760 | } |
---|
| 761 | sia = sia + occupNb*siap; //sum of all contributions; |
---|
| 762 | } |
---|
| 763 | } |
---|
| 764 | XKN = EOEC+ECOE-1+cosTheta*cosTheta; |
---|
| 765 | diffCS = pi*classic_electr_radius*classic_electr_radius*ECOE*ECOE*XKN*sia; |
---|
| 766 | return diffCS; |
---|
| 767 | } |
---|
| 768 | |
---|
| 769 | G4int G4PenelopeCompton::SelectRandomAtomForCompton(const G4Material* material,G4double energy) const |
---|
| 770 | { |
---|
| 771 | G4int nElements = material->GetNumberOfElements(); |
---|
| 772 | //Special case: the material consists of one element |
---|
| 773 | if (nElements == 1) |
---|
| 774 | { |
---|
| 775 | G4int Z = (G4int) material->GetZ(); |
---|
| 776 | return Z; |
---|
| 777 | } |
---|
| 778 | |
---|
| 779 | //Composite material |
---|
| 780 | const G4ElementVector* elementVector = material->GetElementVector(); |
---|
| 781 | size_t materialIndex = material->GetIndex(); |
---|
| 782 | |
---|
| 783 | G4VEMDataSet* materialSet = (*matCrossSections)[materialIndex]; |
---|
| 784 | G4double materialCrossSection0 = 0.0; |
---|
| 785 | G4DataVector cross; |
---|
| 786 | cross.clear(); |
---|
| 787 | G4int i; |
---|
| 788 | for (i=0;i<nElements;i++) |
---|
| 789 | { |
---|
| 790 | G4double cr = (materialSet->GetComponent(i))->FindValue(energy); |
---|
| 791 | materialCrossSection0 += cr; |
---|
| 792 | cross.push_back(materialCrossSection0); //cumulative cross section |
---|
| 793 | } |
---|
| 794 | |
---|
| 795 | G4double random = G4UniformRand()*materialCrossSection0; |
---|
| 796 | for (i=0;i<nElements;i++) |
---|
| 797 | { |
---|
| 798 | if (random <= cross[i]) return (G4int) (*elementVector)[i]->GetZ(); |
---|
| 799 | } |
---|
| 800 | //It should never get here |
---|
| 801 | return 0; |
---|
| 802 | } |
---|
| 803 | |
---|