#include <G4PEEffectModel.hh>
Inheritance diagram for G4PEEffectModel:
Public Member Functions | |
G4PEEffectModel (const G4ParticleDefinition *p=0, const G4String &nam="PhotoElectric") | |
virtual | ~G4PEEffectModel () |
virtual void | Initialise (const G4ParticleDefinition *, const G4DataVector &) |
virtual G4double | ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double, G4double) |
virtual G4double | CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) |
virtual void | SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) |
Definition at line 61 of file G4PEEffectModel.hh.
G4PEEffectModel::G4PEEffectModel | ( | const G4ParticleDefinition * | p = 0 , |
|
const G4String & | nam = "PhotoElectric" | |||
) |
Definition at line 66 of file G4PEEffectModel.cc.
References G4Electron::Electron(), G4cout, G4endl, G4Gamma::Gamma(), and G4VEmModel::SetAngularDistribution().
00068 : G4VEmModel(nam) 00069 { 00070 G4cout << "### G4PEEffectModel is obsolete " 00071 << "and will be removed for the next release." << G4endl; 00072 00073 theGamma = G4Gamma::Gamma(); 00074 theElectron = G4Electron::Electron(); 00075 fminimalEnergy = 1.0*eV; 00076 fParticleChange = 0; 00077 00078 // default generator 00079 SetAngularDistribution(new G4SauterGavrilaAngularDistribution()); 00080 }
G4PEEffectModel::~G4PEEffectModel | ( | ) | [virtual] |
G4double G4PEEffectModel::ComputeCrossSectionPerAtom | ( | const G4ParticleDefinition * | , | |
G4double | kinEnergy, | |||
G4double | Z, | |||
G4double | A, | |||
G4double | , | |||
G4double | ||||
) | [virtual] |
Reimplemented from G4VEmModel.
Definition at line 99 of file G4PEEffectModel.cc.
References G4SandiaTable::GetSandiaCofPerAtom().
Referenced by G4AdjointPhotoElectricModel::AdjointCrossSectionPerAtom().
00103 { 00104 G4double* SandiaCof = G4SandiaTable::GetSandiaCofPerAtom((G4int)Z, energy); 00105 00106 G4double energy2 = energy*energy; 00107 G4double energy3 = energy*energy2; 00108 G4double energy4 = energy2*energy2; 00109 00110 return SandiaCof[0]/energy + SandiaCof[1]/energy2 + 00111 SandiaCof[2]/energy3 + SandiaCof[3]/energy4; 00112 }
G4double G4PEEffectModel::CrossSectionPerVolume | ( | const G4Material * | , | |
const G4ParticleDefinition * | , | |||
G4double | kineticEnergy, | |||
G4double | cutEnergy, | |||
G4double | maxEnergy | |||
) | [virtual] |
Reimplemented from G4VEmModel.
Definition at line 116 of file G4PEEffectModel.cc.
References G4SandiaTable::GetSandiaCofForMaterial(), and G4Material::GetSandiaTable().
00120 { 00121 G4double* SandiaCof = 00122 material->GetSandiaTable()->GetSandiaCofForMaterial(energy); 00123 00124 G4double energy2 = energy*energy; 00125 G4double energy3 = energy*energy2; 00126 G4double energy4 = energy2*energy2; 00127 00128 return SandiaCof[0]/energy + SandiaCof[1]/energy2 + 00129 SandiaCof[2]/energy3 + SandiaCof[3]/energy4; 00130 }
void G4PEEffectModel::Initialise | ( | const G4ParticleDefinition * | , | |
const G4DataVector & | ||||
) | [virtual] |
Implements G4VEmModel.
Definition at line 89 of file G4PEEffectModel.cc.
References G4VEmModel::GetParticleChangeForGamma(), and G4VEmModel::SetDeexcitationFlag().
00091 { 00092 // always false before the run 00093 SetDeexcitationFlag(false); 00094 if(!fParticleChange) { fParticleChange = GetParticleChangeForGamma(); } 00095 }
void G4PEEffectModel::SampleSecondaries | ( | std::vector< G4DynamicParticle * > * | , | |
const G4MaterialCutsCouple * | , | |||
const G4DynamicParticle * | , | |||
G4double | tmin, | |||
G4double | maxEnergy | |||
) | [virtual] |
Implements G4VEmModel.
Definition at line 134 of file G4PEEffectModel.cc.
References G4InuclSpecialFunctions::bindingEnergy(), fStopAndKill, G4VEmModel::GetAngularDistribution(), G4Element::GetAtomicShell(), G4DynamicParticle::GetKineticEnergy(), G4MaterialCutsCouple::GetMaterial(), G4DynamicParticle::GetMomentumDirection(), G4Element::GetNbOfAtomicShells(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4VParticleChange::ProposeTrackStatus(), G4VEmAngularDistribution::SampleDirection(), G4VEmModel::SelectRandomAtom(), and G4ParticleChangeForGamma::SetProposedKineticEnergy().
00139 { 00140 const G4Material* aMaterial = couple->GetMaterial(); 00141 00142 G4double energy = aDynamicPhoton->GetKineticEnergy(); 00143 G4ParticleMomentum PhotonDirection = aDynamicPhoton->GetMomentumDirection(); 00144 00145 // select randomly one element constituing the material. 00146 const G4Element* anElement = SelectRandomAtom(aMaterial,theGamma,energy); 00147 00148 // 00149 // Photo electron 00150 // 00151 00152 // Select atomic shell 00153 G4int nShells = anElement->GetNbOfAtomicShells(); 00154 G4int i = 0; 00155 while ((i<nShells) && (energy<anElement->GetAtomicShell(i))) { ++i; } 00156 00157 G4double edep = energy; 00158 00159 // no shell available 00160 if (i < nShells) { 00161 00162 G4double bindingEnergy = anElement->GetAtomicShell(i); 00163 G4double elecKineEnergy = energy - bindingEnergy; 00164 00165 if (elecKineEnergy > fminimalEnergy) { 00166 00167 edep = bindingEnergy; 00168 G4ThreeVector elecDirection = 00169 GetAngularDistribution()->SampleDirection(aDynamicPhoton, 00170 elecKineEnergy, 00171 i, 00172 couple->GetMaterial()); 00173 00174 G4DynamicParticle* aParticle = 00175 new G4DynamicParticle(theElectron, elecDirection, elecKineEnergy); 00176 fvect->push_back(aParticle); 00177 00178 } 00179 } 00180 00181 fParticleChange->SetProposedKineticEnergy(0.); 00182 fParticleChange->ProposeTrackStatus(fStopAndKill); 00183 if(edep > 0.0) { 00184 fParticleChange->ProposeLocalEnergyDeposit(edep); 00185 } 00186 }