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

#include <GaussianIntegrator.h>

Inheritance diagram for Herwig::GaussianIntegrator:

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
 
std::vector< std::vector< double > > _weights
 The weights for the gaussian quadrature.
 
std::vector< std::vector< double > > _abscissae
 The abscissae.
 
double _abserr
 The parameters controlling the error.
 
double _relerr
 
double _binwidth
 The minimum width of a bin as a fraction of the integration region.
 
int _maxeval
 Maximum number of function evaluations.
 
 GaussianIntegrator ()
 Default Constructor.
 
 GaussianIntegrator (double abserr, double relerr, double binwidth, int maxeval)
 Specify all the parameters.
 
template<class T >
ValT< T > value (const T &, const typename T::ArgType lower, const typename T::ArgType upper) const
 The value of the integral.
 
GaussianIntegratoroperator= (const GaussianIntegrator &)=delete
 The assignment operator is private and must never be called.
 
void Init ()
 Initialise the weights and abscissae.
 

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

Author
Peter Richardson
This class is designed to perform the integral of a function using Gaussian quadrature.The method is adaptive based on using 6th,12th, 24th,48th, or 96th order Gaussian quadrature combined with subdivision of the integral if this is insufficient.

The class is templated on a simple class which should provide a T::operator () (double) const which provides the integrand for the function.

Definition at line 33 of file GaussianIntegrator.h.

Member Typedef Documentation

◆ ValT

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

helper type for the integration result

Definition at line 66 of file GaussianIntegrator.h.

Constructor & Destructor Documentation

◆ GaussianIntegrator() [1/2]

Herwig::GaussianIntegrator::GaussianIntegrator ( )
inline

Default Constructor.

Definition at line 42 of file GaussianIntegrator.h.

References Init().

◆ GaussianIntegrator() [2/2]

Herwig::GaussianIntegrator::GaussianIntegrator ( double  abserr,
double  relerr,
double  binwidth,
int  maxeval 
)
inline

Specify all the parameters.

Parameters
abserrAbsolute error.
relerrRelative error.
binwidthWidth of the bin as a fraction of the integration region.
maxevalMaximum number of function evaluations

Definition at line 55 of file GaussianIntegrator.h.

References Init().

Member Function Documentation

◆ operator=()

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

The assignment operator is private and must never be called.

In fact, it should not even be implemented.

◆ value()

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

The value of the integral.

Parameters
lowerThe lower limit of integration.
upperThe upper limit of integration.

Member Data Documentation

◆ _abscissae

std::vector< std::vector <double> > Herwig::GaussianIntegrator::_abscissae
private

The abscissae.

Definition at line 102 of file GaussianIntegrator.h.

◆ _abserr

double Herwig::GaussianIntegrator::_abserr
private

The parameters controlling the error.

Definition at line 107 of file GaussianIntegrator.h.

◆ _binwidth

double Herwig::GaussianIntegrator::_binwidth
private

The minimum width of a bin as a fraction of the integration region.

Definition at line 112 of file GaussianIntegrator.h.

◆ _maxeval

int Herwig::GaussianIntegrator::_maxeval
private

Maximum number of function evaluations.

Definition at line 117 of file GaussianIntegrator.h.

◆ _relerr

double Herwig::GaussianIntegrator::_relerr
private

Definition at line 107 of file GaussianIntegrator.h.

◆ _weights

std::vector< std::vector<double> > Herwig::GaussianIntegrator::_weights
private

The weights for the gaussian quadrature.

Definition at line 97 of file GaussianIntegrator.h.


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