cepgen is hosted by Hepforge, IPPP Durham
CepGen 1.2.5
Central exclusive processes event generator
Loading...
Searching...
No Matches
BoostTrapAnalyticalIntegrator.cpp
Go to the documentation of this file.
1/*
2 * CepGen: a central exclusive processes event generator
3 * Copyright (C) 2022-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/trapezoidal.hpp>
20
24
25namespace cepgen {
28 public:
30 : AnalyticIntegrator(params),
31 max_refinements_(steerAs<int, size_t>("limit")),
32 tol_(steer<double>("tolerance")) {}
33
36 desc.setDescription("Boost trapezoidal integration algorithm");
37 desc.add<int>("limit", 1000).setDescription("maximum number of subintervals to build");
38 desc.add<double>("tolerance", 1.e-6).setDescription("maximal tolerance");
39 return desc;
40 }
41
42 double integrate(const utils::Function1D& func, void* = nullptr, const Limits& lim = {}) const override {
43 const double xmin = (lim.hasMin() ? lim.min() : range_.min());
44 const double xmax = (lim.hasMax() ? lim.max() : range_.max());
45 return boost::math::quadrature::trapezoidal(func, xmin, xmax, tol_, max_refinements_);
46 }
47
48 private:
49 const size_t max_refinements_;
50 const double tol_;
51 };
52} // namespace cepgen
#define REGISTER_ANALYTIC_INTEGRATOR(name, obj)
Add a generic analytical integrator object builder definition.
Analytic (functional) integration algorithm.
static ParametersDescription description()
BoostAnalyticalIntegrator(const ParametersList &params)
double integrate(const utils::Function1D &func, void *=nullptr, const Limits &lim={}) const override
Evaluate the integral of a function at a given value.
Validity interval for a variable.
Definition Limits.h:28
double min() const
Lower limit to apply on the variable.
Definition Limits.h:52
double max() const
Upper limit to apply on the variable.
Definition Limits.h:54
A description object for parameters collection.
U steerAs(const std::string &key) const
Retrieve a recasted parameters as previously steered.
Definition Steerable.h:44
T steer(const std::string &key) const
Retrieve a parameters as previously steered.
Definition Steerable.h:39
Wrapper to a 1-dimensional function with optional parameters.
Common namespace for this Monte Carlo generator.