SciPy

scipy.interpolate.UnivariateSpline

class scipy.interpolate.UnivariateSpline(x, y, w=None, bbox=[None, None], k=3, s=None)[source]

One-dimensional smoothing spline fit to a given set of data points.

Fits a spline y=s(x) of degree k to the provided x, y data. s specifies the number of knots by specifying a smoothing condition.

Parameters:

x : (N,) array_like

1-D array of independent input data. Must be increasing.

y : (N,) array_like

1-D array of dependent input data, of the same length as x.

w : (N,) array_like, optional

Weights for spline fitting. Must be positive. If None (default), weights are all equal.

bbox : (2,) array_like, optional

2-sequence specifying the boundary of the approximation interval. If None (default), bbox=[x[0], x[-1]].

k : int, optional

Degree of the smoothing spline. Must be <= 5.

s : float or None, optional

Positive smoothing factor used to choose the number of knots. Number of knots will be increased until the smoothing condition is satisfied:

sum((w[i]*(y[i]-s(x[i])))**2,axis=0) <= s

If None (default), s=len(w) which should be a good value if 1/w[i] is an estimate of the standard deviation of y[i]. If 0, spline will interpolate through all data points.

See also

InterpolatedUnivariateSpline
Subclass with smoothing forced to 0
LSQUnivariateSpline
Subclass in which knots are user-selected instead of being set by smoothing condition
splrep
An older, non object-oriented wrapping of FITPACK

splev, sproot, splint, spalde

BivariateSpline
A similar class for two-dimensional spline interpolation

Notes

The number of data points must be larger than the spline degree k.

Examples

>>> from numpy import linspace,exp
>>> from numpy.random import randn
>>> import matplotlib.pyplot as plt
>>> from scipy.interpolate import UnivariateSpline
>>> x = linspace(-3, 3, 100)
>>> y = exp(-x**2) + randn(100)/10
>>> s = UnivariateSpline(x, y, s=1)
>>> xs = linspace(-3, 3, 1000)
>>> ys = s(xs)
>>> plt.plot(x, y, '.-')
>>> plt.plot(xs, ys)
>>> plt.show()

xs,ys is now a smoothed, super-sampled version of the noisy gaussian x,y.

(Source code)

Methods

__call__(x[, nu]) Evaluate spline (or its nu-th derivative) at positions x.
antiderivative([n]) Construct a new spline representing the antiderivative of this spline.
derivative([n]) Construct a new spline representing the derivative of this spline.
derivatives(x) Return all derivatives of the spline at the point x.
get_coeffs() Return spline coefficients.
get_knots() Return positions of (boundary and interior) knots of the spline.
get_residual() Return weighted sum of squared residuals of the spline approximation: sum((w[i] * (y[i]-s(x[i])))**2, axis=0).
integral(a, b) Return definite integral of the spline between two given points.
roots() Return the zeros of the spline.
set_smoothing_factor(s) Continue spline computation with the given smoothing factor s and with the knots found at the last call.