10 #ifndef SRC_INCLUDE_INTEGRATE_H_ 11 #define SRC_INCLUDE_INTEGRATE_H_ 14 #include <gsl/gsl_integration.h> 15 #include <gsl/gsl_monte_plain.h> 16 #include <gsl/gsl_monte_vegas.h> 42 void operator()(gsl_integration_cquad_workspace *ptr)
const {
46 gsl_integration_cquad_workspace_free(ptr);
54 class Result :
public std::pair<double, double> {
56 using Base = std::pair<double, double>;
63 operator double()
const {
return Base::first; }
66 double value()
const {
return Base::first; }
69 double error()
const {
return Base::second; }
80 double relative_tolerance = 5e-4,
81 double absolute_tolerance = 1e-9)
const {
82 const double allowed_error =
83 std::max(absolute_tolerance, value() * relative_tolerance);
84 if (error() > allowed_error) {
85 std::stringstream error_msg;
86 error_msg << integration_name <<
" resulted in I = " << value() <<
" ± " 88 <<
", but the required precision is either absolute error < " 89 << absolute_tolerance <<
" or relative error < " 90 << relative_tolerance << std::endl;
91 throw std::runtime_error(error_msg.str());
123 : workspace_(gsl_integration_cquad_workspace_alloc(workspace_size)) {}
138 template <
typename F>
141 const gsl_function gslfun{
144 [](
double x,
void *type_erased) ->
double {
145 auto &&f = *
static_cast<F *
>(type_erased);
151 const int error_code = gsl_integration_cquad(
152 &gslfun, a, b, accuracy_absolute_, accuracy_relative_, workspace_.get(),
153 &result.first, &result.second,
156 std::stringstream err;
157 err <<
"GSL 1D deterministic integration: " << gsl_strerror(error_code);
158 throw std::runtime_error(err.str());
160 result.
check_error(
"GSL 1D deterministic integration", accuracy_relative_,
167 std::unique_ptr<gsl_integration_cquad_workspace, GslWorkspaceDeleter>
171 const double accuracy_absolute_ = 1.0e-9;
174 const double accuracy_relative_ = 5.0e-4;
205 : state_(gsl_monte_plain_alloc(1)),
206 rng_(gsl_rng_alloc(gsl_rng_mt19937)),
207 number_of_calls_(num_calls) {
208 gsl_monte_plain_init(state_);
211 gsl_rng_set(rng_, seed);
216 gsl_monte_plain_free(state_);
234 template <
typename F>
238 const double lower[1] = {min};
239 const double upper[1] = {max};
244 const gsl_monte_function monte_fun{
246 [](
double *x,
size_t ,
void *params) ->
double {
247 auto &&f = *
static_cast<F *
>(params);
252 const int error_code =
253 gsl_monte_plain_integrate(&monte_fun, lower, upper, 1, number_of_calls_,
254 rng_, state_, &result.first, &result.second);
256 std::stringstream err;
257 err <<
"GSL 1D Monte-Carlo integration: " << gsl_strerror(error_code);
258 throw std::runtime_error(err.str());
261 result.
check_error(
"GSL 1D Monte-Carlo integration");
305 : state_(gsl_monte_plain_alloc(2)),
306 rng_(gsl_rng_alloc(gsl_rng_mt19937)),
307 number_of_calls_(num_calls) {
308 gsl_monte_plain_init(state_);
311 gsl_rng_set(rng_, seed);
316 gsl_monte_plain_free(state_);
336 template <
typename F>
341 const double lower[2] = {min1, min2};
342 const double upper[2] = {max1, max2};
344 if (max1 <= min1 || max2 <= min2)
347 const gsl_monte_function monte_fun{
349 [](
double *x,
size_t ,
void *params) ->
double {
350 auto &&f = *
static_cast<F *
>(params);
351 return f(x[0], x[1]);
355 gsl_monte_plain_integrate(&monte_fun, lower, upper, 2, number_of_calls_,
356 rng_, state_, &result.first, &result.second);
378 template <
typename F>
420 double epsabs = 1e-9)
421 : maxeval_(num_calls), epsrel_(epsrel), epsabs_(epsabs) {}
439 template <
typename F>
443 if (max1 < min1 || max2 < min2) {
446 bool tolerable = (max1 - min1 > -1.e-15) && (max2 - min2 > -1.e-15);
450 std::stringstream err;
451 err <<
"Integrator2dCuhre got wrong integration limits: [" << min1 <<
", " 452 << max1 <<
"], [" << min2 <<
", " << max2 <<
"]";
453 throw std::invalid_argument(err.str());
456 Integrand2d<F> f_with_limits = {min1, max1 - min1, min2, max2 - min2, fun};
458 const integrand_t cuhre_fun{[](
const int * ,
const cubareal xx[],
459 const int * , cubareal ff[],
460 void *userdata) ->
int {
464 ff[0] = (i->f)(i->min1 + i->diff1 * xx[0], i->min2 + i->diff2 * xx[1]) *
471 void *userdata = &f_with_limits;
474 const int mineval = 0;
475 const int maxeval = maxeval_;
477 const char *statefile =
nullptr;
478 void *spin =
nullptr;
480 Cuhre(ndim, ncomp, cuhre_fun, userdata, nvec, epsrel_, epsabs_, flags,
481 mineval, maxeval, key, statefile, spin, &nregions_, &neval_, &fail_,
482 &result.first, &result.second, &prob_);
485 std::stringstream err;
486 err <<
"After " << neval_ <<
" evaluations " 487 <<
"Cuhre integration from Cuba reports error code " << fail_;
488 throw std::runtime_error(err.str());
490 result.
check_error(
"Cuba integration ", epsrel_, epsabs_);
523 #endif // SRC_INCLUDE_INTEGRATE_H_ int maxeval_
The (approximate) maximum number of integrand evaluations allowed.
double min2
the lower bound of the second integrated variable
void check_error(const std::string &integration_name, double relative_tolerance=5e-4, double absolute_tolerance=1e-9) const
Check whether the error is small and alert if it is not.
gsl_monte_plain_state * state_
internal state of the Monte-Carlo integrator
This is a wrapper for the integrand, so we can pass the limits as well for renormalizing to the unit ...
Result operator()(double a, double b, F &&fun)
The function call operator implements the integration functionality.
constexpr GslWorkspaceDeleter()=default
The class has no members, so this is a noop.
A deleter type for std::unique_ptr to be used with gsl_integration_workspace pointers.
~Integrator1dMonte()
Destructor: Clean up internal state and RNG.
Guard type that safely disables floating point traps for the scope in which it is placed...
double value() const
Access the first entry in the pair as the value.
const std::size_t number_of_calls_
number of calls to the integrand
A C++ interface for numerical integration in two dimensions with the GSL Monte-Carlo integration func...
Integrator1dMonte(size_t num_calls=1E6)
Construct an integration functor.
std::pair< double, double > Base
The data type to store the value and the error of the integration.
std::unique_ptr< gsl_integration_cquad_workspace, GslWorkspaceDeleter > workspace_
Holds the workspace pointer.
double min1
the lower bound of the first integrated variable
int neval_
Actual number of integrand evaluations needed.
A C++ interface for numerical integration in two dimensions with the Cuba Cuhre integration function...
A C++ interface for numerical integration in one dimension with the GSL Monte-Carlo integration funct...
The result type returned from integrations, containing the value and an error.
double diff2
the integration range of the second integrated variable
T uniform_int(T min, T max)
Integrator2d(size_t num_calls=1E6)
Construct an integration functor.
Integrator2dCuhre(int num_calls=1e6, double epsrel=5e-4, double epsabs=1e-9)
Construct an integration functor.
gsl_monte_plain_state * state_
internal state of the Monte-Carlo integrator
double diff1
the integration range of the first integrated variable
gsl_rng * rng_
random number generator
Result operator()(double min1, double max1, double min2, double max2, F &&fun)
The function call operator implements the integration functionality.
double prob_
The chi^2 probability that the error is not a reliable estimate of the true integration error...
F f
the integrated function
double epsabs_
Requested absolute accuracy.
gsl_rng * rng_
random number generator
Result operator()(double min, double max, F &&fun)
The function call operator implements the integration functionality.
const std::size_t number_of_calls_
number of calls to the integrand
A C++ interface for numerical integration in one dimension with the GSL CQUAD integration functions...
Result operator()(double min1, double max1, double min2, double max2, F fun)
The function call operator implements the integration functionality.
int nregions_
Actual number of subregions needed.
~Integrator2d()
Destructor: Clean up internal state and RNG.
void operator()(gsl_integration_cquad_workspace *ptr) const
Frees the gsl_integration_cquad_workspace resource if it is non-zero.
Integrator(size_t workspace_size=1000)
Construct an integration functor with the given workspace_size.
double epsrel_
Requested relative accuracy.
double error() const
Access the second entry in the pair as the absolute error.