#include <G4UAtomicDeexcitation.hh>
|
| G4UAtomicDeexcitation () |
|
virtual | ~G4UAtomicDeexcitation () |
|
virtual void | InitialiseForNewRun () |
|
virtual void | InitialiseForExtraAtom (G4int Z) |
|
void | SetCutForSecondaryPhotons (G4double cut) |
|
void | SetCutForAugerElectrons (G4double cut) |
|
virtual const G4AtomicShell * | GetAtomicShell (G4int Z, G4AtomicShellEnumerator shell) |
|
virtual void | GenerateParticles (std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4double gammaCut, G4double eCut) |
|
virtual G4double | GetShellIonisationCrossSectionPerAtom (const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0) |
|
virtual G4double | ComputeShellIonisationCrossSectionPerAtom (const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0) |
|
| G4VAtomDeexcitation (const G4String &modname="Deexcitation", const G4String &pixename="") |
|
virtual | ~G4VAtomDeexcitation () |
|
void | InitialiseAtomicDeexcitation () |
|
void | SetDeexcitationActiveRegion (const G4String &rname, G4bool valDeexcitation, G4bool valAuger, G4bool valPIXE) |
|
void | SetFluo (G4bool) |
|
G4bool | IsFluoActive () const |
|
void | SetAuger (G4bool) |
|
G4bool | IsAugerActive () const |
|
void | SetPIXE (G4bool) |
|
G4bool | IsPIXEActive () const |
|
const G4String & | GetName () const |
|
void | SetPIXECrossSectionModel (const G4String &) |
|
const G4String & | PIXECrossSectionModel () const |
|
void | SetPIXEElectronCrossSectionModel (const G4String &) |
|
const G4String & | PIXEElectronCrossSectionModel () const |
|
const std::vector< G4bool > & | GetListOfActiveAtoms () const |
|
void | SetVerboseLevel (G4int) |
|
G4int | GetVerboseLevel () const |
|
G4bool | CheckDeexcitationActiveRegion (G4int coupleIndex) |
|
G4bool | CheckAugerActiveRegion (G4int coupleIndex) |
|
void | GenerateParticles (std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex) |
|
void | AlongStepDeexcitation (std::vector< G4Track * > &tracks, const G4Step &step, G4double &eLoss, G4int coupleIndex) |
|
Definition at line 61 of file G4UAtomicDeexcitation.hh.
G4UAtomicDeexcitation::G4UAtomicDeexcitation |
( |
| ) |
|
G4UAtomicDeexcitation::~G4UAtomicDeexcitation |
( |
| ) |
|
|
virtual |
Implements G4VAtomDeexcitation.
Definition at line 233 of file G4UAtomicDeexcitation.cc.
References G4Exception(), JustWarning, and G4AtomicShell::ShellId().
244 minGammaEnergy = gammaCut;
245 minElectronEnergy = eCut;
252 G4int provShellId = 0;
267 provShellId = SelectTypeOfTransition(Z, givenShellId);
271 aParticle = GenerateFluorescence(Z,givenShellId,provShellId);
274 else if ( provShellId == -1)
277 aParticle = GenerateAuger(Z, givenShellId);
282 G4Exception(
"G4UAtomicDeexcitation::GenerateParticles()",
"de0002",
JustWarning,
"Energy deposited locally");
289 provShellId = SelectTypeOfTransition(Z,newShellId);
292 aParticle = GenerateFluorescence(Z,newShellId,provShellId);
295 else if ( provShellId == -1)
298 aParticle = GenerateAuger(Z, newShellId);
303 G4Exception(
"G4UAtomicDeexcitation::GenerateParticles()",
"de0002",
JustWarning,
"Energy deposited locally");
309 vectorOfParticles->push_back(aParticle);
312 else {provShellId = -2;}
314 while (provShellId > -2);
318 G4Exception(
"G4UAtomicDeexcitation::GenerateParticles()",
"de0001",
JustWarning,
"Energy deposited locally");
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Implements G4VAtomDeexcitation.
Definition at line 326 of file G4UAtomicDeexcitation.cc.
References G4VhShellCrossSection::CrossSection(), G4EmCorrections::EffectiveChargeSquareRatio(), python.hepunit::eplus, G4AtomicShells::GetNumberOfShells(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), and python.hepunit::proton_mass_c2.
Referenced by ComputeShellIonisationCrossSectionPerAtom().
339 if(Z > 93 || Z < 6 ) {
return xsec; }
344 if(pdef == theElectron || pdef == thePositron) {
345 xsec = ePIXEshellCS->
CrossSection(Z,shellEnum,kineticEnergy,0.0,mat);
357 escaled = kineticEnergy*mass/(pdef->
GetPDGMass());
367 if(PIXEshellCS) { xsec = PIXEshellCS->
CrossSection(Z,shellEnum,escaled,mass,mat); }
370 xsec = anaPIXEshellCS->
CrossSection(Z,shellEnum,escaled,mass,mat);
374 if (q2) {xsec *= q2;}
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
const G4String & GetParticleName() const
G4double GetPDGMass() const
G4double GetPDGCharge() const
virtual G4double CrossSection(G4int Z, G4AtomicShellEnumerator shell, G4double incidentEnergy, G4double mass, const G4Material *mat)=0
static G4int GetNumberOfShells(G4int Z)
void G4UAtomicDeexcitation::InitialiseForExtraAtom |
( |
G4int |
Z | ) |
|
|
virtual |
void G4UAtomicDeexcitation::InitialiseForNewRun |
( |
| ) |
|
|
virtual |
Implements G4VAtomDeexcitation.
Definition at line 97 of file G4UAtomicDeexcitation.cc.
References G4cout, G4endl, G4VhShellCrossSection::GetName(), G4AtomicTransitionManager::Instance(), G4VAtomDeexcitation::IsFluoActive(), G4VAtomDeexcitation::IsPIXEActive(), G4VAtomDeexcitation::PIXECrossSectionModel(), G4VAtomDeexcitation::PIXEElectronCrossSectionModel(), G4VAtomDeexcitation::SetPIXECrossSectionModel(), and G4VAtomDeexcitation::SetPIXEElectronCrossSectionModel().
103 G4cout <<
"### === G4UAtomicDeexcitation::InitialiseForNewRun()" <<
G4endl;
130 G4cout <<
"### G4UAtomicDeexcitation::InitialiseForNewRun WARNING "
133 <<
" is unknown, Analytical cross section will be used" <<
G4endl;
184 G4cout <<
"### G4UAtomicDeexcitation::InitialiseForNewRun WARNING "
187 <<
" is unknown, PIXE is disabled" <<
G4endl;
void SetPIXEElectronCrossSectionModel(const G4String &)
G4bool IsFluoActive() const
G4bool IsPIXEActive() const
const G4String & GetName() const
void SetPIXECrossSectionModel(const G4String &)
const G4String & PIXEElectronCrossSectionModel() const
G4GLOB_DLL std::ostream G4cout
const G4String & PIXECrossSectionModel() const
static G4AtomicTransitionManager * Instance()
void G4UAtomicDeexcitation::SetCutForAugerElectrons |
( |
G4double |
cut | ) |
|
void G4UAtomicDeexcitation::SetCutForSecondaryPhotons |
( |
G4double |
cut | ) |
|
The documentation for this class was generated from the following files: