Version: SMASH-3.1
integrate.h
Go to the documentation of this file.
1 /*
2  *
3  * Copyright (c) 2015-2020,2022
4  * SMASH Team
5  *
6  * GNU General Public License (GPLv3 or later)
7  *
8  */
9 
10 #ifndef SRC_INCLUDE_SMASH_INTEGRATE_H_
11 #define SRC_INCLUDE_SMASH_INTEGRATE_H_
12 
13 #include <cuba.h>
14 
15 #include <algorithm>
16 #include <iostream>
17 #include <memory>
18 #include <sstream>
19 #include <string>
20 #include <tuple>
21 #include <utility>
22 
23 #include "gsl/gsl_integration.h"
24 #include "gsl/gsl_monte_plain.h"
25 #include "gsl/gsl_monte_vegas.h"
26 
27 #include "fpenvironment.h"
28 #include "random.h"
29 
30 namespace smash {
31 
39  constexpr GslWorkspaceDeleter() = default;
40 
42  void operator()(gsl_integration_cquad_workspace *ptr) const {
43  if (ptr == nullptr) {
44  return;
45  }
46  gsl_integration_cquad_workspace_free(ptr);
47  }
48 };
49 
54 class Result : public std::pair<double, double> {
56  using Base = std::pair<double, double>;
57 
58  public:
60  using Base::pair;
61 
63  operator double() const { return Base::first; }
64 
66  double value() const { return Base::first; }
67 
69  double error() const { return Base::second; }
70 
79  void check_error(const std::string &integration_name,
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() << " ± "
87  << error()
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());
92  }
93  }
94 };
95 
106 class Integrator {
107  public:
122  explicit Integrator(size_t workspace_size = 1000)
123  : workspace_(gsl_integration_cquad_workspace_alloc(workspace_size)) {}
124 
138  template <typename F>
139  Result operator()(double a, double b, F &&fun) {
140  Result result;
141  const gsl_function gslfun{
142  /* important! The lambda cannot use captures, otherwise the
143  * conversion to a function pointer type is impossible. */
144  [](double x, void *type_erased) -> double {
145  auto &&f = *static_cast<F *>(type_erased);
146  return f(x);
147  },
148  &fun};
149  // We disable float traps when calling GSL code we cannot control.
150  DisableFloatTraps guard;
151  const int error_code = gsl_integration_cquad(
152  &gslfun, a, b, accuracy_absolute_, accuracy_relative_, workspace_.get(),
153  &result.first, &result.second,
154  nullptr /* Don't store the number of evaluations */);
155  if (error_code) {
156  std::stringstream err;
157  err << "GSL 1D deterministic integration: " << gsl_strerror(error_code);
158  throw std::runtime_error(err.str());
159  }
160  result.check_error("GSL 1D deterministic integration", accuracy_relative_,
162  return result;
163  }
164 
168  void set_precision(double absolute, double relative) {
169  accuracy_absolute_ = absolute;
170  accuracy_relative_ = relative;
171  }
172 
173  private:
175  std::unique_ptr<gsl_integration_cquad_workspace, GslWorkspaceDeleter>
177 
179  double accuracy_absolute_ = 1.0e-9;
180 
182  double accuracy_relative_ = 5.0e-4;
183 };
184 
191 template <typename F>
192 struct Integrand2d {
194  double min1;
196  double diff1;
198  double min2;
200  double diff2;
202  F f;
203 };
204 
220  public:
232  explicit Integrator2d(int num_calls = 1e6, double epsrel = 5e-4,
233  double epsabs = 1e-9)
234  : maxeval_(num_calls), epsrel_(epsrel), epsabs_(epsabs) {}
235 
252  template <typename F>
253  Result operator()(double min1, double max1, double min2, double max2, F fun) {
254  Result result = {0., 0.};
255 
256  if (max1 < min1 || max2 < min2) {
257  /* Tolerance chosen large enough for floating-point arithmetic error in
258  * chained decays, consider increasing if insufficient. */
259  bool tolerable = (max1 - min1 > -1.e-15) && (max2 - min2 > -1.e-15);
260  if (tolerable) {
261  return result;
262  }
263  std::stringstream err;
264  err << "Integrator2d got wrong integration limits: [" << min1 << ", "
265  << max1 << "], [" << min2 << ", " << max2 << "]";
266  throw std::invalid_argument(err.str());
267  }
268 
269  Integrand2d<F> f_with_limits = {min1, max1 - min1, min2, max2 - min2, fun};
270 
271  const integrand_t cuhre_fun{[](const int * /* ndim */, const cubareal xx[],
272  const int * /* ncomp */, cubareal ff[],
273  void *userdata) -> int {
274  auto i = static_cast<Integrand2d<F> *>(userdata);
275  /* We have to transform the integrand to the unit cube.
276  * This is what Cuba expects. */
277  ff[0] = (i->f)(i->min1 + i->diff1 * xx[0], i->min2 + i->diff2 * xx[1]) *
278  i->diff1 * i->diff2;
279  return 0;
280  }};
281 
282  const int ndim = 2;
283  const int ncomp = 1;
284  void *userdata = &f_with_limits;
285  const int nvec = 1;
286  const int flags = 0; // Use the defaults.
287  const int mineval = 0;
288  const int maxeval = maxeval_;
289  const int key = -1; // Use the default.
290  const char *statefile = nullptr;
291  void *spin = nullptr;
292 
293  Cuhre(ndim, ncomp, cuhre_fun, userdata, nvec, epsrel_, epsabs_, flags,
294  mineval, maxeval, key, statefile, spin, &nregions_, &neval_, &fail_,
295  &result.first, &result.second, &prob_);
296 
297  if (fail_) {
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());
302  }
303  result.check_error("Cuba integration ", epsrel_, epsabs_);
304 
305  return result;
306  }
307 
308  private:
310  int maxeval_;
312  double epsrel_;
314  double epsabs_;
318  int neval_;
326  int fail_;
331  double prob_;
332 };
333 
334 } // namespace smash
335 
336 #endif // SRC_INCLUDE_SMASH_INTEGRATE_H_
Guard type that safely disables floating point traps for the scope in which it is placed.
Definition: fpenvironment.h:79
A C++ interface for numerical integration in two dimensions with the Cuba Cuhre integration function.
Definition: integrate.h:219
double prob_
The chi^2 probability that the error is not a reliable estimate of the true integration error.
Definition: integrate.h:331
double epsabs_
Requested absolute accuracy.
Definition: integrate.h:314
int fail_
An error flag.
Definition: integrate.h:326
int neval_
Actual number of integrand evaluations needed.
Definition: integrate.h:318
int nregions_
Actual number of subregions needed.
Definition: integrate.h:316
int maxeval_
The (approximate) maximum number of integrand evaluations allowed.
Definition: integrate.h:310
double epsrel_
Requested relative accuracy.
Definition: integrate.h:312
Integrator2d(int num_calls=1e6, double epsrel=5e-4, double epsabs=1e-9)
Construct an integration functor.
Definition: integrate.h:232
Result operator()(double min1, double max1, double min2, double max2, F fun)
The function call operator implements the integration functionality.
Definition: integrate.h:253
A C++ interface for numerical integration in one dimension with the GSL CQUAD integration functions.
Definition: integrate.h:106
Result operator()(double a, double b, F &&fun)
The function call operator implements the integration functionality.
Definition: integrate.h:139
double accuracy_relative_
Parameter to the GSL integration function: desired relative error limit.
Definition: integrate.h:182
void set_precision(double absolute, double relative)
Set precision for the absolute and relative error of the integration.
Definition: integrate.h:168
Integrator(size_t workspace_size=1000)
Construct an integration functor with the given workspace_size.
Definition: integrate.h:122
std::unique_ptr< gsl_integration_cquad_workspace, GslWorkspaceDeleter > workspace_
Holds the workspace pointer.
Definition: integrate.h:176
double accuracy_absolute_
Parameter to the GSL integration function: desired absolute error limit.
Definition: integrate.h:179
The result type returned from integrations, containing the value and an error.
Definition: integrate.h:54
double value() const
Access the first entry in the pair as the value.
Definition: integrate.h:66
std::pair< double, double > Base
The data type to store the value and the error of the integration.
Definition: integrate.h:56
double error() const
Access the second entry in the pair as the absolute error.
Definition: integrate.h:69
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.
Definition: integrate.h:79
Definition: action.h:24
A deleter type for std::unique_ptr to be used with gsl_integration_workspace pointers.
Definition: integrate.h:37
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.
Definition: integrate.h:42
This is a wrapper for the integrand, so we can pass the limits as well for renormalizing to the unit ...
Definition: integrate.h:192
double diff1
the integration range of the first integrated variable
Definition: integrate.h:196
double min1
the lower bound of the first integrated variable
Definition: integrate.h:194
F f
the integrated function
Definition: integrate.h:202
double min2
the lower bound of the second integrated variable
Definition: integrate.h:198
double diff2
the integration range of the second integrated variable
Definition: integrate.h:200