cepgen is hosted by Hepforge, IPPP Durham
CepGen 1.2.5
Central exclusive processes event generator
Loading...
Searching...
No Matches
Integrator.cpp
Go to the documentation of this file.
1/*
2 * CepGen: a central exclusive processes event generator
3 * Copyright (C) 2021-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
22
23namespace cepgen {
24 namespace cuba {
25 Integrand* Integrator::gIntegrand = nullptr;
26
28 : cepgen::Integrator(params),
29 ncomp_(steer<int>("ncomp")),
30 nvec_(steer<int>("nvec")),
31 epsrel_(steer<double>("epsrel")),
32 epsabs_(steer<double>("epsabs")),
33 mineval_(steer<int>("mineval")),
34 maxeval_(steer<int>("maxeval")),
35 verbose_(steer<int>("verbose")) {}
36
38 gIntegrand = &integr;
39 return integrate();
40 }
41
44 desc.setDescription("Cuba generic integration algorithm");
45 desc.add<int>("ncomp", 1).setDescription("number of components of the integrand");
46 desc.add<int>("nvec", 1).setDescription("number of samples received by the integrand");
47 desc.add<double>("epsrel", 1.e-3).setDescription("requested relative accuracy");
48 desc.add<double>("epsabs", 1.e-12).setDescription("requested absolute accuracy");
49 desc.add<int>("mineval", 0).setDescription("minimum number of integrand evaluations required");
50 desc.add<int>("maxeval", 50'000).setDescription("(approximate) maximum number of integrand evaluations allowed");
51 desc.add<int>("verbose", 0);
52 return desc;
53 }
54
55 int cuba_integrand(const int* ndim, const double xx[], const int* /*ncomp*/, double ff[], void* /*userdata*/) {
57 throw CG_FATAL("cuba_integrand") << "Integrand not set for the Cuba algorithm!";
58 //TODO: handle the non-[0,1] ranges
59 ff[0] = Integrator::gIntegrand->eval(std::vector<double>(xx, xx + *ndim));
60 return 0;
61 }
62 } // namespace cuba
63} // namespace cepgen
#define CG_FATAL(mod)
Definition Exception.h:61
An integrand wrapper placeholder.
Definition Integrand.h:27
virtual double eval(const std::vector< double > &)=0
Compute the integrand for a given coordinates set.
static ParametersDescription description()
A description object for parameters collection.
A scalar value with its uncertainty.
Definition Value.h:26
Cuba integration algorithm.
Definition Integrator.h:29
Integrator(const ParametersList &)
virtual Value integrate()=0
static ParametersDescription description()
static Integrand * gIntegrand
Definition Integrator.h:34
int cuba_integrand(const int *ndim, const double xx[], const int *, double ff[], void *)
Common namespace for this Monte Carlo generator.