[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 | // |
---|
[1055] | 26 | // $Id: G4PEEffectModel.cc,v 1.8 2009/04/09 18:41:18 vnivanch Exp $ |
---|
| 27 | // GEANT4 tag $Name: geant4-09-03-beta-cand-01 $ |
---|
[819] | 28 | // |
---|
| 29 | // ------------------------------------------------------------------- |
---|
| 30 | // |
---|
| 31 | // GEANT4 Class file |
---|
| 32 | // |
---|
| 33 | // |
---|
| 34 | // File name: G4PEEffectModel |
---|
| 35 | // |
---|
| 36 | // Author: Vladimir Ivanchenko on base of Michel Maire code |
---|
| 37 | // |
---|
| 38 | // Creation date: 21.03.2005 |
---|
| 39 | // |
---|
| 40 | // Modifications: |
---|
| 41 | // |
---|
| 42 | // 04.12.05 : SetProposedKineticEnergy(0.) for the killed photon (mma) |
---|
[1055] | 43 | // 20.02.09 : Added initialisation of deexcitation flag and method |
---|
| 44 | // CrossSectionPerVolume instead of mfp (V.Ivanchenko) |
---|
[819] | 45 | // |
---|
| 46 | // Class Description: |
---|
| 47 | // |
---|
| 48 | // ------------------------------------------------------------------- |
---|
| 49 | // |
---|
| 50 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 51 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 52 | |
---|
| 53 | #include "G4PEEffectModel.hh" |
---|
| 54 | #include "G4Electron.hh" |
---|
| 55 | #include "G4Gamma.hh" |
---|
| 56 | #include "Randomize.hh" |
---|
| 57 | #include "G4DataVector.hh" |
---|
| 58 | #include "G4ParticleChangeForGamma.hh" |
---|
| 59 | |
---|
| 60 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 61 | |
---|
| 62 | using namespace std; |
---|
| 63 | |
---|
| 64 | G4PEEffectModel::G4PEEffectModel(const G4ParticleDefinition*, |
---|
| 65 | const G4String& nam) |
---|
| 66 | : G4VEmModel(nam),isInitialized(false) |
---|
| 67 | { |
---|
| 68 | theGamma = G4Gamma::Gamma(); |
---|
| 69 | theElectron = G4Electron::Electron(); |
---|
[1055] | 70 | fminimalEnergy = 1.0*eV; |
---|
[819] | 71 | } |
---|
| 72 | |
---|
| 73 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 74 | |
---|
| 75 | G4PEEffectModel::~G4PEEffectModel() |
---|
[1055] | 76 | {} |
---|
[819] | 77 | |
---|
| 78 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 79 | |
---|
| 80 | void G4PEEffectModel::Initialise(const G4ParticleDefinition*, |
---|
| 81 | const G4DataVector&) |
---|
| 82 | { |
---|
[1055] | 83 | // always false before the run |
---|
| 84 | SetDeexcitationFlag(false); |
---|
[819] | 85 | |
---|
[1055] | 86 | if (isInitialized) return; |
---|
| 87 | fParticleChange = GetParticleChangeForGamma(); |
---|
| 88 | isInitialized = true; |
---|
[819] | 89 | } |
---|
| 90 | |
---|
[1055] | 91 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... |
---|
| 92 | |
---|
| 93 | G4double G4PEEffectModel::ComputeCrossSectionPerAtom(const G4ParticleDefinition*, |
---|
| 94 | G4double energy, |
---|
| 95 | G4double Z, G4double, |
---|
| 96 | G4double, G4double) |
---|
| 97 | { |
---|
| 98 | G4double* SandiaCof = G4SandiaTable::GetSandiaCofPerAtom((G4int)Z, energy); |
---|
| 99 | |
---|
| 100 | G4double energy2 = energy*energy; |
---|
| 101 | G4double energy3 = energy*energy2; |
---|
| 102 | G4double energy4 = energy2*energy2; |
---|
| 103 | |
---|
| 104 | return SandiaCof[0]/energy + SandiaCof[1]/energy2 + |
---|
| 105 | SandiaCof[2]/energy3 + SandiaCof[3]/energy4; |
---|
| 106 | } |
---|
| 107 | |
---|
[819] | 108 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 109 | |
---|
[1055] | 110 | G4double G4PEEffectModel::CrossSectionPerVolume(const G4Material* material, |
---|
| 111 | const G4ParticleDefinition*, |
---|
| 112 | G4double energy, |
---|
| 113 | G4double, G4double) |
---|
| 114 | { |
---|
| 115 | G4double* SandiaCof = |
---|
| 116 | material->GetSandiaTable()->GetSandiaCofForMaterial(energy); |
---|
| 117 | |
---|
| 118 | G4double energy2 = energy*energy; |
---|
| 119 | G4double energy3 = energy*energy2; |
---|
| 120 | G4double energy4 = energy2*energy2; |
---|
| 121 | |
---|
| 122 | return SandiaCof[0]/energy + SandiaCof[1]/energy2 + |
---|
| 123 | SandiaCof[2]/energy3 + SandiaCof[3]/energy4; |
---|
| 124 | } |
---|
| 125 | |
---|
| 126 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 127 | |
---|
[819] | 128 | void G4PEEffectModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect, |
---|
| 129 | const G4MaterialCutsCouple* couple, |
---|
| 130 | const G4DynamicParticle* aDynamicPhoton, |
---|
| 131 | G4double, |
---|
| 132 | G4double) |
---|
| 133 | { |
---|
| 134 | const G4Material* aMaterial = couple->GetMaterial(); |
---|
| 135 | |
---|
| 136 | G4double energy = aDynamicPhoton->GetKineticEnergy(); |
---|
| 137 | G4ParticleMomentum PhotonDirection = aDynamicPhoton->GetMomentumDirection(); |
---|
| 138 | |
---|
| 139 | // select randomly one element constituing the material. |
---|
| 140 | const G4Element* anElement = SelectRandomAtom(aMaterial,theGamma,energy); |
---|
| 141 | |
---|
| 142 | // |
---|
| 143 | // Photo electron |
---|
| 144 | // |
---|
| 145 | |
---|
| 146 | // Select atomic shell |
---|
| 147 | G4int nShells = anElement->GetNbOfAtomicShells(); |
---|
| 148 | G4int i = 0; |
---|
| 149 | while ((i<nShells) && (energy<anElement->GetAtomicShell(i))) i++; |
---|
| 150 | |
---|
| 151 | // no shell available |
---|
| 152 | if (i == nShells) return; |
---|
| 153 | |
---|
| 154 | G4double bindingEnergy = anElement->GetAtomicShell(i); |
---|
| 155 | G4double ElecKineEnergy = energy - bindingEnergy; |
---|
| 156 | |
---|
| 157 | if (ElecKineEnergy > fminimalEnergy) |
---|
| 158 | { |
---|
| 159 | // direction of the photo electron |
---|
| 160 | // |
---|
| 161 | G4double cosTeta = ElecCosThetaDistribution(ElecKineEnergy); |
---|
| 162 | G4double sinTeta = sqrt(1.-cosTeta*cosTeta); |
---|
| 163 | G4double Phi = twopi * G4UniformRand(); |
---|
| 164 | G4double dirx = sinTeta*cos(Phi),diry = sinTeta*sin(Phi),dirz = cosTeta; |
---|
| 165 | G4ThreeVector ElecDirection(dirx,diry,dirz); |
---|
| 166 | ElecDirection.rotateUz(PhotonDirection); |
---|
| 167 | // |
---|
| 168 | G4DynamicParticle* aParticle = new G4DynamicParticle ( |
---|
| 169 | theElectron,ElecDirection, ElecKineEnergy); |
---|
| 170 | fvect->push_back(aParticle); |
---|
| 171 | } |
---|
| 172 | |
---|
| 173 | fParticleChange->SetProposedKineticEnergy(0.); |
---|
| 174 | fParticleChange->ProposeTrackStatus(fStopAndKill); |
---|
| 175 | fParticleChange->ProposeLocalEnergyDeposit(bindingEnergy); |
---|
| 176 | } |
---|
| 177 | |
---|
| 178 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 179 | |
---|
| 180 | G4double G4PEEffectModel::ElecCosThetaDistribution(G4double kineEnergy) |
---|
| 181 | { |
---|
| 182 | // Compute Theta distribution of the emitted electron, with respect to the |
---|
| 183 | // incident Gamma. |
---|
| 184 | // The Sauter-Gavrila distribution for the K-shell is used. |
---|
| 185 | // |
---|
| 186 | G4double costeta = 1.; |
---|
| 187 | G4double gamma = 1. + kineEnergy/electron_mass_c2; |
---|
| 188 | if (gamma > 5.) return costeta; |
---|
| 189 | G4double beta = sqrt(gamma*gamma-1.)/gamma; |
---|
| 190 | G4double b = 0.5*gamma*(gamma-1.)*(gamma-2); |
---|
| 191 | |
---|
| 192 | G4double rndm,term,greject,grejsup; |
---|
| 193 | if (gamma < 2.) grejsup = gamma*gamma*(1.+b-beta*b); |
---|
| 194 | else grejsup = gamma*gamma*(1.+b+beta*b); |
---|
| 195 | |
---|
| 196 | do { rndm = 1.-2*G4UniformRand(); |
---|
| 197 | costeta = (rndm+beta)/(rndm*beta+1.); |
---|
| 198 | term = 1.-beta*costeta; |
---|
| 199 | greject = (1.-costeta*costeta)*(1.+b*term)/(term*term); |
---|
| 200 | } while(greject < G4UniformRand()*grejsup); |
---|
| 201 | |
---|
| 202 | return costeta; |
---|
| 203 | } |
---|
| 204 | |
---|
| 205 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|