chaospy.quadrature.kronrod

chaospy.quadrature.kronrod(order, dist, recurrence_algorithm='stieltjes', rule='clenshaw_curtis', tolerance=1e-10, scaling=3, n_max=5000)[source]

Generate Gauss-Kronrod quadrature abscissas and weights.

Gauss-Kronrod quadrature is an adaptive method for Gaussian quadrature rule. It builds on top of other quadrature rules by extending “pure” Gaussian quadrature rules with extra abscissas and new weights such that already used abscissas can be reused. For more details, see Wikipedia article.

For each order N taken with ordinary Gaussian quadrature, Gauss-Kronrod will create 2N+1 abscissas where all of the N “old” abscissas are all interlaced between the “new” ones.

The algorithm is well suited for any Jacobi scheme, i.e. quadrature involving Uniform or Beta distribution, and might work on others as well. However, it will not work everywhere. For example Kahaner and Monegato showed that higher order Gauss-Kronrod quadrature for Gauss-Hermite and Gauss-Laguerre does not exist.

Args:
order (int):

The order of the quadrature.

dist (chaospy.distributions.baseclass.Distribution):

The distribution which density will be used as weight function.

recurrence_algorithm (str):

Name of the algorithm used to generate abscissas and weights. If omitted, analytical will be tried first, and stieltjes used if that fails.

rule (str):

In the case of lanczos or stieltjes, defines the proxy-integration scheme.

tolerance (float):

The allowed relative error in norm between two quadrature orders before method assumes convergence.

scaling (float):

A multiplier the adaptive order increases with for each step quadrature order is not converged. Use 0 to indicate unit increments.

n_max (int):

The allowed number of quadrature points to use in approximation.

Returns:
(numpy.ndarray, numpy.ndarray):
abscissas:

The quadrature points for where to evaluate the model function with abscissas.shape == (len(dist), N) where N is the number of samples.

weights:

The quadrature weights with weights.shape == (N,).

Raises:
ValueError:

Error raised if Kronrod algorithm results in negative recurrence coefficients.

Notes:

Code is adapted from quadpy, which adapted his code from W. Gautschi. Algorithm for calculating Kronrod-Jacobi matrices was first published Algorithm as proposed by Laurie [5].

Example:
>>> distribution = chaospy.Uniform(-1, 1)
>>> abscissas, weights = chaospy.quadrature.kronrod(4, distribution)
>>> abscissas.round(2)
array([[-0.98, -0.86, -0.64, -0.34,  0.  ,  0.34,  0.64,  0.86,  0.98]])
>>> weights.round(3)
array([0.031, 0.085, 0.133, 0.163, 0.173, 0.163, 0.133, 0.085, 0.031])