chaospy.quadrature.radau

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

Generate the quadrature nodes and weights in Gauss-Radau quadrature.

Gauss-Radau formula for numerical estimation of integrals. It requires \(m+1\) points and fits all Polynomials to degree \(2m\), so it effectively fits exactly all Polynomials of degree \(2m-1\).

It allows for a single abscissas to be user defined, while the others are built around this point.

Canonically, Radau is built around Legendre weight function with the fixed point at the left end. Not all distributions/fixed point combinations allows for the building of a quadrature scheme.

Args:
order (int):

Quadrature order.

dist (chaospy.Distribution):

The distribution weights to be used to create higher order nodes from.

fixed_point (float):

Fixed point abscissas assumed to be included in the quadrature. If omitted, use distribution lower bound.

rule (str):

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

accuracy (int):

In the case rule is used, defines the quadrature order of the scheme used. In practice, must be at least as large as order.

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.

Returns:
abscissas (numpy.ndarray):

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

weights (numpy.ndarray):

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

Example:
>>> distribution = chaospy.Uniform(-1, 1)
>>> abscissas, weights = chaospy.quadrature.radau(4, distribution)
>>> abscissas.round(3)
array([[-1.   , -0.887, -0.64 , -0.295,  0.094,  0.468,  0.771,  0.955]])
>>> weights.round(3)
array([0.016, 0.093, 0.152, 0.188, 0.196, 0.174, 0.125, 0.057])