cepgen is hosted by Hepforge, IPPP Durham
CepGen 1.2.5
Central exclusive processes event generator
Loading...
Searching...
No Matches
BasesIntegrator.cpp
Go to the documentation of this file.
1/*
2 * CepGen: a central exclusive processes event generator
3 * Copyright (C) 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
23#include "CepGen/Utils/String.h"
25
26namespace cepgen {
28 class BasesIntegrator : public Integrator {
29 public:
30 explicit BasesIntegrator(const ParametersList& params) : Integrator(params) {
31 bsinit_();
32 bparm1_.ncall = steer<int>("numFunctionCalls");
33 std::fill(bparm1_.ig.begin(), bparm1_.ig.end(), false);
34 bscntl_.intv = steer<int>("intv");
35 bscntl_.ipnt = steer<int>("verbose");
36 }
37
39 auto desc = Integrator::description();
40 desc.setDescription("Bases integration algorithm");
41 desc.add<int>("numFunctionCalls", 50'000);
42 desc.add<int>("intv", 1);
43 desc.add<int>("verbose", 0);
44 desc.add<std::vector<int> >("wildVars", {}).setDescription("list of 'wild' variables");
45 return desc;
46 }
47
48 void setLimits(const std::vector<Limits>& lims) override {
50 for (size_t i = 0; i < limits_.size(); ++i) {
51 bparm1_.xl[i] = limits_.at(i).min();
52 bparm1_.xu[i] = limits_.at(i).max();
53 }
54 }
55
56 Value integrate(Integrand& integr) override {
57 checkLimits(integr); // check the integration bounds
58 bparm1_.ndim = integr.size();
59 for (const auto& wc : steer<std::vector<int> >("wildVars")) {
60 if (wc < 0 || wc >= bparm1_.ndim)
61 throw CG_FATAL("BasesIntegrator:integrate") << "Invalid 'wild' variable coordinate set: " << wc << ".";
62 bparm1_.ig[wc] = true;
63 ++bparm1_.nwild;
64 }
65 double res, unc, ctime;
66 int it1, it2;
67 if (gIntegrand = &integr; !gIntegrand)
68 throw CG_FATAL("BasesIntegrator") << "Integrand was not specified before integration.";
69 bases_(integrand, res, unc, ctime, it1, it2);
70 CG_DEBUG("BasesIntegrator:integrate")
71 << "Integration performed in " << utils::s("second", ctime, true) << ". " << utils::s("iteration", it1, true)
72 << " for the grid definition, " << utils::s("iteration", it2, true) << " for the integration.";
73 return Value{res, unc};
74 }
75
76 private:
77 static Integrand* gIntegrand;
78 static double integrand(double* in) { return gIntegrand->eval(std::vector<double>(in, in + gIntegrand->size())); }
79 };
80 Integrand* BasesIntegrator::gIntegrand = nullptr;
81} // namespace cepgen
82REGISTER_INTEGRATOR("bases", BasesIntegrator);
void bases_(double(*fxn)(double[]), double &s, double &sigma, double &ctime, int &it1, int &it2)
struct @2 bparm1_
struct @4 bscntl_
void bsinit_()
#define CG_FATAL(mod)
Definition Exception.h:61
#define REGISTER_INTEGRATOR(name, obj)
Add a generic process definition to the list of handled processes.
#define CG_DEBUG(mod)
Definition Message.h:220
Bases integration algorithm.
Value integrate(Integrand &integr) override
Perform the multidimensional Monte Carlo integration.
static ParametersDescription description()
void setLimits(const std::vector< Limits > &lims) override
Specify the variables limits on integration.
BasesIntegrator(const ParametersList &params)
An integrand wrapper placeholder.
Definition Integrand.h:27
virtual double eval(const std::vector< double > &)=0
Compute the integrand for a given coordinates set.
virtual size_t size() const =0
Phase space dimension.
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()
A description object for parameters collection.
A scalar value with its uncertainty.
Definition Value.h:26
std::string s(const std::string &word, float num, bool show_number)
Add a trailing "s" when needed.
Definition String.cpp:228
Common namespace for this Monte Carlo generator.