chaospy.quadrature.legendre_proxy

chaospy.quadrature.legendre_proxy(order, domain=(0, 1), segments=1)[source]

Generate proxy abscissas and weights from Legendre quadrature.

Legendre provides optimal abscissas \(X_i\) and weights \(W_i\) to solve the integration problem:

\[\int f(x) dx \approx \sum W_i f(X_i)\]

over a function \(f\), where the probability density function \(p\) is uniform.

Since the weight function is constant, it can in principle be used to integrate any density function by considering it a part of the function. In other words:

\[\int p(x) f(x) \approx \sum W_i p(X_i) f(X_i) = \sum W_i^' f(X_i)\]

So when providing non-uniform distribution as domain, the weights will be adjusted with:

\[W_i^' = W_i p(X_i)\]

Bounds of the Legendre schemes is chosen to be the same as the distribution provided. This makes it a bad choice for unbound distributions.

To get optimal abscissas and weights directly from a density, use chaospy.quadrature.gaussian() instead.

Args:
order (int, numpy.ndarray):

Quadrature order.

domain (chaospy.Distribution, numpy.ndarray):

Either distribution or bounding of interval to integrate over.

segments (int):

Split intervals into steps subintervals and create a patched quadrature based on the segmented quadrature. Can not be lower than order. If 0 is provided, default to square root of order. Nested samples only appear when the number of segments are fixed.

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:
>>> abscissas, weights = chaospy.quadrature.legendre_proxy(3)
>>> abscissas.round(4)
array([[0.0694, 0.33  , 0.67  , 0.9306]])
>>> weights.round(4)
array([0.1739, 0.3261, 0.3261, 0.1739])
>>> abscissas, weights = chaospy.quadrature.legendre_proxy(3, chaospy.Uniform(0, 1))
>>> abscissas.round(4)
array([[0.0694, 0.33  , 0.67  , 0.9306]])
>>> weights.round(4)
array([0.1739, 0.3261, 0.3261, 0.1739])
>>> abscissas, weights = chaospy.quadrature.legendre_proxy(3, chaospy.Beta(2, 2))
>>> abscissas.round(4)
array([[0.0694, 0.33  , 0.67  , 0.9306]])
>>> weights.round(4)
array([0.0674, 0.4326, 0.4326, 0.0674])
See also:

chaospy.quadrature.legendre(), chaospy.quadrature.gaussian()