{ } Raw JSON

bundles / scipy latest / scipy / interpolate / _bsplines / make_smoothing_spline

function

scipy.interpolate._bsplines:make_smoothing_spline

source: /scipy/interpolate/_bsplines.py :2248

Signature

def   make_smoothing_spline ( x y w = None lam = None * axis = 0 )

Summary

Create a smoothing B-spline satisfying the Generalized Cross Validation (GCV) criterion.

Extended Summary

Compute the (coefficients of) smoothing cubic spline function using lam to control the tradeoff between the amount of smoothness of the curve and its proximity to the data. In case lam is None, using the GCV criteria [1] to find it.

A smoothing spline is found as a solution to the regularized weighted linear regression problem:

where is a spline function, is a vector of weights and is a regularization parameter.

If lam is None, we use the GCV criteria to find an optimal regularization parameter, otherwise we solve the regularized weighted linear regression problem with given parameter. The parameter controls the tradeoff in the following way: the larger the parameter becomes, the smoother the function gets.

Parameters

x : array_like, shape (n,)

Abscissas. n must be at least 5.

y : array_like, shape (n, ...)

Ordinates. n must be at least 5.

w : array_like, shape (n,), optional

Vector of weights. Default is np.ones_like(x).

lam : float, (:math:`\lambda \geq 0`), optional

Regularization parameter. If lam is None, then it is found from the GCV criteria. Default is None.

axis : int, optional

The data axis. Default is zero. The assumption is that y.shape[axis] == n, and all other axes of y are batching axes.

Returns

func : `BSpline` object

An object representing a spline in the B-spline basis as a solution of the problem of smoothing splines using the GCV criteria [1] in case lam is None, otherwise using the given parameter lam.

Notes

This algorithm is a clean room reimplementation of the algorithm introduced by Woltring in FORTRAN [2]. The original version cannot be used in SciPy source code because of the license issues. The details of the reimplementation are discussed here (available only in Russian) [4].

If the vector of weights w is None, we assume that all the points are equal in terms of weights, and vector of weights is vector of ones.

Note that in weighted residual sum of squares, weights are not squared: while in splrep the sum is built from the squared weights.

In cases when the initial problem is ill-posed (for example, the product where is a design matrix is not a positive defined matrix) a ValueError is raised.

Array API Standard Support

make_smoothing_spline has experimental support for Python Array API Standard compatible backends in addition to NumPy. Please consider testing these features by setting an environment variable SCIPY_ARRAY_API=1 and providing CuPy, PyTorch, JAX, or Dask arrays as array arguments. The following combinations of backend and device (or other capability) are supported.

====================  ====================  ====================
Library               CPU                   GPU
====================  ====================  ====================
NumPy                 ✅                     n/a                 
CuPy                  n/a                   ⛔                   
PyTorch               ✅                     ⛔                   
JAX                   ⚠️ no JIT
Dask                  ⚠️ computes graph     n/a                 
====================  ====================  ====================

See dev-arrayapi for more information.

Examples

Generate some noisy data
import numpy as np
np.random.seed(1234)
n = 200
def func(x):
   return x**3 + x**2 * np.sin(4 * x)
x = np.sort(np.random.random_sample(n) * 4 - 2)
y = func(x) + np.random.normal(scale=1.5, size=n)
Make a smoothing spline function
from scipy.interpolate import make_smoothing_spline
spl = make_smoothing_spline(x, y)
Plot both
import matplotlib.pyplot as plt
grid = np.linspace(x[0], x[-1], 400)
plt.plot(x, y, '.')
plt.plot(grid, spl(grid), label='Spline')
plt.plot(grid, func(grid), label='Original function')
plt.legend(loc='best')
plt.show()
fig-b77f48254baa5213.png

Aliases

  • scipy.interpolate.make_smoothing_spline

Referenced by