cepgen is hosted by Hepforge, IPPP Durham
CepGen 1.2.5
Central exclusive processes event generator
Loading...
Searching...
No Matches
Integrator.h
Go to the documentation of this file.
1/*
2 * CepGen: a central exclusive processes event generator
3 * Copyright (C) 2013-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#ifndef CepGen_Integration_Integrator_h
20#define CepGen_Integration_Integrator_h
21
23#include "CepGen/Utils/Value.h"
24
25namespace cepgen {
26 class Integrand;
28 class Integrator : public NamedModule<Integrator> {
29 public:
31 explicit Integrator(const ParametersList& params);
32
34
36 void checkLimits(const Integrand&);
38 virtual void setLimits(const std::vector<Limits>& limits) { limits_ = limits; }
39
41 virtual double eval(Integrand&, const std::vector<double>&) const;
43 virtual double uniform(const Limits& = {0., 1.}) const;
44
47 virtual Value integrate(Integrand& result) = 0;
49 static Value integrate(const std::function<double(const std::vector<double>&)>&, const ParametersList&, size_t);
51 static Value integrate(const std::function<double(const std::vector<double>&)>&,
52 const ParametersList&,
53 const std::vector<Limits>&);
54
55 protected:
56 const std::unique_ptr<utils::RandomGenerator> rnd_gen_;
58 std::vector<Limits> limits_;
59 };
60} // namespace cepgen
61
62#endif
An integrand wrapper placeholder.
Definition Integrand.h:27
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
int verbosity_
Integrator verbosity.
Definition Integrator.h:57
virtual double uniform(const Limits &={0., 1.}) const
Generate a uniformly distributed (between 0 and 1) random number.
virtual double eval(Integrand &, const std::vector< double > &) const
Compute the function value at the given phase space point.
std::vector< Limits > limits_
List of per-variable integration limits.
Definition Integrator.h:58
const std::unique_ptr< utils::RandomGenerator > rnd_gen_
Definition Integrator.h:56
static ParametersDescription description()
virtual Value integrate(Integrand &result)=0
Perform the multidimensional Monte Carlo integration.
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.