ClenshawCurtis

class probnum.quad.ClenshawCurtis(npts_per_dim, ndim, bounds)[source]

Bases: probnum.quad.PolynomialQuadrature

Clenshaw-Curtis quadrature rule.

Method of numerical integration based on an expansion of the integrand in terms of a discrete cosine transform.

The nodes of the Clenshaw-Curtis rule are the roots of the Chebyshev polynomials. The \(i^\text{th}\) root is

\[x_i = \frac{1}{2} \left(1 - \cos\left( \frac{i \pi}{n+1} \right) \right)\]

for \(i=1, ..., n\). The \(i^\text{th}\) weight is given by

\[w_i = \frac{2}{n+1} \sin\left(\frac{i \pi}{n+1}\right)\sum_{j=1}^{(n+1)/2} \frac{1}{2j-1}\sin\left(\frac{(2j-1)i \pi}{n+1}\right).\]

These formulas can be found in [1]. For an \(r\)-times differentiable integrand, the Clenshaw-Curtis approximation error is proportional to \(\mathcal{O}(n^{-r})\). It integrates polynomials of degree \(\leq n+1\) exactly.

Parameters:
  • npts_per_dim (int) – Number of evaluation points per dimension. The resulting mesh will have npts_per_dim**ndim elements.
  • ndim (int) – Number of dimensions.
  • bounds (ndarray, shape=(n, 2)) – Integration bounds.

See also

PolynomialQuadrature
Quadrature rule based on polynomial functions.

References

[1]Holtz, M., Sparse Grid Quadrature in High Dimensions with Applications in Finance and Insurance, Springer, 2010

Examples

>>> cc = ClenshawCurtis(npts_per_dim=3, ndim=2, bounds=np.array([[0, 1], [0, 0.1]]))
>>> print(cc.nodes)
[[0.14644661 0.01464466]
 [0.14644661 0.05      ]
 [0.14644661 0.08535534]
 [0.5        0.01464466]
 [0.5        0.05      ]
 [0.5        0.08535534]
 [0.85355339 0.01464466]
 [0.85355339 0.05      ]
 [0.85355339 0.08535534]]
>>> print(cc.weights)
[0.01111111 0.01111111 0.01111111 0.01111111 0.01111111 0.01111111
 0.01111111 0.01111111 0.01111111]
>>> print(cc.integrate(lambda x: x[0] + x[1]))
0.05500000000000001
>>> cc = ClenshawCurtis(npts_per_dim=7, ndim=1, bounds=np.array([[0, 1]]))
>>> print(cc.weights)
[0.08898234 0.12380952 0.19673195 0.18095238 0.19673195 0.12380952
 0.08898234]
>>> print(cc.nodes)
[[0.03806023]
 [0.14644661]
 [0.30865828]
 [0.5       ]
 [0.69134172]
 [0.85355339]
 [0.96193977]]
>>> print(cc.integrate(lambda x: np.sin(x)))
[0.45969769]

Methods Summary

integrate(fun[, isvectorized]) Numerically integrate the function fun.

Methods Documentation

integrate(fun, isvectorized=False, **kwargs)

Numerically integrate the function fun.

Parameters:
  • fun (callable) – Function to be integrated. Signature is fun(x, **kwargs) where x is either a float or an ndarray with shape (d,). If fun can be evaluated vectorized, the implementation expects signature fun(X, **kwargs) where X is is an ndarray of shape (n, d). Making use of vectorization is recommended wherever possible for improved speed of computation.
  • isvectorized (bool) – Whether integrand allows vectorised evaluation (i.e. evaluation of all nodes at once).
  • kwargs – Key-word arguments being passed down to fun at each evaluation. For example (hyper)parameters.