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
19
#include "
CepGen/Core/Exception.h
"
20
#include "
CepGen/Integration/Integrand.h
"
21
#include "
CepGen/Integration/Integrator.h
"
22
#include "
CepGen/Modules/IntegratorFactory.h
"
23
#include "
CepGen/Utils/String.h
"
24
#include "
CepGenAddOns/BasesWrapper/BasesCommonBlocks.h
"
25
26
namespace
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
38
static
ParametersDescription
description
() {
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
{
49
Integrator::setLimits
(lims);
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
82
REGISTER_INTEGRATOR
(
"bases"
, BasesIntegrator);
BasesCommonBlocks.h
bases_
void bases_(double(*fxn)(double[]), double &s, double &sigma, double &ctime, int &it1, int &it2)
bparm1_
struct @2 bparm1_
bscntl_
struct @4 bscntl_
bsinit_
void bsinit_()
Exception.h
CG_FATAL
#define CG_FATAL(mod)
Definition
Exception.h:61
Integrand.h
Integrator.h
IntegratorFactory.h
REGISTER_INTEGRATOR
#define REGISTER_INTEGRATOR(name, obj)
Add a generic process definition to the list of handled processes.
Definition
IntegratorFactory.h:25
CG_DEBUG
#define CG_DEBUG(mod)
Definition
Message.h:220
String.h
cepgen::BasesIntegrator
Bases integration algorithm.
Definition
BasesIntegrator.cpp:28
cepgen::BasesIntegrator::integrate
Value integrate(Integrand &integr) override
Perform the multidimensional Monte Carlo integration.
Definition
BasesIntegrator.cpp:56
cepgen::BasesIntegrator::description
static ParametersDescription description()
Definition
BasesIntegrator.cpp:38
cepgen::BasesIntegrator::setLimits
void setLimits(const std::vector< Limits > &lims) override
Specify the variables limits on integration.
Definition
BasesIntegrator.cpp:48
cepgen::BasesIntegrator::BasesIntegrator
BasesIntegrator(const ParametersList ¶ms)
Definition
BasesIntegrator.cpp:30
cepgen::Integrand
An integrand wrapper placeholder.
Definition
Integrand.h:27
cepgen::Integrand::eval
virtual double eval(const std::vector< double > &)=0
Compute the integrand for a given coordinates set.
cepgen::Integrand::size
virtual size_t size() const =0
Phase space dimension.
cepgen::Integrator
Monte-Carlo integration algorithm.
Definition
Integrator.h:28
cepgen::Integrator::checkLimits
void checkLimits(const Integrand &)
Ensure the integration bounds are properly set.
Definition
Integrator.cpp:32
cepgen::Integrator::setLimits
virtual void setLimits(const std::vector< Limits > &limits)
Specify the variables limits on integration.
Definition
Integrator.h:38
cepgen::Integrator::limits_
std::vector< Limits > limits_
List of per-variable integration limits.
Definition
Integrator.h:58
cepgen::Integrator::description
static ParametersDescription description()
Definition
Integrator.cpp:71
cepgen::ParametersDescription
A description object for parameters collection.
Definition
ParametersDescription.h:26
cepgen::ParametersList
Definition
ParametersList.h:52
cepgen::Value
A scalar value with its uncertainty.
Definition
Value.h:26
cepgen::utils::s
std::string s(const std::string &word, float num, bool show_number)
Add a trailing "s" when needed.
Definition
String.cpp:228
cepgen
Common namespace for this Monte Carlo generator.
Definition
CommandLineHandler.cpp:36
CepGenAddOns
BasesWrapper
BasesIntegrator.cpp
Generated on Mon Jul 29 2024 for CepGen by
1.9.7