#include <integrate.h>
A C++ interface for numerical integration in one dimension with the GSL Monte-Carlo integration functions.
Example:
[](double x) { return x * x; });
Definition at line 188 of file integrate.h.
◆ Integrator1dMonte()
smash::Integrator1dMonte::Integrator1dMonte |
( |
size_t |
num_calls = 1E6 | ) |
|
|
inlineexplicit |
Construct an integration functor.
- Parameters
-
[in] | num_calls | The desired number of calls to the integrand function (defaults to 1E6 if omitted), i.e. how often the integrand is sampled in the Monte-Carlo integration. Larger numbers lead to a more precise result, but also to increased runtime. |
- Note
- Since the workspace is allocated in the constructor and deallocated on destruction, you should not recreate Integrator objects unless required. Thus, if you want to calculate multiple integrals with the same
workspace_size
, keep the Integrator object around.
Definition at line 204 of file integrate.h.
205 :
state_(gsl_monte_plain_alloc(1)),
206 rng_(gsl_rng_alloc(gsl_rng_mt19937)),
208 gsl_monte_plain_init(
state_);
211 gsl_rng_set(
rng_, seed);
◆ ~Integrator1dMonte()
smash::Integrator1dMonte::~Integrator1dMonte |
( |
| ) |
|
|
inline |
Destructor: Clean up internal state and RNG.
Definition at line 215 of file integrate.h.
216 gsl_monte_plain_free(
state_);
◆ operator()()
template<typename F >
Result smash::Integrator1dMonte::operator() |
( |
double |
min, |
|
|
double |
max, |
|
|
F && |
fun |
|
) |
| |
|
inline |
The function call operator implements the integration functionality.
- Parameters
-
[in] | min | The lower limit of the integration. |
[in] | max | The upper limit of the integration. |
- Template Parameters
-
F | Type of the integrand function. |
- Parameters
-
[in] | fun | The callable to integrate over. This callable may be a function pointer, lambda, or a functor object. In any case, the callable must return a double and take two double arguments. If you want to pass additional data to the callable you can e.g. use lambda captures. |
- Returns
- Pair of integral value and absolute error estimate.
Definition at line 235 of file integrate.h.
236 Result result = {0, 0};
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 =
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");
◆ state_
gsl_monte_plain_state* smash::Integrator1dMonte::state_ |
|
private |
internal state of the Monte-Carlo integrator
Definition at line 268 of file integrate.h.
◆ rng_
gsl_rng* smash::Integrator1dMonte::rng_ |
|
private |
◆ number_of_calls_
const std::size_t smash::Integrator1dMonte::number_of_calls_ |
|
private |
number of calls to the integrand
Definition at line 274 of file integrate.h.
The documentation for this class was generated from the following file: