17 return {-0.577350269189625765, 0.577350269189625765};
23 return {-0.774596669241483377, 0., 0.774596669241483377};
26 return {5. / 9., 8. / 9., 5. / 9.};
29 return {-0.861136311594052575, -0.339981043584856265,
30 0.339981043584856265, 0.861136311594052575};
33 return {0.347854845137453857, 0.652145154862546143,
34 0.652145154862546143, 0.347854845137453857};
37 return {-0.906179845938663993, -0.538469310105683091, 0.,
38 0.538469310105683091, 0.906179845938663993};
41 return {0.236926885056189088, 0.478628670499366468, 128. / 225.,
42 0.478628670499366468, 0.236926885056189088};
45 return {-0.932469514203152028, -0.661209386466264514,
46 -0.238619186083196909, 0.238619186083196909,
47 0.661209386466264514, 0.932469514203152028};
50 return {0.171324492379170345, 0.360761573048138608,
51 0.467913934572691047, 0.467913934572691047,
52 0.360761573048138608, 0.171324492379170345};
101void qagi(std::function<
double(
double)> f,
double bound,
const int inf,
102 const double epsabs,
const double epsrel,
103 double& result,
double& abserr,
unsigned int& status);
116void qk15i(std::function<
double(
double)> f,
double bound,
const int inf,
117 const double a,
const double b,
double& result,
double& abserr,
118 double& resabs,
double& resasc);
121void qk15(std::function<
double(
double)> f,
const double a,
const double b,
122 double& result,
double& abserr,
double& resabs,
double& resasc);
134int deqn(
const int n, std::vector<std::vector<double> >& a,
135 std::vector<double>& b);
137int deqinv(
const int n, std::vector<std::vector<double> >& a,
138 std::vector<double>& b);
140void dfact(
const int n, std::vector<std::vector<double> >& a,
141 std::vector<int>& ir,
int& ifail,
double& det,
int& jfail);
142void dfeqn(
const int n, std::vector<std::vector<double> >& a,
143 std::vector<int>& ir, std::vector<double>& b);
144void dfinv(
const int n, std::vector<std::vector<double> >& a,
145 std::vector<int>& ir);
147int dinv(
const int n, std::vector<std::vector<double> >& a);
149void cfact(
const int n, std::vector<std::vector<std::complex<double> > >& a,
150 std::vector<int>& ir,
int& ifail, std::complex<double>& det,
152void cfinv(
const int n, std::vector<std::vector<std::complex<double> > >& a,
153 std::vector<int>& ir);
156int cinv(
const int n, std::vector<std::vector<std::complex<double> > >& a);
158void cfft(std::vector<std::complex<double> >& a,
const int msign);
162inline double Legendre(
const unsigned int n,
const double x) {
164 if (std::abs(x) > 1.)
return 0.;
167 if (n == 0)
return p0;
168 if (n == 1)
return p1;
169 for (
unsigned int k = 1; k < n; ++k) {
170 p0 = ((2 * k + 1) * x * p1 - k * p0) / (k + 1);
162inline double Legendre(
const unsigned int n,
const double x) {
…}
179 const double y = xx / 3.75;
180 const double y2 = y * y;
181 return 1. + 3.5156229 * y2 + 3.0899424 * y2 * y2 +
182 1.2067492 * pow(y2, 3) + 0.2659732 * pow(y2, 4) +
183 0.0360768 * pow(y2, 5) + 0.0045813 * pow(y2, 6);
187 const double y = xx / 3.75;
188 const double y2 = y * y;
190 (0.5 + 0.87890594 * y2 + 0.51498869 * y2 * y2 +
191 0.15084934 * pow(y2, 3) + 0.02658733 * pow(y2, 4) +
192 0.00301532 * pow(y2, 5) + 0.00032411 * pow(y2, 6));
196 const double y = xx / 2.;
197 const double y2 = y * y;
199 0.42278420 * y2 + 0.23069756 * y2 * y2 +
200 0.03488590 * pow(y2, 3) + 0.00262698 * pow(y2, 4) +
201 0.00010750 * pow(y2, 5) + 0.00000740 * pow(y2, 6);
205 const double y = 2. / xx;
206 return (exp(-xx) / sqrt(xx)) *
207 (1.25331414 - 0.07832358 * y + 0.02189568 * y * y -
208 0.01062446 * pow(y, 3) + 0.00587872 * pow(y, 4) -
209 0.00251540 * pow(y, 5) + 0.00053208 * pow(y, 6));
213 const double y = xx / 2.;
214 const double y2 = y * y;
217 (1. + 0.15443144 * y2 - 0.67278579 * y2 * y2 -
218 0.18156897 * pow(y2, 3) - 0.01919402 * pow(y2, 4) -
219 0.00110404 * pow(y2, 5) - 0.00004686 * pow(y2, 6));
223 const double y = 2. / xx;
224 return (exp(-xx) / sqrt(xx)) *
225 (1.25331414 + 0.23498619 * y - 0.03655620 * y * y +
226 0.01504268 * pow(y, 3) - 0.00780353 * pow(y, 4) +
227 0.00325614 * pow(y, 5) - 0.00068245 * pow(y, 6));
232double Divdif(
const std::vector<double>& f,
const std::vector<double>& a,
233 int nn,
double x,
int mm);
237bool Boxin2(
const std::vector<std::vector<double> >& value,
238 const std::vector<double>& xAxis,
const std::vector<double>& yAxis,
239 const int nx,
const int ny,
const double xx,
const double yy,
240 double& f,
const int iOrder);
243bool Boxin3(
const std::vector<std::vector<std::vector<double> > >& value,
244 const std::vector<double>& xAxis,
const std::vector<double>& yAxis,
245 const std::vector<double>& zAxis,
const int nx,
const int ny,
246 const int nz,
const double xx,
const double yy,
const double zz,
247 double& f,
const int iOrder);
251 std::function<
double(
double,
const std::vector<double>&)> f,
252 std::vector<double>& par, std::vector<double>& epar,
253 const std::vector<double>& x,
const std::vector<double>& y,
254 const std::vector<double>& ey,
const unsigned int nMaxIter,
255 const double diff,
double& chi2,
const double eps,
256 const bool debug,
const bool verbose);
Linear algebra routines from CERNLIB.
void dfinv(const int n, std::vector< std::vector< double > > &a, std::vector< int > &ir)
void cfact(const int n, std::vector< std::vector< std::complex< double > > > &a, std::vector< int > &ir, int &ifail, std::complex< double > &det, int &jfail)
int deqn(const int n, std::vector< std::vector< double > > &a, std::vector< double > &b)
void dfact(const int n, std::vector< std::vector< double > > &a, std::vector< int > &ir, int &ifail, double &det, int &jfail)
void cfinv(const int n, std::vector< std::vector< std::complex< double > > > &a, std::vector< int > &ir)
int dinv(const int n, std::vector< std::vector< double > > &a)
Replace square matrix A by its inverse.
void dfeqn(const int n, std::vector< std::vector< double > > &a, std::vector< int > &ir, std::vector< double > &b)
int cinv(const int n, std::vector< std::vector< std::complex< double > > > &a)
Replace square matrix A by its inverse.
int deqinv(const int n, std::vector< std::vector< double > > &a, std::vector< double > &b)
Replaces b by the solution x of Ax = b, and replace A by its inverse.
void cfft(std::vector< std::complex< double > > &a, const int msign)
void qk15(std::function< double(double)> f, const double a, const double b, double &result, double &abserr, double &resabs, double &resasc)
15-point Gauss-Kronrod integration with finite integration range.
void qagi(std::function< double(double)> f, double bound, const int inf, const double epsabs, const double epsrel, double &result, double &abserr, unsigned int &status)
void qk15i(std::function< double(double)> f, double bound, const int inf, const double a, const double b, double &result, double &abserr, double &resabs, double &resasc)
Collection of numerical routines.
constexpr std::array< double, 4 > GaussLegendreWeights4()
constexpr std::array< double, 2 > GaussLegendreNodes2()
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)
double BesselK0L(const double xx)
double BesselI0S(const double xx)
double BesselI1S(const double xx)
constexpr std::array< double, 5 > GaussLegendreWeights5()
double Legendre(const unsigned int n, const double x)
Legendre polynomials.
constexpr std::array< double, 6 > GaussLegendreWeights6()
double BesselK1S(const double xx)
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.
constexpr std::array< double, 2 > GaussLegendreWeights2()
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)
constexpr std::array< double, 3 > GaussLegendreNodes3()
constexpr std::array< double, 4 > GaussLegendreNodes4()
double Divdif(const std::vector< double > &f, const std::vector< double > &a, int nn, double x, int mm)
constexpr std::array< double, 5 > GaussLegendreNodes5()
constexpr std::array< double, 3 > GaussLegendreWeights3()
double BesselK0S(const double xx)
double BesselK1L(const double xx)
constexpr std::array< double, 6 > GaussLegendreNodes6()