Open CASCADE Technology Reference Manual 8.0.0
Loading...
Searching...
No Matches
Data Structures | Functions
MathInteg Namespace Reference

Numerical integration algorithms. More...

Data Structures

struct  DoubleExpConfig
 Configuration for double exponential integration. More...
 
struct  KronrodConfig
 Configuration for Gauss-Kronrod integration. More...
 
struct  MultipleConfig
 Configuration for multi-dimensional Gauss integration. More...
 
struct  SetResult
 Result for vector function integration. More...
 

Functions

template<typename Function >
IntegResult Gauss (Function &theFunc, double theLower, double theUpper, int theNbPoints=15)
 Gauss-Legendre quadrature for definite integrals. Computes integral of f(x) from theLower to theUpper using n-point Gauss-Legendre rule.
 
template<typename Function >
IntegResult GaussAdaptive (Function &theFunc, double theLower, double theUpper, const IntegConfig &theConfig=IntegConfig())
 Adaptive Gauss-Legendre integration. Recursively subdivides interval until error estimate is below tolerance.
 
template<typename Function >
IntegResult GaussComposite (Function &theFunc, double theLower, double theUpper, int theNbIntervals, int theNbPoints=7)
 Composite Gauss-Legendre integration. Divides interval into subintervals and applies Gauss-Legendre to each. Simple alternative to adaptive integration.
 
template<typename Function >
IntegResult KronrodRule (Function &theFunc, double theLower, double theUpper, int theNbGauss=7)
 Apply Gauss-Kronrod rule to a single interval.
 
template<typename Function >
IntegResult Kronrod (Function &theFunc, double theLower, double theUpper, const KronrodConfig &theConfig=KronrodConfig())
 Gauss-Kronrod adaptive integration.
 
template<typename Function >
IntegResult KronrodAuto (Function &theFunc, double theLower, double theUpper, double theTolerance=1.0e-10, int theMaxOrder=30)
 Gauss-Kronrod integration with automatic order selection.
 
template<typename Function >
IntegResult KronrodSemiInfinite (Function &theFunc, double theLower, const KronrodConfig &theConfig=KronrodConfig())
 Gauss-Kronrod integration over semi-infinite interval [a, +infinity).
 
template<typename Function >
IntegResult KronrodInfinite (Function &theFunc, const KronrodConfig &theConfig=KronrodConfig())
 Gauss-Kronrod integration over infinite interval (-infinity, +infinity).
 
template<typename Function >
IntegResult TanhSinh (Function &theFunc, double theLower, double theUpper, const DoubleExpConfig &theConfig=DoubleExpConfig())
 Tanh-Sinh (Double Exponential) quadrature for finite interval [a,b].
 
template<typename Function >
IntegResult ExpSinh (Function &theFunc, double theLower, const DoubleExpConfig &theConfig=DoubleExpConfig())
 Exp-Sinh quadrature for semi-infinite interval [a, +infinity).
 
template<typename Function >
IntegResult SinhSinh (Function &theFunc, const DoubleExpConfig &theConfig=DoubleExpConfig())
 Sinh-Sinh quadrature for infinite interval (-infinity, +infinity).
 
template<typename Function >
IntegResult DoubleExponential (Function &theFunc, double theLower, double theUpper, const DoubleExpConfig &theConfig=DoubleExpConfig())
 Double exponential integration with automatic interval detection.
 
template<typename Function >
IntegResult TanhSinhSingular (Function &theFunc, double theLower, double theUpper, double theTolerance=1.0e-10)
 Tanh-Sinh quadrature optimized for endpoint singularities.
 
template<typename Function >
IntegResult TanhSinhWithSingularity (Function &theFunc, double theLower, double theUpper, double theSingularity, const DoubleExpConfig &theConfig=DoubleExpConfig())
 Integrates a function with a known singularity at an interior point.
 
template<typename Func >
IntegResult GaussMultiple (Func &theFunc, int theNVars, const math_Vector &theLower, const math_Vector &theUpper, const math_IntegerVector &theOrder, const MultipleConfig &theConfig=MultipleConfig())
 Gauss-Legendre integration of a multi-variable function.
 
template<typename Func >
IntegResult GaussMultipleUniform (Func &theFunc, int theNVars, const math_Vector &theLower, const math_Vector &theUpper, int theOrder)
 Gauss-Legendre integration with uniform order for all variables.
 
template<typename Func >
SetResult GaussSet (Func &theFunc, double theLower, double theUpper, int theOrder)
 Gauss-Legendre integration of a vector-valued function.
 
template<typename Func >
SetResult GaussSet (Func &theFunc, const math_Vector &theLower, const math_Vector &theUpper, int theOrder)
 Gauss-Legendre integration of vector function using math_Vector bounds.
 

Detailed Description

Numerical integration algorithms.

Function Documentation

◆ DoubleExponential()

template<typename Function >
IntegResult MathInteg::DoubleExponential ( Function & theFunc,
double theLower,
double theUpper,
const DoubleExpConfig & theConfig = DoubleExpConfig() )

Double exponential integration with automatic interval detection.

Selects the appropriate DE quadrature based on the interval:

  • Finite [a,b]: TanhSinh
  • Semi-infinite [a, +infinity): ExpSinh
  • Infinite (-infinity, +infinity): SinhSinh
Template Parameters
Functiontype with Value(double theX, double& theF) method
Parameters
theFuncfunction to integrate
theLowerlower bound (use -HUGE_VAL for -infinity)
theUpperupper bound (use HUGE_VAL for +infinity)
theConfigintegration configuration
Returns
integration result

◆ ExpSinh()

template<typename Function >
IntegResult MathInteg::ExpSinh ( Function & theFunc,
double theLower,
const DoubleExpConfig & theConfig = DoubleExpConfig() )

Exp-Sinh quadrature for semi-infinite interval [a, +infinity).

The transformation maps [0, +infinity) to (-infinity, +infinity): x = a + exp(pi/2 * sinh(t))

Template Parameters
Functiontype with Value(double theX, double& theF) method
Parameters
theFuncfunction to integrate
theLowerlower bound a
theConfigintegration configuration
Returns
integration result

◆ Gauss()

template<typename Function >
IntegResult MathInteg::Gauss ( Function & theFunc,
double theLower,
double theUpper,
int theNbPoints = 15 )

Gauss-Legendre quadrature for definite integrals. Computes integral of f(x) from theLower to theUpper using n-point Gauss-Legendre rule.

Algorithm:

  1. Transform interval [theLower, theUpper] to [-1, 1]
  2. Evaluate f at Gauss-Legendre points
  3. Sum weighted function values

Exact for polynomials of degree up to 2n-1.

Template Parameters
Functiontype with Value(double theX, double& theF) method
Parameters
theFuncfunction to integrate
theLowerlower integration bound
theUpperupper integration bound
theNbPointsnumber of quadrature points (>= 1)
Returns
result containing integral value

◆ GaussAdaptive()

template<typename Function >
IntegResult MathInteg::GaussAdaptive ( Function & theFunc,
double theLower,
double theUpper,
const IntegConfig & theConfig = IntegConfig() )

Adaptive Gauss-Legendre integration. Recursively subdivides interval until error estimate is below tolerance.

Algorithm:

  1. Compute integral using n and 2n points
  2. Estimate error as difference between the two
  3. If error > tolerance, subdivide and recurse
Template Parameters
Functiontype with Value(double theX, double& theF) method
Parameters
theFuncfunction to integrate
theLowerlower integration bound
theUpperupper integration bound
theConfigintegration configuration
Returns
result containing integral value and error estimate

◆ GaussComposite()

template<typename Function >
IntegResult MathInteg::GaussComposite ( Function & theFunc,
double theLower,
double theUpper,
int theNbIntervals,
int theNbPoints = 7 )

Composite Gauss-Legendre integration. Divides interval into subintervals and applies Gauss-Legendre to each. Simple alternative to adaptive integration.

Template Parameters
Functiontype with Value(double theX, double& theF) method
Parameters
theFuncfunction to integrate
theLowerlower integration bound
theUpperupper integration bound
theNbIntervalsnumber of subintervals
theNbPointsGauss points per interval (>= 1)
Returns
result containing integral value

◆ GaussMultiple()

template<typename Func >
IntegResult MathInteg::GaussMultiple ( Func & theFunc,
int theNVars,
const math_Vector & theLower,
const math_Vector & theUpper,
const math_IntegerVector & theOrder,
const MultipleConfig & theConfig = MultipleConfig() )

Gauss-Legendre integration of a multi-variable function.

Computes the N-dimensional integral using tensor product of 1D Gauss-Legendre quadrature: I = integral_{Lower}^{Upper} F(x1,...,xN) dx1...dxN

Uses recursive summation over all Gauss points in each dimension.

Template Parameters
Funcfunction type with bool Value(const math_Vector&, double&)
Parameters
theFuncN-dimensional function to integrate
theNVarsnumber of variables
theLowerlower bounds for each variable
theUpperupper bounds for each variable
theOrderintegration order for each variable (max 61)
theConfigoptional configuration
Returns
IntegResult containing the integral value

◆ GaussMultipleUniform()

template<typename Func >
IntegResult MathInteg::GaussMultipleUniform ( Func & theFunc,
int theNVars,
const math_Vector & theLower,
const math_Vector & theUpper,
int theOrder )

Gauss-Legendre integration with uniform order for all variables.

Template Parameters
Funcfunction type with bool Value(const math_Vector&, double&)
Parameters
theFuncN-dimensional function to integrate
theNVarsnumber of variables
theLowerlower bounds for each variable
theUpperupper bounds for each variable
theOrderintegration order for all variables
Returns
IntegResult containing the integral value

◆ GaussSet() [1/2]

template<typename Func >
SetResult MathInteg::GaussSet ( Func & theFunc,
const math_Vector & theLower,
const math_Vector & theUpper,
int theOrder )

Gauss-Legendre integration of vector function using math_Vector bounds.

Convenience overload that extracts scalar bounds from math_Vector.

Template Parameters
Funcfunction type
Parameters
theFuncvector-valued function
theLowerlower bound vector (only first element used)
theUpperupper bound vector (only first element used)
theOrderintegration order
Returns
SetResult containing vector of integrals

◆ GaussSet() [2/2]

template<typename Func >
SetResult MathInteg::GaussSet ( Func & theFunc,
double theLower,
double theUpper,
int theOrder )

Gauss-Legendre integration of a vector-valued function.

Integrates F: R -> R^N over interval [Lower, Upper]: Result[i] = integral_{Lower}^{Upper} F_i(x) dx

Uses Gauss-Legendre quadrature applied to each component.

Note: Only 1D input is supported (M=1). Multi-dimensional input is not implemented in the legacy API.

Template Parameters
Funcfunction type with:
  • int NbEquations() - number of output components
  • bool Value(const math_Vector& theX, math_Vector& theF)
Parameters
theFuncvector-valued function to integrate
theLowerlower bound
theUpperupper bound
theOrderintegration order (max 61)
Returns
SetResult containing vector of integrals

◆ Kronrod()

template<typename Function >
IntegResult MathInteg::Kronrod ( Function & theFunc,
double theLower,
double theUpper,
const KronrodConfig & theConfig = KronrodConfig() )

Gauss-Kronrod adaptive integration.

Uses adaptive bisection to achieve the requested tolerance. At each subdivision, the interval with the largest error estimate is bisected, and both halves are reintegrated.

Template Parameters
Functiontype with Value(double theX, double& theF) method
Parameters
theFuncfunction to integrate
theLowerlower integration bound
theUpperupper integration bound
theConfigintegration configuration
Returns
integration result with error estimate

◆ KronrodAuto()

template<typename Function >
IntegResult MathInteg::KronrodAuto ( Function & theFunc,
double theLower,
double theUpper,
double theTolerance = 1.0e-10,
int theMaxOrder = 30 )

Gauss-Kronrod integration with automatic order selection.

Starts with a low-order rule and increases the order until the tolerance is met or the maximum order is reached.

Template Parameters
Functiontype with Value(double theX, double& theF) method
Parameters
theFuncfunction to integrate
theLowerlower integration bound
theUpperupper integration bound
theTolerancerelative tolerance
theMaxOrdermaximum Gauss order to try
Returns
integration result

◆ KronrodInfinite()

template<typename Function >
IntegResult MathInteg::KronrodInfinite ( Function & theFunc,
const KronrodConfig & theConfig = KronrodConfig() )

Gauss-Kronrod integration over infinite interval (-infinity, +infinity).

Uses the substitution x = t / (1 - t^2) to map (-infinity, +infinity) to (-1, 1). The function must decay sufficiently fast at infinity.

Template Parameters
Functiontype with Value(double theX, double& theF) method
Parameters
theFuncfunction to integrate
theConfigintegration configuration
Returns
integration result

◆ KronrodRule()

template<typename Function >
IntegResult MathInteg::KronrodRule ( Function & theFunc,
double theLower,
double theUpper,
int theNbGauss = 7 )

Apply Gauss-Kronrod rule to a single interval.

The Gauss-Kronrod rule uses n Gauss points embedded in 2n+1 Kronrod points. The difference between the Gauss and Kronrod estimates provides an error estimate without additional function evaluations.

Template Parameters
Functiontype with Value(double theX, double& theF) method
Parameters
theFuncfunction to integrate
theLowerlower integration bound
theUpperupper integration bound
theNbGaussnumber of Gauss points (determines rule order)
Returns
integration result with error estimate

◆ KronrodSemiInfinite()

template<typename Function >
IntegResult MathInteg::KronrodSemiInfinite ( Function & theFunc,
double theLower,
const KronrodConfig & theConfig = KronrodConfig() )

Gauss-Kronrod integration over semi-infinite interval [a, +infinity).

Uses the substitution x = a + t / (1 - t) to map [a, +infinity) to [0, 1).

Template Parameters
Functiontype with Value(double theX, double& theF) method
Parameters
theFuncfunction to integrate
theLowerlower bound a
theConfigintegration configuration
Returns
integration result

◆ SinhSinh()

template<typename Function >
IntegResult MathInteg::SinhSinh ( Function & theFunc,
const DoubleExpConfig & theConfig = DoubleExpConfig() )

Sinh-Sinh quadrature for infinite interval (-infinity, +infinity).

The transformation maps (-infinity, +infinity) to (-infinity, +infinity): x = sinh(pi/2 * sinh(t))

Template Parameters
Functiontype with Value(double theX, double& theF) method
Parameters
theFuncfunction to integrate
theConfigintegration configuration
Returns
integration result

◆ TanhSinh()

template<typename Function >
IntegResult MathInteg::TanhSinh ( Function & theFunc,
double theLower,
double theUpper,
const DoubleExpConfig & theConfig = DoubleExpConfig() )

Tanh-Sinh (Double Exponential) quadrature for finite interval [a,b].

The double exponential transformation maps [a,b] to (-inf, +inf): x = (b+a)/2 + (b-a)/2 * tanh(pi/2 * sinh(t))

The integrand is multiplied by the Jacobian: dx/dt = (b-a)/2 * (pi/2) * cosh(t) / cosh^2(pi/2 * sinh(t))

Properties:

  • Excellent for functions with endpoint singularities
  • Exponential convergence for analytic functions
  • Self-adaptive through level refinement
  • Handles algebraic and logarithmic singularities
Template Parameters
Functiontype with Value(double theX, double& theF) method
Parameters
theFuncfunction to integrate
theLowerlower integration bound
theUpperupper integration bound
theConfigintegration configuration
Returns
integration result with error estimate

◆ TanhSinhSingular()

template<typename Function >
IntegResult MathInteg::TanhSinhSingular ( Function & theFunc,
double theLower,
double theUpper,
double theTolerance = 1.0e-10 )

Tanh-Sinh quadrature optimized for endpoint singularities.

This variant uses precomputed nodes and weights for better performance when integrating functions with algebraic singularities at the endpoints.

Template Parameters
Functiontype with Value(double theX, double& theF) method
Parameters
theFuncfunction to integrate
theLowerlower bound
theUpperupper bound
theTolerancerelative tolerance
Returns
integration result

◆ TanhSinhWithSingularity()

template<typename Function >
IntegResult MathInteg::TanhSinhWithSingularity ( Function & theFunc,
double theLower,
double theUpper,
double theSingularity,
const DoubleExpConfig & theConfig = DoubleExpConfig() )

Integrates a function with a known singularity at an interior point.

Splits the interval at the singularity and integrates each part.

Template Parameters
Functiontype with Value(double theX, double& theF) method
Parameters
theFuncfunction to integrate
theLowerlower bound
theUpperupper bound
theSingularitylocation of interior singularity
theConfigintegration configuration
Returns
combined integration result