Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4INCL::Math Namespace Reference

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 >
max (const T t1, const T t2)
 brief Return the largest of the two arguments
 
template<typename 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.)
 

Function Documentation

◆ arcCos()

G4double G4INCL::Math::arcCos ( const G4double x)

Calculates arccos with some tolerance on illegal arguments.

Definition at line 103 of file G4INCLGlobals.cc.

103 {
104// assert(x>-1.000001 && x<1.000001);
105 return ((x > 1.) ? 0. : ((x<-1.) ? pi : std::acos(x)));
106 }

Referenced by G4INCL::CoulombNonRelativistic::distortOut(), G4INCL::EventInfo::fillInverseKinematics(), G4INCL::SurfaceAvatar::getTransmissionProbability(), and G4INCL::EventInfo::remnantToParticle().

◆ arcSin()

G4double G4INCL::Math::arcSin ( const G4double x)

Calculates arcsin with some tolerance on illegal arguments.

Definition at line 98 of file G4INCLGlobals.cc.

98 {
99// assert(x>-1.000001 && x<1.000001);
100 return ((x > 1.) ? 0. : ((x<-1.) ? pi : std::asin(x)));
101 }

◆ aSinH()

G4double G4INCL::Math::aSinH ( G4double x)
inline

Definition at line 100 of file G4INCLGlobals.hh.

100 {
101 return std::log(x + std::sqrt(x*x+1.));
102 }

◆ gaussianCDF() [1/2]

G4double G4INCL::Math::gaussianCDF ( const G4double x)

Cumulative distribution function for Gaussian.

A public-domain approximation taken from Abramowitz and Stegun. Applies to a Gaussian with mean=0 and sigma=1.

Parameters
xa Gaussian variable

Definition at line 74 of file G4INCLGlobals.cc.

75 {
76 // Save the sign of x
77 const G4double sgn = sign(x);
78 const G4double z = std::fabs(x) * oneOverSqrtTwo;
79
80 // A&S formula 7.1.26
81 G4double t = 1.0/(1.0 + gcdfp*z);
82 G4double y = 1.0 - (((((gcdfa5*t + gcdfa4)*t) + gcdfa3)*t + gcdfa2)*t + gcdfa1)*t*std::exp(-z*z);
83
84 return 0.5*(1.0 + sgn*y);
85 }
double G4double
Definition G4Types.hh:83
const G4double oneOverSqrtTwo
G4int sign(const T t)

Referenced by G4INCL::Random::correlatedUniform(), and gaussianCDF().

◆ gaussianCDF() [2/2]

G4double G4INCL::Math::gaussianCDF ( const G4double x,
const G4double x0,
const G4double sigma )

Generic cumulative distribution function for Gaussian.

A public-domain approximation taken from Abramowitz and Stegun. Applies to a generic Gaussian.

Parameters
xa Gaussian variable
x0mean of the Gaussian
sigmastandard deviation of the Gaussian

Definition at line 87 of file G4INCLGlobals.cc.

87 {
88 return gaussianCDF((x-x0)/sigma);
89 }
G4double gaussianCDF(const G4double x)
Cumulative distribution function for Gaussian.

◆ heaviside()

G4int G4INCL::Math::heaviside ( G4int n)
inline

Definition at line 83 of file G4INCLGlobals.hh.

83 {
84 if(n < 0) return 0;
85 else return 1;
86 }

Referenced by G4INCL::Nucleus::insertParticle().

◆ inverseGaussianCDF()

G4double G4INCL::Math::inverseGaussianCDF ( const G4double x)

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.

Parameters
xa uniform variate
Returns
a Gaussian variate

Definition at line 91 of file G4INCLGlobals.cc.

91 {
92 if (x < 0.5)
93 return -inverseGaussianCDFRational( std::sqrt(-2.0*std::log(x)) );
94 else
95 return inverseGaussianCDFRational( std::sqrt(-2.0*std::log(1.-x)) );
96 }

◆ max()

template<typename T >
T G4INCL::Math::max ( const T t1,
const T t2 )
inline

brief Return the largest of the two arguments

Definition at line 112 of file G4INCLGlobals.hh.

112 {
113 return t1 > t2 ? t1 : t2;
114 }

Referenced by G4INCL::ParticleTable::getLargestNuclearRadius().

◆ min()

template<typename T >
T G4INCL::Math::min ( const T t1,
const T t2 )
inline

brief Return the smallest of the two arguments

Definition at line 117 of file G4INCLGlobals.hh.

117 {
118 return t1 < t2 ? t1 : t2;
119 }

◆ pow13()

G4double G4INCL::Math::pow13 ( G4double x)
inline

◆ pow23()

G4double G4INCL::Math::pow23 ( G4double x)
inline

Definition at line 96 of file G4INCLGlobals.hh.

96 {
97 return std::pow(x, twoThirds);
98 }

◆ powMinus13()

G4double G4INCL::Math::powMinus13 ( G4double x)
inline

Definition at line 92 of file G4INCLGlobals.hh.

92 {
93 return std::pow(x, -oneThird);
94 }

◆ sign()

template<typename T >
G4int G4INCL::Math::sign ( const T t)
inline

A simple sign function that allows us to port fortran code to c++ more easily.

Definition at line 107 of file G4INCLGlobals.hh.

107 {
108 return t > 0 ? 1: t < 0 ? -1 : 0;
109 }

Referenced by G4INCL::DeltaProductionChannel::fillFinalState(), G4INCL::ElasticChannel::fillFinalState(), G4INCL::NNToMissingStrangenessChannel::fillFinalState(), G4INCL::NpiToMissingStrangenessChannel::fillFinalState(), gaussianCDF(), and G4INCL::RootFinder::solve().

◆ toDegrees()

G4double G4INCL::Math::toDegrees ( G4double radians)
inline

Definition at line 79 of file G4INCLGlobals.hh.

79 {
80 return radians * (180.0 / pi);
81 }

Referenced by G4INCL::Nucleus::fillEventInfo(), G4INCL::EventInfo::fillInverseKinematics(), and G4INCL::EventInfo::remnantToParticle().

Variable Documentation

◆ oneOverSqrtThree

◆ oneOverSqrtTwo

const G4double G4INCL::Math::oneOverSqrtTwo = 1./std::sqrt((G4double)2.)

Definition at line 72 of file G4INCLGlobals.hh.

Referenced by gaussianCDF().

◆ oneThird

const G4double G4INCL::Math::oneThird = 1./3.

Definition at line 74 of file G4INCLGlobals.hh.

Referenced by pow13(), and powMinus13().

◆ pi

◆ piOverTwo

const G4double G4INCL::Math::piOverTwo = 0.5 * pi

Definition at line 71 of file G4INCLGlobals.hh.

Referenced by G4INCL::CoulombNonRelativistic::distortOut().

◆ sqrtFiveThirds

const G4double G4INCL::Math::sqrtFiveThirds = std::sqrt(5./3.)

◆ sqrtThreeFifths

const G4double G4INCL::Math::sqrtThreeFifths = std::sqrt(3./5.)

Definition at line 77 of file G4INCLGlobals.hh.

Referenced by G4INCL::ParticleTable::getMomentumRMS().

◆ tenPi

◆ twoPi

◆ twoThirds

const G4double G4INCL::Math::twoThirds = 2./3.

Definition at line 75 of file G4INCLGlobals.hh.

Referenced by G4INCL::InteractionAvatar::postInteraction(), and pow23().