cepgen is hosted by Hepforge, IPPP Durham

Monte Carlo numerical integrators#

The modularity of CepGen also allows for multiple integration algorithms to be steered and profit for some increased numerical stability of a wise choice of parameters. Several interfaces to external algorithms are provided in the core and CepGenAddOns libraries.

In the Python cards parsing these can be steered through the integrator keyword. All modules are derived from a common cepgen::Integrator object, described below:

A full list of the algorithms and their parameters can be found here.

class Integrator : public NamedModule<Integrator>#

Subclassed by BasesIntegrator, FoamIntegrator, GSLIntegrator, NaiveBoostIntegrator, Integrator, Integrator, Integrator

Detailed description

class Integrator : public NamedModule<Integrator>

Monte-Carlo integration algorithm.

Subclassed by BasesIntegrator, FoamIntegrator, GSLIntegrator, NaiveBoostIntegrator, Integrator, Integrator, Integrator

Public Functions

explicit Integrator(const ParametersList &params)

Integrator algorithm constructor.

void checkLimits(const Integrand&)

Ensure the integration bounds are properly set.

inline virtual void setLimits(const std::vector<Limits> &limits)

Specify the variables limits on integration.

virtual double eval(Integrand&, const std::vector<double>&) const

Compute the function value at the given phase space point.

virtual double uniform(const Limits& = {0., 1.}) const

Generate a uniformly distributed (between 0 and 1) random number.

virtual Value integrate(Integrand &result) = 0

Perform the multidimensional Monte Carlo integration.

Parameters:

result[out] integral computed over the full phase space

inline const std::string &name() const

Module unique indexing name.

inline bool operator==(const SteeredObject &oth) const

Equality operator.

inline bool operator!=(const SteeredObject &oth) const

Inequality operator.

inline virtual const ParametersList &parameters() const override

Module user-defined parameters.

inline virtual void setParameters(const ParametersList &params) override

Set module parameters.

inline void setDescribedParameters(const ParametersList &params_orig)

Set (documented) module parameters.

Public Static Functions

static Value integrate(const std::function<double(const std::vector<double>&)>&, const ParametersList&, size_t)

Perform an integration with a given functional and a given set of parameters.

static Value integrate(const std::function<double(const std::vector<double>&)>&, const ParametersList&, const std::vector<Limits>&)

Perform an integration with a given functional and a given set of parameters.


GSL Monte Carlo integrator algorithms#

Historical integrator algorithms, relying on the GSL implementation of gsl_monte_xxx integration routines (see the GSL manual for a description of all relevent parameters).

class GSLIntegrator : public Integrator#

Subclassed by MISERIntegrator, PlainIntegrator, VegasIntegrator

Vegas integration algorithm#

class VegasIntegrator : public GSLIntegrator#

Vegas integration algorithm developed by P. Lepage, as documented in [18].

Private Functions

inline virtual double eval(Integrand &integrand, const std::vector<double> &x) const override#

Compute the function value at the given phase space point.

Private Members

const bool treat_#

Is the integrand to be smoothed for events generation?

std::unique_ptr<gsl_monte_vegas_state, gsl_monte_vegas_deleter> vegas_state_#

A Vegas integrator state for integration (optional) and/or “treated” event generation.

struct gsl_monte_vegas_deleter#

A trivial deleter for the Vegas integrator.

Miser integration algorithm#

class MISERIntegrator : public GSLIntegrator#

MISER integration algorithm developed by W.H. Press and G.R. Farrar, as documented in [23].

Plain integration algorithm#

class PlainIntegrator : public GSLIntegrator#

Plain integration algorithm randomly sampling points in the phase space.

ROOT integration algorithms#

Interfacing to the two ROOT numerical integration algorithms:

According to the user’s request, either of the two objects is populated and configured. The following parameters are to be steered by the end user:

  • type: integrator algorithm type: : - gauss, legendre, adaptive, adaptiveSingular, nonAdaptive, for one-dimensional integration,

    • adaptive, plain, miser, vegas for multidimensional integration. The last three are one-to-one equivalent to the GSL Monte Carlo algorithms described above (except for the interface and parameters definition).

  • absToL: desired absolute error ;

  • relToL: desired relative error ;

  • size: maximum number of sub-intervals.

Added in version 0.9.10.

Warning

doxygenclass: Cannot find class “cepgen::ROOTIntegrator” in doxygen xml output for project “CepGen” from directory: /Documentation/build/xml

Python integration algorithms#

Among the newest features is the interfacing between CepGen and any Python numerical integrator, provided the definition of a wrapper function (in Python) called by an interfacing object:

Added in version 1.2.0.

Warning

doxygenclass: Cannot find class “cepgen::PythonIntegrator” in doxygen xml output for project “CepGen” from directory: /Documentation/build/xml

The Python interfacing/wrapper function for CepGen numerical integration has the form:

def integrate(f,  # CepGen function to integrate (a [](vector<double>) -> double function)
              num_dim: int,  # number of dimensions to integrate
              num_iter: int,  # number of iterations for integration (if used)
              num_warmup: int,  # number of calls to perform for an integrator warmup (if any)
              num_calls: int,  # number of function calls for each iteration
              limits: list[tuple[float]]=[]  # list of (min, max) limits for each dimension
              ):
    # definition of integration procedure
    # ...
    return (average, standard_deviation)