scipy.interpolate.UnivariateSpline#
- class scipy.interpolate.UnivariateSpline(x, y, w=None, bbox=[None, None], k=3, s=None, ext=0, check_finite=False)[source]#
- 1-D smoothing spline fit to a given set of data points. - Fits a spline y = spl(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; must be strictly increasing if s is 0. 
- 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 w is None, weights are all 1. Default is None. 
- bbox(2,) array_like, optional
- 2-sequence specifying the boundary of the approximation interval. If bbox is None, - bbox=[x[0], x[-1]]. Default is None.
- kint, optional
- Degree of the smoothing spline. Must be 1 <= k <= 5. - k = 3is a cubic spline. Default is 3.
- sfloat 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]-spl(x[i])))**2, axis=0) <= s - However, because of numerical issues, the actual condition is: - abs(sum((w[i] * (y[i]-spl(x[i])))**2, axis=0) - s) < 0.001 * s - If s is None, s will be set as len(w) for a smoothing spline that uses all data points. If 0, spline will interpolate through all data points. This is equivalent to - InterpolatedUnivariateSpline. Default is None. The user can use the s to control the tradeoff between closeness and smoothness of fit. Larger s means more smoothing while smaller values of s indicate less smoothing. Recommended values of s depend on the weights, w. If the weights represent the inverse of the standard-deviation of y, then a good s value should be found in the range (m-sqrt(2*m),m+sqrt(2*m)) where m is the number of datapoints in x, y, and w. This means- s = len(w)should be a good value if- 1/w[i]is an estimate of the standard deviation of- y[i].
- extint or str, optional
- Controls the extrapolation mode for elements not in the interval defined by the knot sequence. - if ext=0 or ‘extrapolate’, return the extrapolated value. 
- if ext=1 or ‘zeros’, return 0 
- if ext=2 or ‘raise’, raise a ValueError 
- if ext=3 or ‘const’, return the boundary value. 
 - Default is 0. 
- check_finitebool, optional
- Whether to check that the input arrays contain only finite numbers. Disabling may give a performance gain, but may result in problems (crashes, non-termination or non-sensical results) if the inputs do contain infinities or NaNs. Default is False. 
 
 - See also - BivariateSpline
- a base class for bivariate splines. 
- SmoothBivariateSpline
- a smoothing bivariate spline through the given points 
- LSQBivariateSpline
- a bivariate spline using weighted least-squares fitting 
- RectSphereBivariateSpline
- a bivariate spline over a rectangular mesh on a sphere 
- SmoothSphereBivariateSpline
- a smoothing bivariate spline in spherical coordinates 
- LSQSphereBivariateSpline
- a bivariate spline in spherical coordinates using weighted least-squares fitting 
- RectBivariateSpline
- a bivariate spline over a rectangular mesh 
- InterpolatedUnivariateSpline
- a interpolating univariate spline for a given set of data points. 
- bisplrep
- a function to find a bivariate B-spline representation of a surface 
- bisplev
- a function to evaluate a bivariate B-spline and its derivatives 
- splrep
- a function to find the B-spline representation of a 1-D curve 
- splev
- a function to evaluate a B-spline or its derivatives 
- sproot
- a function to find the roots of a cubic B-spline 
- splint
- a function to evaluate the definite integral of a B-spline between two given points 
- spalde
- a function to evaluate all derivatives of a B-spline 
 - Notes - The number of data points must be larger than the spline degree k. - NaN handling: If the input arrays contain - nanvalues, the result is not useful, since the underlying spline fitting routines cannot deal with- nan. A workaround is to use zero weights for not-a-number data points:- >>> import numpy as np >>> from scipy.interpolate import UnivariateSpline >>> x, y = np.array([1, 2, 3, 4]), np.array([1, np.nan, 3, 4]) >>> w = np.isnan(y) >>> y[w] = 0. >>> spl = UnivariateSpline(x, y, w=~w) - Notice the need to replace a - nanby a numerical value (precise value does not matter as long as the corresponding weight is zero.)- References - Based on algorithms described in [1], [2], [3], and [4]: [1]- P. Dierckx, “An algorithm for smoothing, differentiation and integration of experimental data using spline functions”, J.Comp.Appl.Maths 1 (1975) 165-184. [2]- P. Dierckx, “A fast algorithm for smoothing data on a rectangular grid while using spline functions”, SIAM J.Numer.Anal. 19 (1982) 1286-1304. [3]- P. Dierckx, “An improved algorithm for curve fitting with spline functions”, report tw54, Dept. Computer Science,K.U. Leuven, 1981. [4]- P. Dierckx, “Curve and surface fitting with splines”, Monographs on Numerical Analysis, Oxford University Press, 1993. - Examples - >>> import numpy as np >>> import matplotlib.pyplot as plt >>> from scipy.interpolate import UnivariateSpline >>> rng = np.random.default_rng() >>> x = np.linspace(-3, 3, 50) >>> y = np.exp(-x**2) + 0.1 * rng.standard_normal(50) >>> plt.plot(x, y, 'ro', ms=5) - Use the default value for the smoothing parameter: - >>> spl = UnivariateSpline(x, y) >>> xs = np.linspace(-3, 3, 1000) >>> plt.plot(xs, spl(xs), 'g', lw=3) - Manually change the amount of smoothing: - >>> spl.set_smoothing_factor(0.5) >>> plt.plot(xs, spl(xs), 'b', lw=3) >>> plt.show()   - Methods - __call__(x[, nu, ext])- 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. - Return spline coefficients. - Return positions of interior knots of the spline. - Return weighted sum of squared residuals of the spline approximation. - integral(a, b)- Return definite integral of the spline between two given points. - roots()- Return the zeros of the spline. - Continue spline computation with the given smoothing factor s and with the knots found at the last call. - validate_input