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