cepgen is hosted by Hepforge, IPPP Durham
CepGen 1.2.5
Central exclusive processes event generator
Loading...
Searching...
No Matches
NaiveBoostIntegrator.cpp
Go to the documentation of this file.
1/*
2 * CepGen: a central exclusive processes event generator
3 * Copyright (C) 2020-2024 Laurent Forthomme
4 *
5 * This program is free software: you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation, either version 3 of the License, or
8 * any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program. If not, see <http://www.gnu.org/licenses/>.
17 */
18
19#include <boost/math/quadrature/naive_monte_carlo.hpp>
20
24
25namespace cepgen {
27 class NaiveBoostIntegrator final : public Integrator {
28 public:
29 explicit NaiveBoostIntegrator(const ParametersList& params) : Integrator(params) {}
30
32 auto desc = Integrator::description();
33 desc.setDescription("\"Naive\" Boost integrator");
34 return desc;
35 }
36
37 void setLimits(const std::vector<Limits>& lims) override {
39 bounds_.clear();
40 std::transform(
41 limits_.begin(), limits_.end(), std::back_inserter(bounds_), [](const auto& lim) { return lim.raw(); });
42 }
43 Value integrate(Integrand& integrand) override {
44 checkLimits(integrand);
45
46 auto funct = [&](const std::vector<double>& coord) -> double { return integrand.eval(coord); };
47 boost::math::quadrature::naive_monte_carlo<double, decltype(funct)> mc(funct, bounds_, 1.e-2, true, 1);
48 auto task = mc.integrate();
49
50 return Value{task.get(), mc.current_error_estimate()};
51 }
52
53 private:
54 std::vector<std::pair<double, double> > bounds_;
55 };
56} // namespace cepgen
57
58REGISTER_INTEGRATOR("Naive", NaiveBoostIntegrator);
#define REGISTER_INTEGRATOR(name, obj)
Add a generic process definition to the list of handled processes.
An integrand wrapper placeholder.
Definition Integrand.h:27
virtual double eval(const std::vector< double > &)=0
Compute the integrand for a given coordinates set.
Monte-Carlo integration algorithm.
Definition Integrator.h:28
void checkLimits(const Integrand &)
Ensure the integration bounds are properly set.
virtual void setLimits(const std::vector< Limits > &limits)
Specify the variables limits on integration.
Definition Integrator.h:38
std::vector< Limits > limits_
List of per-variable integration limits.
Definition Integrator.h:58
static ParametersDescription description()
Boost's Naive integration algorithm.
Value integrate(Integrand &integrand) override
Perform the multidimensional Monte Carlo integration.
static ParametersDescription description()
void setLimits(const std::vector< Limits > &lims) override
Specify the variables limits on integration.
NaiveBoostIntegrator(const ParametersList &params)
A description object for parameters collection.
A scalar value with its uncertainty.
Definition Value.h:26
Common namespace for this Monte Carlo generator.