cepgen is hosted by Hepforge, IPPP Durham
CepGen N/A
Central exclusive processes event generator
Integrator.h
Go to the documentation of this file.
1/*
2 * CepGen: a central exclusive processes event generator
3 * Copyright (C) 2022-2025 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#ifndef CepGen_Integration_Integrator_h
20#define CepGen_Integration_Integrator_h
21
22#include <functional>
23
25#include "CepGen/Utils/Limits.h"
26#include "CepGen/Utils/Value.h"
27
28namespace cepgen {
29 class Integrand;
30} // namespace cepgen
31
32namespace cepgen {
34 class Integrator : public NamedModule<Integrator> {
35 public:
36 explicit Integrator(const ParametersList&);
37
39
40 virtual bool oneDimensional() const { return false; }
41 virtual double eval(Integrand&, const std::vector<double>&) const;
42
44 Value integrate(Integrand& integrand, const std::vector<Limits>& = {});
48 Value integrate(const std::function<double(double)>& integrand, const Limits& range_1d = {0., 1.});
52 Value integrate(const std::function<double(const std::vector<double>&)>& integrand,
53 const std::vector<Limits>& range);
54
55 protected:
59 virtual Value run(Integrand& integrand, const std::vector<Limits>& range) = 0;
60
61 const int verbosity_;
62 };
63} // namespace cepgen
64
65#endif
An integrand wrapper placeholder.
Definition Integrand.h:26
Integration algorithm.
Definition Integrator.h:34
Value integrate(Integrand &integrand, const std::vector< Limits > &={})
Evaluate the integral for a given range.
Value integrate(const std::function< double(double)> &integrand, const Limits &range_1d={0., 1.})
Evaluate the integral of a function for a given range.
virtual double eval(Integrand &, const std::vector< double > &) const
Compute function value at one point.
Integrator(const ParametersList &)
virtual Value run(Integrand &integrand, const std::vector< Limits > &range)=0
Evaluate the integral of a function for a given range.
Value integrate(const std::function< double(const std::vector< double > &)> &integrand, const std::vector< Limits > &range)
Evaluate the integral of a function for a given range.
static ParametersDescription description()
virtual bool oneDimensional() const
Is the integrator designed for one-dimensional case?
Definition Integrator.h:40
const int verbosity_
Integrator verbosity.
Definition Integrator.h:61
Validity interval for a variable.
Definition Limits.h:28
Base runtime module object.
Definition NamedModule.h:28
A description object for parameters collection.
A scalar value with its uncertainty.
Definition Value.h:26
Common namespace for this Monte Carlo generator.
Definition Handler.h:26