Geant4-11
|
#include <G4PenelopeCrossSection.hh>
Public Member Functions | |
void | AddCrossSectionPoint (size_t binNumber, G4double energy, G4double XH0, G4double XH1, G4double XH2, G4double XS0, G4double XS1, G4double XS2) |
void | AddShellCrossSectionPoint (size_t binNumber, size_t shellID, G4double energy, G4double xs) |
G4PenelopeCrossSection (const G4PenelopeCrossSection &)=delete | |
G4PenelopeCrossSection (size_t nOfEnergyPoints, size_t nOfShells=0) | |
G4double | GetHardCrossSection (G4double energy) const |
Returns hard cross section at the given energy. More... | |
G4double | GetNormalizedShellCrossSection (size_t shellID, G4double energy) const |
Returns the hard cross section for the given shell (normalized to 1) More... | |
size_t | GetNumberOfShells () const |
G4double | GetShellCrossSection (size_t shellID, G4double energy) const |
Returns the hard cross section for the given shell (per molecule) More... | |
G4double | GetSoftStoppingPower (G4double energy) const |
Returns the total stopping power due to soft collisions. More... | |
G4double | GetTotalCrossSection (G4double energy) const |
Returns total cross section at the given energy. More... | |
void | NormalizeShellCrossSections () |
G4PenelopeCrossSection & | operator= (const G4PenelopeCrossSection &right)=delete |
~G4PenelopeCrossSection () | |
Private Attributes | |
G4PhysicsTable * | fHardCrossSections |
G4bool | fIsNormalized |
size_t | fNumberOfEnergyPoints |
size_t | fNumberOfShells |
G4PhysicsTable * | fShellCrossSections |
G4PhysicsTable * | fShellNormalizedCrossSections |
G4PhysicsTable * | fSoftCrossSections |
Definition at line 71 of file G4PenelopeCrossSection.hh.
|
explicit |
Definition at line 45 of file G4PenelopeCrossSection.cc.
References FatalException, fHardCrossSections, fIsNormalized, fNumberOfEnergyPoints, fNumberOfShells, fShellCrossSections, fShellNormalizedCrossSections, fSoftCrossSections, G4endl, G4Exception(), and G4PhysicsTable::push_back().
G4PenelopeCrossSection::~G4PenelopeCrossSection | ( | ) |
Definition at line 100 of file G4PenelopeCrossSection.cc.
References G4PhysicsTable::clearAndDestroy(), fHardCrossSections, fShellCrossSections, fShellNormalizedCrossSections, and fSoftCrossSections.
|
delete |
void G4PenelopeCrossSection::AddCrossSectionPoint | ( | size_t | binNumber, |
G4double | energy, | ||
G4double | XH0, | ||
G4double | XH1, | ||
G4double | XH2, | ||
G4double | XS0, | ||
G4double | XS1, | ||
G4double | XS2 | ||
) |
Public interface for the master thread
Definition at line 126 of file G4PenelopeCrossSection.cc.
References cm2, G4INCL::KinematicsUtils::energy(), eV, fHardCrossSections, fNumberOfEnergyPoints, fSoftCrossSections, G4cout, G4endl, G4Log(), G4INCL::Math::max(), and G4PhysicsFreeVector::PutValues().
Referenced by G4PenelopeIonisationXSHandler::BuildXSTable(), and G4PenelopeBremsstrahlungModel::BuildXSTable().
void G4PenelopeCrossSection::AddShellCrossSectionPoint | ( | size_t | binNumber, |
size_t | shellID, | ||
G4double | energy, | ||
G4double | xs | ||
) |
Definition at line 186 of file G4PenelopeCrossSection.cc.
References cm2, G4INCL::KinematicsUtils::energy(), fNumberOfEnergyPoints, fNumberOfShells, fShellCrossSections, G4cout, G4endl, G4Log(), G4INCL::Math::max(), and G4PhysicsFreeVector::PutValues().
Referenced by G4PenelopeIonisationXSHandler::BuildXSTable().
Returns hard cross section at the given energy.
Definition at line 270 of file G4PenelopeCrossSection.cc.
References G4INCL::KinematicsUtils::energy(), fHardCrossSections, fNumberOfEnergyPoints, G4cout, G4endl, G4Exp(), G4Log(), G4PhysicsVector::GetVectorLength(), and G4PhysicsVector::Value().
Referenced by G4PenelopeBremsstrahlungModel::CrossSectionPerVolume(), and G4PenelopeIonisationModel::CrossSectionPerVolume().
G4double G4PenelopeCrossSection::GetNormalizedShellCrossSection | ( | size_t | shellID, |
G4double | energy | ||
) | const |
Returns the hard cross section for the given shell (normalized to 1)
Definition at line 364 of file G4PenelopeCrossSection.cc.
References G4INCL::KinematicsUtils::energy(), fIsNormalized, fNumberOfEnergyPoints, fNumberOfShells, fShellNormalizedCrossSections, G4cout, G4endl, G4Exp(), G4Log(), G4PhysicsVector::GetVectorLength(), and G4PhysicsVector::Value().
Referenced by G4PenelopeIonisationModel::SampleFinalStateElectron(), and G4PenelopeIonisationModel::SampleFinalStatePositron().
|
inline |
Returns the hard cross section for the given shell (per molecule)
Definition at line 328 of file G4PenelopeCrossSection.cc.
References G4INCL::KinematicsUtils::energy(), fNumberOfEnergyPoints, fNumberOfShells, fShellCrossSections, G4cout, G4endl, G4Exp(), G4Log(), G4PhysicsVector::GetVectorLength(), and G4PhysicsVector::Value().
Referenced by G4PenelopeIonisationCrossSection::CrossSection().
Returns the total stopping power due to soft collisions.
Definition at line 299 of file G4PenelopeCrossSection.cc.
References G4INCL::KinematicsUtils::energy(), fNumberOfEnergyPoints, fSoftCrossSections, G4cout, G4endl, G4Exp(), G4Log(), G4PhysicsVector::GetVectorLength(), and G4PhysicsVector::Value().
Referenced by G4PenelopeBremsstrahlungModel::ComputeDEDXPerVolume(), and G4PenelopeIonisationModel::ComputeDEDXPerVolume().
Returns total cross section at the given energy.
Definition at line 227 of file G4PenelopeCrossSection.cc.
References G4INCL::KinematicsUtils::energy(), fHardCrossSections, fNumberOfEnergyPoints, fSoftCrossSections, G4cout, G4endl, G4Exp(), G4Log(), G4PhysicsVector::GetVectorLength(), and G4PhysicsVector::Value().
Referenced by G4PenelopeIonisationModel::CrossSectionPerVolume().
void G4PenelopeCrossSection::NormalizeShellCrossSections | ( | ) |
Definition at line 409 of file G4PenelopeCrossSection.cc.
References fIsNormalized, fNumberOfEnergyPoints, fNumberOfShells, fShellCrossSections, fShellNormalizedCrossSections, G4cout, G4endl, G4Exp(), G4Log(), G4PhysicsVector::GetLowEdgeEnergy(), and G4PhysicsFreeVector::PutValues().
Referenced by G4PenelopeIonisationXSHandler::BuildXSTable().
|
delete |
|
private |
Definition at line 116 of file G4PenelopeCrossSection.hh.
Referenced by AddCrossSectionPoint(), G4PenelopeCrossSection(), GetHardCrossSection(), GetTotalCrossSection(), and ~G4PenelopeCrossSection().
|
private |
Definition at line 124 of file G4PenelopeCrossSection.hh.
Referenced by G4PenelopeCrossSection(), GetNormalizedShellCrossSection(), and NormalizeShellCrossSections().
|
private |
Definition at line 122 of file G4PenelopeCrossSection.hh.
Referenced by AddCrossSectionPoint(), AddShellCrossSectionPoint(), G4PenelopeCrossSection(), GetHardCrossSection(), GetNormalizedShellCrossSection(), GetShellCrossSection(), GetSoftStoppingPower(), GetTotalCrossSection(), and NormalizeShellCrossSections().
|
private |
Definition at line 123 of file G4PenelopeCrossSection.hh.
Referenced by AddShellCrossSectionPoint(), G4PenelopeCrossSection(), GetNormalizedShellCrossSection(), GetNumberOfShells(), GetShellCrossSection(), and NormalizeShellCrossSections().
|
private |
Definition at line 119 of file G4PenelopeCrossSection.hh.
Referenced by AddShellCrossSectionPoint(), G4PenelopeCrossSection(), GetShellCrossSection(), NormalizeShellCrossSections(), and ~G4PenelopeCrossSection().
|
private |
Definition at line 120 of file G4PenelopeCrossSection.hh.
Referenced by G4PenelopeCrossSection(), GetNormalizedShellCrossSection(), NormalizeShellCrossSections(), and ~G4PenelopeCrossSection().
|
private |
Definition at line 113 of file G4PenelopeCrossSection.hh.
Referenced by AddCrossSectionPoint(), G4PenelopeCrossSection(), GetSoftStoppingPower(), GetTotalCrossSection(), and ~G4PenelopeCrossSection().