10 #ifndef SRC_INCLUDE_SMASH_INTEGRATE_H_ 
   11 #define SRC_INCLUDE_SMASH_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(
 
  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());
 
  175   std::unique_ptr<gsl_integration_cquad_workspace, GslWorkspaceDeleter>
 
  191 template <
typename F>
 
  233                         double epsabs = 1e-9)
 
  252   template <
typename F>
 
  256     if (max1 < min1 || max2 < min2) {
 
  259       bool tolerable = (max1 - min1 > -1.e-15) && (max2 - min2 > -1.e-15);
 
  263       std::stringstream err;
 
  264       err << 
"Integrator2d got wrong integration limits: [" << min1 << 
", " 
  265           << max1 << 
"], [" << min2 << 
", " << max2 << 
"]";
 
  266       throw std::invalid_argument(err.str());
 
  269     Integrand2d<F> f_with_limits = {min1, max1 - min1, min2, max2 - min2, fun};
 
  271     const integrand_t cuhre_fun{[](
const int * , 
const cubareal xx[],
 
  272                                    const int * , cubareal ff[],
 
  273                                    void *userdata) -> 
int {
 
  277       ff[0] = (i->f)(i->min1 + i->diff1 * xx[0], i->min2 + i->diff2 * xx[1]) *
 
  284     void *userdata = &f_with_limits;
 
  287     const int mineval = 0;
 
  290     const char *statefile = 
nullptr;
 
  291     void *spin = 
nullptr;
 
  293     Cuhre(ndim, ncomp, cuhre_fun, userdata, nvec, 
epsrel_, 
epsabs_, flags,
 
  295           &result.first, &result.second, &
prob_);
 
  298       std::stringstream err;
 
  299       err << 
"After " << 
neval_ << 
" evaluations " 
  300           << 
"Cuhre integration from Cuba reports error code " << 
fail_;
 
  301       throw std::runtime_error(err.str());
 
Guard type that safely disables floating point traps for the scope in which it is placed.
 
A C++ interface for numerical integration in two dimensions with the Cuba Cuhre integration function.
 
double prob_
The chi^2 probability that the error is not a reliable estimate of the true integration error.
 
double epsabs_
Requested absolute accuracy.
 
int neval_
Actual number of integrand evaluations needed.
 
int nregions_
Actual number of subregions needed.
 
int maxeval_
The (approximate) maximum number of integrand evaluations allowed.
 
double epsrel_
Requested relative accuracy.
 
Integrator2d(int num_calls=1e6, double epsrel=5e-4, double epsabs=1e-9)
Construct an integration functor.
 
Result operator()(double min1, double max1, double min2, double max2, F fun)
The function call operator implements the integration functionality.
 
A C++ interface for numerical integration in one dimension with the GSL CQUAD integration functions.
 
Result operator()(double a, double b, F &&fun)
The function call operator implements the integration functionality.
 
double accuracy_relative_
Parameter to the GSL integration function: desired relative error limit.
 
void set_precision(double absolute, double relative)
Set precision for the absolute and relative error of the integration.
 
Integrator(size_t workspace_size=1000)
Construct an integration functor with the given workspace_size.
 
std::unique_ptr< gsl_integration_cquad_workspace, GslWorkspaceDeleter > workspace_
Holds the workspace pointer.
 
double accuracy_absolute_
Parameter to the GSL integration function: desired absolute error limit.
 
The result type returned from integrations, containing the value and an error.
 
double value() const
Access the first entry in the pair as the value.
 
std::pair< double, double > Base
The data type to store the value and the error of the integration.
 
double error() const
Access the second entry in the pair as the absolute error.
 
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.
 
A deleter type for std::unique_ptr to be used with gsl_integration_workspace pointers.
 
constexpr GslWorkspaceDeleter()=default
The class has no members, so this is a noop.
 
void operator()(gsl_integration_cquad_workspace *ptr) const
Frees the gsl_integration_cquad_workspace resource if it is non-zero.
 
This is a wrapper for the integrand, so we can pass the limits as well for renormalizing to the unit ...
 
double diff1
the integration range of the first integrated variable
 
double min1
the lower bound of the first integrated variable
 
F f
the integrated function
 
double min2
the lower bound of the second integrated variable
 
double diff2
the integration range of the second integrated variable