herwig is hosted by Hepforge, IPPP Durham
Herwig 7.3.0
Herwig::GSLIntegrator Class Reference

This class is designed to integrate a given function between 2 limits using the gsl QAGS integration subroutine. More...

#include <GSLIntegrator.h>

Inheritance diagram for Herwig::GSLIntegrator:

Standard constructors and destructors.

template<class T >
using ValT = decltype(std::declval< typename T::ValType >() *std::declval< typename T::ArgType >())
 helper type for the integration result
 
double _abserr
 The parameters controlling the absolute error.
 
double _relerr
 The parameters controlling the relative error.
 
int _nbins
 The maximum number of intervals to use.
 
 GSLIntegrator ()
 Default Constructor uses values in GSL manual as parameters.
 
 GSLIntegrator (double abserr, double relerr, int nbins)
 Specify all the parameters.
 
template<class T >
ValT< T > value (const T &function, const typename T::ArgType lower, const typename T::ArgType upper) const
 The value of the integral.
 
template<class T >
ValT< T > value (const T &function, const typename T::ArgType lower, const typename T::ArgType upper, ValT< T > &error) const
 The value of the integral.
 
GSLIntegratoroperator= (const GSLIntegrator &)=delete
 The assignment operator is private and must never be called.
 

Additional Inherited Members

- Public Types inherited from ThePEG::Pointer::ReferenceCounted
typedef unsigned int CounterType
 
- Public Member Functions inherited from ThePEG::Pointer::ReferenceCounted
CounterType referenceCount () const
 
- Public Attributes inherited from ThePEG::Pointer::ReferenceCounted
const unsigned long uniqueId
 
- Protected Member Functions inherited from ThePEG::Pointer::ReferenceCounted
 ReferenceCounted (const ReferenceCounted &)
 
ReferenceCountedoperator= (const ReferenceCounted &)
 

Detailed Description

This class is designed to integrate a given function between 2 limits using the gsl QAGS integration subroutine.

The function is supplied using a templated class that must define operator(argument). The units of the argument ArgType and return type ValType must be supplied in the integrand class using a typedef i.e.
struct integrand {
...
Energy operator(double arg) const;
typedef double ArgType
typedef Energy ValType
...
}

Definition at line 40 of file GSLIntegrator.h.

Member Typedef Documentation

◆ ValT

template<class T >
using Herwig::GSLIntegrator::ValT = decltype(std::declval<typename T::ValType>() * std::declval<typename T::ArgType>())

helper type for the integration result

Definition at line 63 of file GSLIntegrator.h.

Constructor & Destructor Documentation

◆ GSLIntegrator() [1/2]

Herwig::GSLIntegrator::GSLIntegrator ( )
inline

Default Constructor uses values in GSL manual as parameters.

Definition at line 49 of file GSLIntegrator.h.

◆ GSLIntegrator() [2/2]

Herwig::GSLIntegrator::GSLIntegrator ( double  abserr,
double  relerr,
int  nbins 
)
inline

Specify all the parameters.

Parameters
abserrAbsolute error.
relerrRelative error.
nbinsNumber of bins

Definition at line 57 of file GSLIntegrator.h.

Member Function Documentation

◆ operator=()

GSLIntegrator & Herwig::GSLIntegrator::operator= ( const GSLIntegrator )
privatedelete

The assignment operator is private and must never be called.

In fact, it should not even be implemented.

◆ value() [1/2]

template<class T >
ValT< T > Herwig::GSLIntegrator::value ( const T &  function,
const typename T::ArgType  lower,
const typename T::ArgType  upper 
) const
inline

The value of the integral.

Parameters
functionThe integrand class that defines operator()
lowerThe lower limit of integration.
upperThe upper limit of integration.

Referenced by Herwig::ThreeBodyAllOnCalculator< T >::Outer::operator()().

◆ value() [2/2]

template<class T >
ValT< T > Herwig::GSLIntegrator::value ( const T &  function,
const typename T::ArgType  lower,
const typename T::ArgType  upper,
ValT< T > &  error 
) const
inline

The value of the integral.

Parameters
functionThe integrand class that defines operator()
lowerThe lower limit of integration.
upperThe upper limit of integration.
errorReturns the estimated error of the integral

Member Data Documentation

◆ _abserr

double Herwig::GSLIntegrator::_abserr
private

The parameters controlling the absolute error.

Definition at line 105 of file GSLIntegrator.h.

◆ _nbins

int Herwig::GSLIntegrator::_nbins
private

The maximum number of intervals to use.

Definition at line 115 of file GSLIntegrator.h.

◆ _relerr

double Herwig::GSLIntegrator::_relerr
private

The parameters controlling the relative error.

Definition at line 110 of file GSLIntegrator.h.


The documentation for this class was generated from the following file: