scipy.interpolate.KroghInterpolator#

class scipy.interpolate.KroghInterpolator(xi, yi, axis=0)[source]#

Interpolating polynomial for a set of points.

The polynomial passes through all the pairs (xi, yi). One may additionally specify a number of derivatives at each point xi; this is done by repeating the value xi and specifying the derivatives as successive yi values.

Allows evaluation of the polynomial and all its derivatives. For reasons of numerical stability, this function does not compute the coefficients of the polynomial, although they can be obtained by evaluating all the derivatives.

Parameters:
xiarray_like, shape (npoints, )

Known x-coordinates. Must be sorted in increasing order.

yiarray_like, shape (…, npoints, …)

Known y-coordinates. When an xi occurs two or more times in a row, the corresponding yi’s represent derivative values. The length of yi along the interpolation axis must be equal to the length of xi. Use the axis parameter to select the correct axis.

axisint, optional

Axis in the yi array corresponding to the x-coordinate values. Defaults to axis=0.

Notes

Be aware that the algorithms implemented here are not necessarily the most numerically stable known. Moreover, even in a world of exact computation, unless the x coordinates are chosen very carefully - Chebyshev zeros (e.g., cos(i*pi/n)) are a good choice - polynomial interpolation itself is a very ill-conditioned process due to the Runge phenomenon. In general, even with well-chosen x values, degrees higher than about thirty cause problems with numerical instability in this code.

Based on [1].

References

[1]

Krogh, “Efficient Algorithms for Polynomial Interpolation and Numerical Differentiation”, 1970.

Examples

To produce a polynomial that is zero at 0 and 1 and has derivative 2 at 0, call

>>> from scipy.interpolate import KroghInterpolator
>>> KroghInterpolator([0,0,1],[0,2,0])

This constructs the quadratic \(2x^2-2x\). The derivative condition is indicated by the repeated zero in the xi array; the corresponding yi values are 0, the function value, and 2, the derivative value.

For another example, given xi, yi, and a derivative ypi for each point, appropriate arrays can be constructed as:

>>> import numpy as np
>>> rng = np.random.default_rng()
>>> xi = np.linspace(0, 1, 5)
>>> yi, ypi = rng.random((2, 5))
>>> xi_k, yi_k = np.repeat(xi, 2), np.ravel(np.dstack((yi,ypi)))
>>> KroghInterpolator(xi_k, yi_k)

To produce a vector-valued polynomial, supply a higher-dimensional array for yi:

>>> KroghInterpolator([0,1],[[2,3],[4,5]])

This constructs a linear polynomial giving (2,3) at 0 and (4,5) at 1.

Attributes:
dtype

Methods

__call__(x)

Evaluate the interpolant

derivative(x[, der])

Evaluate a single derivative of the polynomial at the point x.

derivatives(x[, der])

Evaluate several derivatives of the polynomial at the point x