Adaptive Clenshaw Curtis Quadrature#
-
template<typename MemorySpace = Kokkos::HostSpace>
class AdaptiveClenshawCurtis : public mpart::RecursiveQuadratureBase<Kokkos::HostSpace># Adaptive quadrature based on applying a Clenshaw-Curtis recursively on subintervals.
The points in a Clenshaw-Curtis quadrature rule with \(2^{L}+1\) points are a subset of the points in a rule with \(2^{L+1}+1\) points. This “nesting” allows us to approximate the integral at two different levels with minimal additional integrand evaluations. Comparing the integral at these values gives an indication of the error in the integral approximation. If the error is unacceptably large, we can subdivide the integration domain and apply Clenshaw-Curtis rules on each subinterval. This class implements a recursive version of this process. Two nested Clenshaw-Curtis rules are used to estimate the integral, and its error. If the error is too large, then the integration domain is split into two equal halves, where the nested quadrature rules can again be applied. This recursive subdivision repeats until an acceptable error level is reached or until a maximum number of subdivisions has occured.
Public Functions
-
inline AdaptiveClenshawCurtis(unsigned int level, unsigned int maxSub, unsigned int maxDim, double absTol, double relTol, QuadError::Type errorMetric, unsigned int minSub = 0)#
Construct a new adaptive quadrature class with specified stopping criteria.
- Parameters:
maxSub – The maximum number of subintervals allowed.
maxDim – The maximum dimension of the integrand.
absTol – [in] An absolute error tolerance used to stop the adaptive integration.
relTol – [in] A relative error tolerance used to stop te adaptive integration.
errorMetric – [in] A flag specifying the type of error metric to use.
level – [in] The nesting level \(L\). A coarse rule with \(2^{L}+1\) points will be used and a fine rule with \(w^{L+1}+1\) points will be used.
-
inline AdaptiveClenshawCurtis(unsigned int level, unsigned int maxSub, unsigned int maxDim, double *workspace, double absTol, double relTol, QuadError::Type errorMetric, unsigned int minSub = 0)#
-
inline AdaptiveClenshawCurtis(Kokkos::View<double*, MemorySpace> coarsePts, Kokkos::View<double*, MemorySpace> coarseWts, Kokkos::View<double*, MemorySpace> finePts, Kokkos::View<double*, MemorySpace> fineWts, unsigned int maxSub, unsigned int maxDim, double *workspace, double absTol, double relTol, QuadError::Type errorMetric, unsigned int minSub = 0)#
-
inline void SetDim(unsigned int fdim)#
-
template<class FunctionType>
inline void Integrate(FunctionType const &f, double lb, double ub, double *res) const# Approximates the integral \(\int_{x_L}^{x_U} f(x) dx\).
- Parameters:
f – [in] The integrand. Can return any type that overloads multiplication with a double and the += operator. doubles and Eigen::VectorXd are examples.
lb – [in] The lower bound \(x_L\) in the integration.
lb – [in] The upper bound \(x_U\) in the integration.
- Template Parameters:
FunctionType – The type of the integrand. Must have an operator()(double x) function.
- Returns:
An approximation of \(\int_{x_L}^{x_U} f(x) dx\). The return type will be the same as the type returned by the integrand function.
-
template<class FunctionType>
inline void Integrate(double *workspace, FunctionType const &f, double lb, double ub, double *res) const#
-
inline AdaptiveClenshawCurtis(unsigned int level, unsigned int maxSub, unsigned int maxDim, double absTol, double relTol, QuadError::Type errorMetric, unsigned int minSub)
-
inline AdaptiveClenshawCurtis(unsigned int level, unsigned int maxSub, unsigned int maxDim, double *workspace, double absTol, double relTol, QuadError::Type errorMetric, unsigned int minSub)
Public Static Functions
-
static inline unsigned int GetWorkspaceSize(unsigned int maxSub, unsigned int fdim)#
-
inline AdaptiveClenshawCurtis(unsigned int level, unsigned int maxSub, unsigned int maxDim, double absTol, double relTol, QuadError::Type errorMetric, unsigned int minSub = 0)#