Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization measurement in gases and semiconductors
Loading...
Searching...
No Matches
Garfield::Numerics Namespace Reference

Collection of numerical routines. More...

Namespaces

namespace  CERNLIB
 Linear algebra routines from CERNLIB.
namespace  QUADPACK
 Functions for performing numerical integration (quadrature).

Functions

constexpr std::array< double, 2 > GaussLegendreNodes2 ()
constexpr std::array< double, 2 > GaussLegendreWeights2 ()
constexpr std::array< double, 3 > GaussLegendreNodes3 ()
constexpr std::array< double, 3 > GaussLegendreWeights3 ()
constexpr std::array< double, 4 > GaussLegendreNodes4 ()
constexpr std::array< double, 4 > GaussLegendreWeights4 ()
constexpr std::array< double, 5 > GaussLegendreNodes5 ()
constexpr std::array< double, 5 > GaussLegendreWeights5 ()
constexpr std::array< double, 6 > GaussLegendreNodes6 ()
constexpr std::array< double, 6 > GaussLegendreWeights6 ()
double Legendre (const unsigned int n, const double x)
 Legendre polynomials.
double BesselI0S (const double xx)
 Modified Bessel functions.
double BesselI1S (const double xx)
double BesselK0S (const double xx)
double BesselK0L (const double xx)
double BesselK1S (const double xx)
double BesselK1L (const double xx)
double Divdif (const std::vector< double > &f, const std::vector< double > &a, int nn, double x, int mm)
 C++ version of DIVDIF (CERN program library E105) which performs tabular interpolation using symmetrically placed argument points.
double LinearInterpolation (const std::vector< double > &ytab, const std::vector< double > &xtab, const double x)
bool Boxin2 (const std::vector< std::vector< double > > &value, const std::vector< double > &xAxis, const std::vector< double > &yAxis, const int nx, const int ny, const double xx, const double yy, double &f, const int iOrder)
 Interpolation of order 1 and 2 in an irregular rectangular two-dimensional grid.
bool Boxin3 (const std::vector< std::vector< std::vector< double > > > &value, const std::vector< double > &xAxis, const std::vector< double > &yAxis, const std::vector< double > &zAxis, const int nx, const int ny, const int nz, const double xx, const double yy, const double zz, double &f, const int iOrder)
 Interpolation of order 1 and 2 in an irregular rectangular three-dimensional grid.
bool LeastSquaresFit (std::function< double(double, const std::vector< double > &)> f, std::vector< double > &par, std::vector< double > &epar, const std::vector< double > &x, const std::vector< double > &y, const std::vector< double > &ey, const unsigned int nMaxIter, const double diff, double &chi2, const double eps, const bool debug, const bool verbose)
 Least-squares minimisation.

Detailed Description

Collection of numerical routines.

Function Documentation

◆ BesselI0S()

double Garfield::Numerics::BesselI0S ( const double xx)
inline

Modified Bessel functions.

Series expansions from Abramowitz and Stegun.

Definition at line 174 of file Numerics.hh.

174 {
175 const double y = xx / 3.75;
176 const double y2 = y * y;
177 return 1. + 3.5156229 * y2 + 3.0899424 * y2 * y2 + 1.2067492 * pow(y2, 3) +
178 0.2659732 * pow(y2, 4) + 0.0360768 * pow(y2, 5) +
179 0.0045813 * pow(y2, 6);
180}

◆ BesselI1S()

double Garfield::Numerics::BesselI1S ( const double xx)
inline

Definition at line 182 of file Numerics.hh.

182 {
183 const double y = xx / 3.75;
184 const double y2 = y * y;
185 return xx * (0.5 + 0.87890594 * y2 + 0.51498869 * y2 * y2 +
186 0.15084934 * pow(y2, 3) + 0.02658733 * pow(y2, 4) +
187 0.00301532 * pow(y2, 5) + 0.00032411 * pow(y2, 6));
188}

◆ BesselK0L()

double Garfield::Numerics::BesselK0L ( const double xx)
inline

Definition at line 199 of file Numerics.hh.

199 {
200 const double y = 2. / xx;
201 return (exp(-xx) / sqrt(xx)) *
202 (1.25331414 - 0.07832358 * y + 0.02189568 * y * y -
203 0.01062446 * pow(y, 3) + 0.00587872 * pow(y, 4) -
204 0.00251540 * pow(y, 5) + 0.00053208 * pow(y, 6));
205}

◆ BesselK0S()

double Garfield::Numerics::BesselK0S ( const double xx)
inline

Definition at line 190 of file Numerics.hh.

190 {
191 const double y = xx / 2.;
192 const double y2 = y * y;
193 return -log(y) * BesselI0S(xx) - Gamma + 0.42278420 * y2 +
194 0.23069756 * y2 * y2 + 0.03488590 * pow(y2, 3) +
195 0.00262698 * pow(y2, 4) + 0.00010750 * pow(y2, 5) +
196 0.00000740 * pow(y2, 6);
197}
double BesselI0S(const double xx)
Modified Bessel functions.
Definition Numerics.hh:174

◆ BesselK1L()

double Garfield::Numerics::BesselK1L ( const double xx)
inline

Definition at line 216 of file Numerics.hh.

216 {
217 const double y = 2. / xx;
218 return (exp(-xx) / sqrt(xx)) *
219 (1.25331414 + 0.23498619 * y - 0.03655620 * y * y +
220 0.01504268 * pow(y, 3) - 0.00780353 * pow(y, 4) +
221 0.00325614 * pow(y, 5) - 0.00068245 * pow(y, 6));
222}

◆ BesselK1S()

double Garfield::Numerics::BesselK1S ( const double xx)
inline

Definition at line 207 of file Numerics.hh.

207 {
208 const double y = xx / 2.;
209 const double y2 = y * y;
210 return log(y) * BesselI1S(xx) +
211 (1. / xx) * (1. + 0.15443144 * y2 - 0.67278579 * y2 * y2 -
212 0.18156897 * pow(y2, 3) - 0.01919402 * pow(y2, 4) -
213 0.00110404 * pow(y2, 5) - 0.00004686 * pow(y2, 6));
214}
double BesselI1S(const double xx)
Definition Numerics.hh:182

◆ Boxin2()

bool Garfield::Numerics::Boxin2 ( const std::vector< std::vector< double > > & value,
const std::vector< double > & xAxis,
const std::vector< double > & yAxis,
const int nx,
const int ny,
const double xx,
const double yy,
double & f,
const int iOrder )

Interpolation of order 1 and 2 in an irregular rectangular two-dimensional grid.

◆ Boxin3()

bool Garfield::Numerics::Boxin3 ( const std::vector< std::vector< std::vector< double > > > & value,
const std::vector< double > & xAxis,
const std::vector< double > & yAxis,
const std::vector< double > & zAxis,
const int nx,
const int ny,
const int nz,
const double xx,
const double yy,
const double zz,
double & f,
const int iOrder )

Interpolation of order 1 and 2 in an irregular rectangular three-dimensional grid.

◆ Divdif()

double Garfield::Numerics::Divdif ( const std::vector< double > & f,
const std::vector< double > & a,
int nn,
double x,
int mm )

C++ version of DIVDIF (CERN program library E105) which performs tabular interpolation using symmetrically placed argument points.

◆ GaussLegendreNodes2()

std::array< double, 2 > Garfield::Numerics::GaussLegendreNodes2 ( )
constexpr

Definition at line 16 of file Numerics.hh.

16 {
17 return {-0.577350269189625765, 0.577350269189625765};
18}

◆ GaussLegendreNodes3()

std::array< double, 3 > Garfield::Numerics::GaussLegendreNodes3 ( )
constexpr

Definition at line 20 of file Numerics.hh.

20 {
21 return {-0.774596669241483377, 0., 0.774596669241483377};
22}

◆ GaussLegendreNodes4()

std::array< double, 4 > Garfield::Numerics::GaussLegendreNodes4 ( )
constexpr

Definition at line 26 of file Numerics.hh.

26 {
27 return {-0.861136311594052575, -0.339981043584856265, 0.339981043584856265,
28 0.861136311594052575};
29}

◆ GaussLegendreNodes5()

std::array< double, 5 > Garfield::Numerics::GaussLegendreNodes5 ( )
constexpr

Definition at line 34 of file Numerics.hh.

34 {
35 return {-0.906179845938663993, -0.538469310105683091, 0.,
36 0.538469310105683091, 0.906179845938663993};
37}

◆ GaussLegendreNodes6()

std::array< double, 6 > Garfield::Numerics::GaussLegendreNodes6 ( )
constexpr

Definition at line 42 of file Numerics.hh.

42 {
43 return {-0.932469514203152028, -0.661209386466264514, -0.238619186083196909,
44 0.238619186083196909, 0.661209386466264514, 0.932469514203152028};
45}

◆ GaussLegendreWeights2()

std::array< double, 2 > Garfield::Numerics::GaussLegendreWeights2 ( )
constexpr

Definition at line 19 of file Numerics.hh.

19{ return {1., 1.}; }

◆ GaussLegendreWeights3()

std::array< double, 3 > Garfield::Numerics::GaussLegendreWeights3 ( )
constexpr

Definition at line 23 of file Numerics.hh.

23 {
24 return {5. / 9., 8. / 9., 5. / 9.};
25}

◆ GaussLegendreWeights4()

std::array< double, 4 > Garfield::Numerics::GaussLegendreWeights4 ( )
constexpr

Definition at line 30 of file Numerics.hh.

30 {
31 return {0.347854845137453857, 0.652145154862546143, 0.652145154862546143,
32 0.347854845137453857};
33}

◆ GaussLegendreWeights5()

std::array< double, 5 > Garfield::Numerics::GaussLegendreWeights5 ( )
constexpr

Definition at line 38 of file Numerics.hh.

38 {
39 return {0.236926885056189088, 0.478628670499366468, 128. / 225.,
40 0.478628670499366468, 0.236926885056189088};
41}

◆ GaussLegendreWeights6()

std::array< double, 6 > Garfield::Numerics::GaussLegendreWeights6 ( )
constexpr

Definition at line 46 of file Numerics.hh.

46 {
47 return {0.171324492379170345, 0.360761573048138608, 0.467913934572691047,
48 0.467913934572691047, 0.360761573048138608, 0.171324492379170345};
49}

◆ LeastSquaresFit()

bool Garfield::Numerics::LeastSquaresFit ( std::function< double(double, const std::vector< double > &)> f,
std::vector< double > & par,
std::vector< double > & epar,
const std::vector< double > & x,
const std::vector< double > & y,
const std::vector< double > & ey,
const unsigned int nMaxIter,
const double diff,
double & chi2,
const double eps,
const bool debug,
const bool verbose )

Least-squares minimisation.

◆ Legendre()

double Garfield::Numerics::Legendre ( const unsigned int n,
const double x )
inline

Legendre polynomials.

Definition at line 159 of file Numerics.hh.

159 {
160 if (std::abs(x) > 1.) return 0.;
161 double p0 = 1.;
162 double p1 = x;
163 if (n == 0) return p0;
164 if (n == 1) return p1;
165 for (unsigned int k = 1; k < n; ++k) {
166 p0 = ((2 * k + 1) * x * p1 - k * p0) / (k + 1);
167 std::swap(p0, p1);
168 }
169 return p1;
170}

◆ LinearInterpolation()

double Garfield::Numerics::LinearInterpolation ( const std::vector< double > & ytab,
const std::vector< double > & xtab,
const double x )