Geant4
10.01.p02
|
#include <G4PenelopeBremsstrahlungFS.hh>
Public Member Functions | |
G4PenelopeBremsstrahlungFS (G4int verbosity=0) | |
Only master models are supposed to create instances. More... | |
~G4PenelopeBremsstrahlungFS () | |
G4double | GetEffectiveZSquared (const G4Material *mat) const |
Master and workers (do not touch tables) All of them are const. More... | |
size_t | GetNBinsX () const |
G4double | GetMomentumIntegral (G4double *y, G4double up, G4int momOrder) const |
const G4PhysicsTable * | GetScaledXSTable (const G4Material *, const G4double cut) const |
G4double | SampleGammaEnergy (G4double energy, const G4Material *, const G4double cut) const |
void | ClearTables (G4bool isMaster=true) |
Reserved for the master model: they build and handle tables. More... | |
void | BuildScaledXSTable (const G4Material *material, G4double cut, G4bool isMaster) |
void | SetVerbosity (G4int ver) |
G4int | GetVerbosity () |
Private Member Functions | |
G4PenelopeBremsstrahlungFS & | operator= (const G4PenelopeBremsstrahlungFS &right) |
G4PenelopeBremsstrahlungFS (const G4PenelopeBremsstrahlungFS &) | |
void | ReadDataFile (G4int Z) |
void | InitializeEnergySampling (const G4Material *, G4double cut) |
Private Attributes | |
std::map< std::pair< const G4Material *, G4double > , G4PhysicsTable * > * | theReducedXSTable |
std::map< const G4Material *, G4double > * | theEffectiveZSq |
G4double | theXGrid [nBinsX] |
G4double | theEGrid [nBinsE] |
std::map< G4int, G4DataVector * > * | theElementData |
std::map< std::pair< const G4Material *, G4double > , G4PhysicsTable * > * | theSamplingTable |
std::map< std::pair< const G4Material *, G4double > , G4PhysicsFreeVector * > * | thePBcut |
G4Cache< G4PhysicsFreeVector * > | fCache |
G4int | fVerbosity |
Static Private Attributes | |
static const size_t | nBinsE = 57 |
static const size_t | nBinsX = 32 |
Definition at line 61 of file G4PenelopeBremsstrahlungFS.hh.
G4PenelopeBremsstrahlungFS::G4PenelopeBremsstrahlungFS | ( | G4int | verbosity = 0 | ) |
Only master models are supposed to create instances.
Definition at line 51 of file G4PenelopeBremsstrahlungFS.cc.
References fCache, nBinsE, nBinsX, G4Cache< VALTYPE >::Put(), theEGrid, theElementData, and theXGrid.
G4PenelopeBremsstrahlungFS::~G4PenelopeBremsstrahlungFS | ( | ) |
Definition at line 73 of file G4PenelopeBremsstrahlungFS.cc.
References ClearTables(), and theElementData.
|
private |
void G4PenelopeBremsstrahlungFS::BuildScaledXSTable | ( | const G4Material * | material, |
G4double | cut, | ||
G4bool | isMaster | ||
) |
Definition at line 177 of file G4PenelopeBremsstrahlungFS.cc.
References barn, G4PhysicsVector::Energy(), eV, FatalException, fVerbosity, g, G4cout, G4endl, G4Exception(), G4Material::GetElementVector(), G4Material::GetFractionVector(), GetMomentumIntegral(), G4Material::GetName(), G4Material::GetNumberOfElements(), InitializeEnergySampling(), keV, millibarn, mole, nBinsE, nBinsX, G4PhysicsTable::push_back(), G4PhysicsFreeVector::PutValue(), ReadDataFile(), theEffectiveZSq, theEGrid, theElementData, thePBcut, theReducedXSTable, and theSamplingTable.
Referenced by G4PenelopeBremsstrahlungModel::GetCrossSectionTableForCouple(), and G4PenelopeBremsstrahlungModel::Initialise().
void G4PenelopeBremsstrahlungFS::ClearTables | ( | G4bool | isMaster = true | ) |
Reserved for the master model: they build and handle tables.
Definition at line 95 of file G4PenelopeBremsstrahlungFS.cc.
References FatalException, G4Exception(), tab(), theEffectiveZSq, thePBcut, theReducedXSTable, and theSamplingTable.
Referenced by G4PenelopeBremsstrahlungModel::ClearTables(), and ~G4PenelopeBremsstrahlungFS().
G4double G4PenelopeBremsstrahlungFS::GetEffectiveZSquared | ( | const G4Material * | mat | ) | const |
Master and workers (do not touch tables) All of them are const.
Definition at line 150 of file G4PenelopeBremsstrahlungFS.cc.
References FatalException, G4endl, G4Exception(), G4Material::GetName(), and theEffectiveZSq.
Referenced by G4PenelopeBremsstrahlungModel::BuildXSTable(), and G4PenelopeBremsstrahlungModel::GetPositronXSCorrection().
G4double G4PenelopeBremsstrahlungFS::GetMomentumIntegral | ( | G4double * | y, |
G4double | up, | ||
G4int | momOrder | ||
) | const |
Definition at line 441 of file G4PenelopeBremsstrahlungFS.cc.
References a, eps, FatalException, G4endl, G4Exception(), G4INCL::Math::max(), and G4INCL::Math::min().
Referenced by BuildScaledXSTable(), G4PenelopeBremsstrahlungModel::BuildXSTable(), and InitializeEnergySampling().
|
inline |
Definition at line 73 of file G4PenelopeBremsstrahlungFS.hh.
References nBinsX.
Referenced by G4PenelopeBremsstrahlungModel::BuildXSTable().
const G4PhysicsTable * G4PenelopeBremsstrahlungFS::GetScaledXSTable | ( | const G4Material * | mat, |
const G4double | cut | ||
) | const |
Definition at line 517 of file G4PenelopeBremsstrahlungFS.cc.
References FatalException, G4Exception(), and theReducedXSTable.
Referenced by G4PenelopeBremsstrahlungModel::BuildXSTable().
|
inline |
Definition at line 87 of file G4PenelopeBremsstrahlungFS.hh.
References fVerbosity.
|
private |
Definition at line 534 of file G4PenelopeBremsstrahlungFS.cc.
References A, FatalException, fVerbosity, G4cout, G4endl, G4Exception(), GetMomentumIntegral(), G4Material::GetName(), G4INCL::Math::max(), nBinsE, nBinsX, G4PhysicsTable::push_back(), G4PhysicsFreeVector::PutValue(), theEGrid, thePBcut, theReducedXSTable, theSamplingTable, and theXGrid.
Referenced by BuildScaledXSTable().
|
private |
|
private |
Definition at line 372 of file G4PenelopeBremsstrahlungFS.cc.
References eV, FatalException, G4endl, G4Exception(), millibarn, nBinsE, nBinsX, theEGrid, and theElementData.
Referenced by BuildScaledXSTable().
G4double G4PenelopeBremsstrahlungFS::SampleGammaEnergy | ( | G4double | energy, |
const G4Material * | mat, | ||
const G4double | cut | ||
) | const |
Definition at line 610 of file G4PenelopeBremsstrahlungFS.cc.
References A, G4INCL::KinematicsUtils::energy(), FatalException, fCache, fVerbosity, G4cout, G4endl, G4Exception(), G4UniformRand, G4Cache< VALTYPE >::Get(), G4Material::GetName(), JustWarning, keV, G4INCL::Math::max(), nBinsE, nBinsX, G4Cache< VALTYPE >::Put(), G4PhysicsFreeVector::PutValue(), theEGrid, thePBcut, theReducedXSTable, theSamplingTable, and theXGrid.
Referenced by G4PenelopeBremsstrahlungModel::SampleSecondaries().
|
inline |
Definition at line 86 of file G4PenelopeBremsstrahlungFS.hh.
References fVerbosity.
|
private |
Definition at line 133 of file G4PenelopeBremsstrahlungFS.hh.
Referenced by G4PenelopeBremsstrahlungFS(), and SampleGammaEnergy().
|
private |
Definition at line 135 of file G4PenelopeBremsstrahlungFS.hh.
Referenced by BuildScaledXSTable(), GetVerbosity(), InitializeEnergySampling(), SampleGammaEnergy(), and SetVerbosity().
|
staticprivate |
Definition at line 105 of file G4PenelopeBremsstrahlungFS.hh.
Referenced by BuildScaledXSTable(), G4PenelopeBremsstrahlungFS(), InitializeEnergySampling(), ReadDataFile(), and SampleGammaEnergy().
|
staticprivate |
Definition at line 106 of file G4PenelopeBremsstrahlungFS.hh.
Referenced by BuildScaledXSTable(), G4PenelopeBremsstrahlungFS(), GetNBinsX(), InitializeEnergySampling(), ReadDataFile(), and SampleGammaEnergy().
|
private |
Definition at line 101 of file G4PenelopeBremsstrahlungFS.hh.
Referenced by BuildScaledXSTable(), ClearTables(), and GetEffectiveZSquared().
Definition at line 109 of file G4PenelopeBremsstrahlungFS.hh.
Referenced by BuildScaledXSTable(), G4PenelopeBremsstrahlungFS(), InitializeEnergySampling(), ReadDataFile(), and SampleGammaEnergy().
|
private |
Definition at line 117 of file G4PenelopeBremsstrahlungFS.hh.
Referenced by BuildScaledXSTable(), G4PenelopeBremsstrahlungFS(), ReadDataFile(), and ~G4PenelopeBremsstrahlungFS().
|
private |
Definition at line 127 of file G4PenelopeBremsstrahlungFS.hh.
Referenced by BuildScaledXSTable(), ClearTables(), InitializeEnergySampling(), and SampleGammaEnergy().
|
private |
Definition at line 99 of file G4PenelopeBremsstrahlungFS.hh.
Referenced by BuildScaledXSTable(), ClearTables(), GetScaledXSTable(), InitializeEnergySampling(), and SampleGammaEnergy().
|
private |
Definition at line 125 of file G4PenelopeBremsstrahlungFS.hh.
Referenced by BuildScaledXSTable(), ClearTables(), InitializeEnergySampling(), and SampleGammaEnergy().
Definition at line 108 of file G4PenelopeBremsstrahlungFS.hh.
Referenced by G4PenelopeBremsstrahlungFS(), InitializeEnergySampling(), and SampleGammaEnergy().