PARTONS/NumA++  
Numerical Analysis C++ routines
Public Member Functions | Protected Member Functions | Private Attributes | List of all members
NumA::DExpIntegrator1D Class Reference

This is an implementation of the double exponential rule. More...

Inheritance diagram for NumA::DExpIntegrator1D:
NumA::Integrator1D

Public Member Functions

 DExpIntegrator1D ()
 Default constructor. More...
 
virtual ~DExpIntegrator1D ()
 Default destructor. More...
 
DExpIntegrator1Dclone () const
 
virtual double integrate (FunctionType1D *pFunction, double a, double b, std::vector< double > &parameters)
 Integration routine. More...
 
- Public Member Functions inherited from NumA::Integrator1D
 Integrator1D ()
 Default constructor. More...
 
virtual ~Integrator1D ()
 Default destructor. More...
 
virtual void configure (const ElemUtils::Parameters &parameters)
 Provides a generic method to configure all types of integrations by passing a Parameters object. More...
 
const TolerancesgetTolerances () const
 
void setTolerances (const Tolerances &tolerances)
 
const ErrorsgetErrors () const
 
void setErrors (const Errors &errors)
 

Protected Member Functions

 DExpIntegrator1D (const DExpIntegrator1D &other)
 Copy constructor. More...
 
- Protected Member Functions inherited from NumA::Integrator1D
 Integrator1D (const Integrator1D &other)
 Copy constructor. More...
 

Private Attributes

double m_errorEstimate
 Absolute error estimate. More...
 
int m_numFunctionEvaluations
 Number of evaluations. More...
 
std::vector< double > m_nodes
 Nodes of the quadrature. More...
 
std::vector< double > m_weights
 Weights of the quadrature. More...
 
std::vector< int > m_offsets
 

Additional Inherited Members

- Static Public Member Functions inherited from NumA::Integrator1D
template<typename PointerToObj , typename PointerToMemFn >
static Functor1D< PointerToObj, PointerToMemFn > * newIntegrationFunctor (PointerToObj *object, PointerToMemFn function)
 Use FunctorUtils::newFunctor1D instead. More...
 
static Integrator1DnewIntegrator (const IntegratorType1D::Type &oneDimIntegratorType)
 Instantiates an Integrator object. More...
 
- Static Public Attributes inherited from NumA::Integrator1D
static const std::string PARAM_NAME_ABSOLUTE_TOLERANCE
 Parameter used in the configure() to set the absolute tolerance. More...
 
static const std::string PARAM_NAME_RELATIVE_TOLERANCE
 Parameter used in the configure() to set the relative tolerance. More...
 
- Protected Attributes inherited from NumA::Integrator1D
Tolerances m_tolerances
 Absolute and relative tolerances. More...
 
Errors m_errors
 Absolute and relative errors estimations. More...
 

Detailed Description

This is an implementation of the double exponential rule.

The functions being integrated are required to be smooth (no discontinuities in the function or any of its derivatives), with no infinities in the middle of the integration interval (but it does allow functions that become infinite at the ends of the integration interval).

The method is based on the observation that the trapezoid rule converges extremely fast for functions that go to zero like exp(-exp(t)). In practice, it implements a change of variables to implement this property in the function to be integrated.

The change of variable is x = tanh( pi sinh(t) /2) which transforms an integral over [-1, 1] into an integral with integrand suited to the double exponential rule. The transformed integral is infinite, but one can truncate the domain of integration to [-3, 3].

The limit '3' is chosen for two reasons:

  1. The transformed x values are equal to 1 for 12 or more significant figures;
  2. The smallest weights are 12 orders of magnitude smaller than the largest weights (setting the cutoff larger than 3 would not have a significant impact on the integral value unless there is a strong singularity at one of the end points).

The change of variables x(t) is an odd function, whereas its derivative w(t) is even; thus in the cpp file we will store only positive values of nodes and weights.

The integration first applies the trapezoid rule to [-3, 3] in steps of size 1, subsequently cutting the step size in half each time, and comparing the results. The routine stops when subsequent iterations are close enough together or the maximum integration points have been used. Notice that by cutting h in half, the previous integral can be reused; we only need evaluate the integrand at the newly added points.

As we assume that values at the ends of the interval hardly matter due to the rapid decay of the integrand we don't treat the end points differently.

See Integrator1D documentation for an example.

Constructor & Destructor Documentation

◆ DExpIntegrator1D() [1/2]

NumA::DExpIntegrator1D::DExpIntegrator1D ( )

Default constructor.

◆ ~DExpIntegrator1D()

NumA::DExpIntegrator1D::~DExpIntegrator1D ( )
virtual

Default destructor.

◆ DExpIntegrator1D() [2/2]

NumA::DExpIntegrator1D::DExpIntegrator1D ( const DExpIntegrator1D other)
protected

Copy constructor.

Called by clone().

Parameters
otherDExpIntegrator1D object to copy.

Member Function Documentation

◆ clone()

DExpIntegrator1D * NumA::DExpIntegrator1D::clone ( ) const
virtual

Implements NumA::Integrator1D.

◆ integrate()

double NumA::DExpIntegrator1D::integrate ( FunctionType1D pFunction,
double  a,
double  b,
std::vector< double > &  parameters 
)
virtual

Integration routine.

Parameters
pFunctionFunctor representing the one-dimensional function to integrate.
aLower bound.
bUpper bound.
parametersParameters that can be passed to the function.
Returns
Integral.

Implements NumA::Integrator1D.

Member Data Documentation

◆ m_errorEstimate

double NumA::DExpIntegrator1D::m_errorEstimate
private

Absolute error estimate.

◆ m_nodes

std::vector<double> NumA::DExpIntegrator1D::m_nodes
private

Nodes of the quadrature.

◆ m_numFunctionEvaluations

int NumA::DExpIntegrator1D::m_numFunctionEvaluations
private

Number of evaluations.

◆ m_offsets

std::vector<int> NumA::DExpIntegrator1D::m_offsets
private

◆ m_weights

std::vector<double> NumA::DExpIntegrator1D::m_weights
private

Weights of the quadrature.


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