{ } Raw JSON

bundles / scipy latest / scipy / interpolate / _bsplines / BSpline

class

scipy.interpolate._bsplines:BSpline

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

Signature

class   BSpline ( t c k extrapolate = True axis = 0 )

Members

Summary

Univariate spline in the B-spline basis.

Extended Summary

where are B-spline basis functions of degree k and knots t.

Parameters

t : ndarray, shape (n+k+1,)

knots

c : ndarray, shape (>=n, ...)

spline coefficients

k : int

B-spline degree

extrapolate : bool or 'periodic', optional

whether to extrapolate beyond the base interval, t[k] .. t[n], or to return nans. If True, extrapolates the first and last polynomial pieces of b-spline functions active on the base interval. If 'periodic', periodic extrapolation is used. Default is True.

axis : int, optional

Interpolation axis. Default is zero.

Attributes

t : ndarray

knot vector

c : ndarray

spline coefficients

k : int

spline degree

extrapolate : bool

If True, extrapolates the first and last polynomial pieces of b-spline functions active on the base interval.

axis : int

Interpolation axis.

tck : tuple

A read-only equivalent of (self.t, self.c, self.k)

Methods

__call__
basis_element
derivative
antiderivative
integrate
insert_knot
construct_fast
design_matrix
from_power_basis

Notes

B-spline basis elements are defined via

Implementation details

  • At least k+1 coefficients are required for a spline of degree k, so that n >= k+1. Additional coefficients, c[j] with j > n, are ignored.

  • B-spline basis elements of degree k form a partition of unity on the base interval, t[k] <= x <= t[n].

Array API Standard Support

BSpline 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                  ⛔                     n/a                 
====================  ====================  ====================

See dev-arrayapi for more information.

Examples

Translating the recursive definition of B-splines into Python code, we have:
def B(x, k, i, t):
   if k == 0:
      return 1.0 if t[i] <= x < t[i+1] else 0.0
   if t[i+k] == t[i]:
      c1 = 0.0
   else:
      c1 = (x - t[i])/(t[i+k] - t[i]) * B(x, k-1, i, t)
   if t[i+k+1] == t[i+1]:
      c2 = 0.0
   else:
      c2 = (t[i+k+1] - x)/(t[i+k+1] - t[i+1]) * B(x, k-1, i+1, t)
   return c1 + c2
def bspline(x, t, c, k):
   n = len(t) - k - 1
   assert (n >= k+1) and (len(c) >= n)
   return sum(c[i] * B(x, k, i, t) for i in range(n))
Note that this is an inefficient (if straightforward) way to evaluate B-splines --- this spline class does it in an equivalent, but much more efficient way. Here we construct a quadratic spline function on the base interval ``2 <= x <= 4`` and compare with the naive way of evaluating the spline:
from scipy.interpolate import BSpline
k = 2
t = [0, 1, 2, 3, 4, 5, 6]
c = [-1, 2, 0, -1]
spl = BSpline(t, c, k)
spl(2.5)
bspline(2.5, t, c, k)
Note that outside of the base interval results differ. This is because `BSpline` extrapolates the first and last polynomial pieces of B-spline functions active on the base interval.
import matplotlib.pyplot as plt
import numpy as np
fig, ax = plt.subplots()
xx = np.linspace(1.5, 4.5, 50)
ax.plot(xx, [bspline(x, t, c ,k) for x in xx], 'r-', lw=3, label='naive')
ax.plot(xx, spl(xx), 'b-', lw=4, alpha=0.7, label='BSpline')
ax.grid(True)
ax.legend(loc='best')
plt.show()
fig-a22e510536011789.png

Aliases

  • scipy.interpolate.BSpline

Referenced by