chaospy.quadrature.lobatto¶
- chaospy.quadrature.lobatto(order, dist, recurrence_algorithm='stieltjes', rule='fejer_2', tolerance=1e-10, scaling=3, n_max=5000)[source]¶
Generate the abscissas and weights in Gauss-Loboto quadrature.
Also known as Lobatto quadrature, named after Dutch mathematician Rehuel Lobatto. It is similar to Gaussian quadrature with the following differences:
The integration points include the end points of the integration interval.
It is accurate for polynomials up to degree \(2n–3\), where \(n\) is the number of integration points.
- Args:
- order (int):
Quadrature order.
- dist (chaospy.distributions.baseclass.Distribution):
The distribution weights to be used to create higher order nodes from.
- recurrence_algorithm (str):
Name of the algorithm used to generate abscissas and weights.
- rule (str):
In the case of
lanczos
orstieltjes
, 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)
whereN
is the number of samples.- weights:
The quadrature weights with
weights.shape == (N,)
.
- Example:
>>> distribution = chaospy.Uniform(-1, 1) >>> abscissas, weights = chaospy.quadrature.lobatto(4, distribution) >>> abscissas.round(3) array([[-1. , -0.872, -0.592, -0.209, 0.209, 0.592, 0.872, 1. ]]) >>> weights.round(3) array([0.018, 0.105, 0.171, 0.206, 0.206, 0.171, 0.105, 0.018])