Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
|
Functions | |
G4double | toDegrees (G4double radians) |
G4int | heaviside (G4int n) |
G4double | pow13 (G4double x) |
G4double | powMinus13 (G4double x) |
G4double | pow23 (G4double x) |
G4double | aSinH (G4double x) |
template<typename T > | |
G4int | sign (const T t) |
template<typename T > | |
T | max (const T t1, const T t2) |
brief Return the largest of the two arguments | |
template<typename T > | |
T | min (const T t1, const T t2) |
brief Return the smallest of the two arguments | |
G4double | gaussianCDF (const G4double x) |
Cumulative distribution function for Gaussian. | |
G4double | gaussianCDF (const G4double x, const G4double x0, const G4double sigma) |
Generic cumulative distribution function for Gaussian. | |
G4double | inverseGaussianCDF (const G4double x) |
Inverse cumulative distribution function for Gaussian. | |
G4double | arcSin (const G4double x) |
Calculates arcsin with some tolerance on illegal arguments. | |
G4double | arcCos (const G4double x) |
Calculates arccos with some tolerance on illegal arguments. | |
Variables | |
const G4double | pi = 3.14159265358979323846264338328 |
const G4double | twoPi = 2.0 * pi |
const G4double | tenPi = 10.0 * pi |
const G4double | piOverTwo = 0.5 * pi |
const G4double | oneOverSqrtTwo = 1./std::sqrt((G4double)2.) |
const G4double | oneOverSqrtThree = 1./std::sqrt((G4double)3.) |
const G4double | oneThird = 1./3. |
const G4double | twoThirds = 2./3. |
const G4double | sqrtFiveThirds = std::sqrt(5./3.) |
const G4double | sqrtThreeFifths = std::sqrt(3./5.) |
Calculates arccos with some tolerance on illegal arguments.
Definition at line 103 of file G4INCLGlobals.cc.
Referenced by G4INCL::CoulombNonRelativistic::distortOut(), G4INCL::EventInfo::fillInverseKinematics(), G4INCL::SurfaceAvatar::getTransmissionProbability(), and G4INCL::EventInfo::remnantToParticle().
Calculates arcsin with some tolerance on illegal arguments.
Definition at line 98 of file G4INCLGlobals.cc.
Definition at line 100 of file G4INCLGlobals.hh.
Cumulative distribution function for Gaussian.
A public-domain approximation taken from Abramowitz and Stegun. Applies to a Gaussian with mean=0 and sigma=1.
x | a Gaussian variable |
Definition at line 74 of file G4INCLGlobals.cc.
Referenced by G4INCL::Random::correlatedUniform(), and gaussianCDF().
Generic cumulative distribution function for Gaussian.
A public-domain approximation taken from Abramowitz and Stegun. Applies to a generic Gaussian.
x | a Gaussian variable |
x0 | mean of the Gaussian |
sigma | standard deviation of the Gaussian |
Definition at line 87 of file G4INCLGlobals.cc.
Definition at line 83 of file G4INCLGlobals.hh.
Referenced by G4INCL::Nucleus::insertParticle().
Inverse cumulative distribution function for Gaussian.
A public-domain approximation taken from Abramowitz and Stegun. Applies to a Gaussian with mean=0 and sigma=1.
x | a uniform variate |
Definition at line 91 of file G4INCLGlobals.cc.
|
inline |
brief Return the largest of the two arguments
Definition at line 112 of file G4INCLGlobals.hh.
Referenced by G4INCL::ParticleTable::getLargestNuclearRadius().
|
inline |
brief Return the smallest of the two arguments
Definition at line 117 of file G4INCLGlobals.hh.
Definition at line 88 of file G4INCLGlobals.hh.
Referenced by G4INCL::NuclearDensityFactory::createRPCorrelationTable(), G4INCL::NuclearPotential::INuclearPotential::INuclearPotential(), and G4INCL::Random::sphereVector().
Definition at line 96 of file G4INCLGlobals.hh.
Definition at line 92 of file G4INCLGlobals.hh.
|
inline |
A simple sign function that allows us to port fortran code to c++ more easily.
Definition at line 107 of file G4INCLGlobals.hh.
Referenced by G4INCL::DeltaProductionChannel::fillFinalState(), G4INCL::ElasticChannel::fillFinalState(), G4INCL::NNToMissingStrangenessChannel::fillFinalState(), G4INCL::NpiToMissingStrangenessChannel::fillFinalState(), gaussianCDF(), and G4INCL::RootFinder::solve().
Definition at line 79 of file G4INCLGlobals.hh.
Referenced by G4INCL::Nucleus::fillEventInfo(), G4INCL::EventInfo::fillInverseKinematics(), and G4INCL::EventInfo::remnantToParticle().
Definition at line 73 of file G4INCLGlobals.hh.
Referenced by G4INCL::NuclearDensityFactory::createPCDFTable(), G4INCL::NuclearDensityFactory::createRCDFTable(), G4INCL::NuclearDensityFactory::createRPCorrelationTable(), and G4INCL::Random::gaussVector().
Definition at line 72 of file G4INCLGlobals.hh.
Referenced by gaussianCDF().
const G4double G4INCL::Math::oneThird = 1./3. |
Definition at line 74 of file G4INCLGlobals.hh.
Referenced by pow13(), and powMinus13().
const G4double G4INCL::Math::pi = 3.14159265358979323846264338328 |
Definition at line 68 of file G4INCLGlobals.hh.
Referenced by arcCos(), arcSin(), G4INCL::NNToNSKpiChannel::fillFinalState(), G4INCL::NpiToLK2piChannel::fillFinalState(), G4INCL::NpiToLKpiChannel::fillFinalState(), G4INCL::NpiToNKKbChannel::fillFinalState(), G4INCL::NpiToSK2piChannel::fillFinalState(), G4INCL::NpiToSKpiChannel::fillFinalState(), G4INCL::PauliStandard::getBlockingProbability(), and toDegrees().
Definition at line 71 of file G4INCLGlobals.hh.
Referenced by G4INCL::CoulombNonRelativistic::distortOut().
const G4double G4INCL::Math::sqrtFiveThirds = std::sqrt(5./3.) |
Definition at line 76 of file G4INCLGlobals.hh.
Referenced by G4INCL::ParticleTable::getFermiMomentumConstantLight().
const G4double G4INCL::Math::sqrtThreeFifths = std::sqrt(3./5.) |
Definition at line 77 of file G4INCLGlobals.hh.
Referenced by G4INCL::ParticleTable::getMomentumRMS().
Definition at line 70 of file G4INCLGlobals.hh.
Referenced by G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar(), G4INCL::BinaryCollisionAvatar::getChannel(), G4INCL::CrossSections::interactionDistanceKbarN(), G4INCL::CrossSections::interactionDistanceKN(), G4INCL::CrossSections::interactionDistanceNN(), G4INCL::CrossSections::interactionDistancePiN(), G4INCL::CrossSections::interactionDistanceYN(), and G4INCL::INCL::prepareReaction().
Definition at line 69 of file G4INCLGlobals.hh.
Referenced by G4INCL::DeuteronDensity::derivDensityR(), G4INCL::DeltaProductionChannel::fillFinalState(), G4INCL::ElasticChannel::fillFinalState(), G4INCL::Random::gaussWithMemory(), G4INCL::PauliStandard::getBlockingProbability(), G4INCL::NKbElasticChannel::KaonMomentum(), G4INCL::NKbToLpiChannel::KaonMomentum(), G4INCL::NKbToNKbChannel::KaonMomentum(), G4INCL::NKbToSpiChannel::KaonMomentum(), G4INCL::NpiToLKChannel::KaonMomentum(), G4INCL::NpiToSKChannel::KaonMomentum(), and G4INCL::Random::normVector().
const G4double G4INCL::Math::twoThirds = 2./3. |
Definition at line 75 of file G4INCLGlobals.hh.
Referenced by G4INCL::InteractionAvatar::postInteraction(), and pow23().