![]() |
Open CASCADE Technology Reference Manual 8.0.0
|
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. | |
Numerical integration algorithms.
| 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:
| Function | type with Value(double theX, double& theF) method |
| theFunc | function to integrate |
| theLower | lower bound (use -HUGE_VAL for -infinity) |
| theUpper | upper bound (use HUGE_VAL for +infinity) |
| theConfig | integration configuration |
| 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))
| Function | type with Value(double theX, double& theF) method |
| theFunc | function to integrate |
| theLower | lower bound a |
| theConfig | integration configuration |
| 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:
Exact for polynomials of degree up to 2n-1.
| Function | type with Value(double theX, double& theF) method |
| theFunc | function to integrate |
| theLower | lower integration bound |
| theUpper | upper integration bound |
| theNbPoints | number of quadrature points (>= 1) |
| 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:
| Function | type with Value(double theX, double& theF) method |
| theFunc | function to integrate |
| theLower | lower integration bound |
| theUpper | upper integration bound |
| theConfig | integration configuration |
| 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.
| Function | type with Value(double theX, double& theF) method |
| theFunc | function to integrate |
| theLower | lower integration bound |
| theUpper | upper integration bound |
| theNbIntervals | number of subintervals |
| theNbPoints | Gauss points per interval (>= 1) |
| 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.
| Func | function type with bool Value(const math_Vector&, double&) |
| theFunc | N-dimensional function to integrate |
| theNVars | number of variables |
| theLower | lower bounds for each variable |
| theUpper | upper bounds for each variable |
| theOrder | integration order for each variable (max 61) |
| theConfig | optional configuration |
| 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.
| Func | function type with bool Value(const math_Vector&, double&) |
| theFunc | N-dimensional function to integrate |
| theNVars | number of variables |
| theLower | lower bounds for each variable |
| theUpper | upper bounds for each variable |
| theOrder | integration order for all variables |
| 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.
| Func | function type |
| theFunc | vector-valued function |
| theLower | lower bound vector (only first element used) |
| theUpper | upper bound vector (only first element used) |
| theOrder | integration order |
| 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.
| Func | function type with:
|
| theFunc | vector-valued function to integrate |
| theLower | lower bound |
| theUpper | upper bound |
| theOrder | integration order (max 61) |
| 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.
| Function | type with Value(double theX, double& theF) method |
| theFunc | function to integrate |
| theLower | lower integration bound |
| theUpper | upper integration bound |
| theConfig | integration configuration |
| 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.
| Function | type with Value(double theX, double& theF) method |
| theFunc | function to integrate |
| theLower | lower integration bound |
| theUpper | upper integration bound |
| theTolerance | relative tolerance |
| theMaxOrder | maximum Gauss order to try |
| 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.
| Function | type with Value(double theX, double& theF) method |
| theFunc | function to integrate |
| theConfig | integration configuration |
| 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.
| Function | type with Value(double theX, double& theF) method |
| theFunc | function to integrate |
| theLower | lower integration bound |
| theUpper | upper integration bound |
| theNbGauss | number of Gauss points (determines rule order) |
| 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).
| Function | type with Value(double theX, double& theF) method |
| theFunc | function to integrate |
| theLower | lower bound a |
| theConfig | integration configuration |
| 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))
| Function | type with Value(double theX, double& theF) method |
| theFunc | function to integrate |
| theConfig | integration configuration |
| 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:
| Function | type with Value(double theX, double& theF) method |
| theFunc | function to integrate |
| theLower | lower integration bound |
| theUpper | upper integration bound |
| theConfig | integration configuration |
| 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.
| Function | type with Value(double theX, double& theF) method |
| theFunc | function to integrate |
| theLower | lower bound |
| theUpper | upper bound |
| theTolerance | relative tolerance |
| 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.
| Function | type with Value(double theX, double& theF) method |
| theFunc | function to integrate |
| theLower | lower bound |
| theUpper | upper bound |
| theSingularity | location of interior singularity |
| theConfig | integration configuration |