Skip to content

Basis Functions API

The basis module provides functions for constructing B-spline basis matrices and their derivatives.

Functions

Compute B-spline regression basis on [xl, xr].

Parameters:

  • x ((array - like, shape(n))) –
  • xl (floats) –

    Domain endpoints (xl < xr).

  • xr (floats) –

    Domain endpoints (xl < xr).

  • nseg (int) –

    Number of equal-length segments.

  • degree (int, default: 3 ) –

    Degree of the spline.

Returns:

  • B ( (csr_matrix, shape(n, nseg + degree)) ) –

    Sparse design matrix, each row has at most degree+1 nonzeros.

  • knots ( ndarray ) –

    Full knot vector (length nseg + degree + 1 + degree).

Source code in src/psplines/basis.py
def b_spline_basis(
    x: ArrayLike, xl: float, xr: float, nseg: int, degree: int = 3
) -> tuple[sp.csr_matrix, NDArray]:
    """
    Compute B-spline regression basis on [xl, xr].

    Parameters
    ----------
    x : array-like, shape (n,)
    xl, xr : floats
        Domain endpoints (xl < xr).
    nseg : int
        Number of equal-length segments.
    degree : int
        Degree of the spline.

    Returns
    -------
    B : csr_matrix, shape (n, nseg+degree)
        Sparse design matrix, each row has at most degree+1 nonzeros.
    knots : ndarray
        Full knot vector (length nseg + degree + 1 + degree).
    """
    x = np.asarray(x, float)
    knots = _make_knots(xl, xr, nseg, degree)
    # number of basis functions
    nb = nseg + degree
    # clip to domain
    x_clipped = np.clip(x, xl, xr)
    # build BSpline object: vectorized multi-coeff basis
    # coefficients as identity to extract each basis
    coeffs = np.eye(nb)
    spline = BSpline(knots, coeffs, degree, extrapolate=False)
    # evaluate basis at all x in one call
    B_full = spline(x_clipped)
    # convert to CSR for sparsity
    B = sp.csr_matrix(B_full)
    return B, knots

Compute k-th derivative of B-spline basis on [xl, xr].

Parameters:

  • x (array - like) –
  • xl (floats) –
  • xr (floats) –
  • nseg (int) –
  • degree (int, default: 3 ) –
  • deriv_order (int, default: 1 ) –

    Order of derivative (0 returns same as b_spline_basis).

  • knots (ndarray or None, default: None ) –

    Precomputed knot vector.

Returns:

  • B_deriv ( (csr_matrix, shape(n, nseg + degree)) ) –

    Sparse derivative basis.

  • knots ( ndarray ) –
Source code in src/psplines/basis.py
def b_spline_derivative_basis(
    x: ArrayLike,
    xl: float,
    xr: float,
    nseg: int,
    degree: int = 3,
    deriv_order: int = 1,
    knots: NDArray | None = None,
) -> tuple[sp.csr_matrix, NDArray]:
    """
    Compute k-th derivative of B-spline basis on [xl, xr].

    Parameters
    ----------
    x : array-like
    xl, xr : floats
    nseg : int
    degree : int
    deriv_order : int
        Order of derivative (0 returns same as b_spline_basis).
    knots : ndarray or None
        Precomputed knot vector.

    Returns
    -------
    B_deriv : csr_matrix, shape (n, nseg+degree)
        Sparse derivative basis.
    knots : ndarray
    """
    if deriv_order < 0:
        raise ValueError("deriv_order must be non-negative")
    if knots is None:
        knots = _make_knots(xl, xr, nseg, degree)
    # nothing to do for zero-order
    if deriv_order == 0:
        return b_spline_basis(x, xl, xr, nseg, degree)
    # if derivative exceeds degree, result is zero
    if deriv_order > degree:
        n = np.atleast_1d(x).shape[0]
        return sp.csr_matrix((n, nseg + degree)), knots

    x = np.asarray(x, float)
    x_clipped = np.clip(x, xl, xr)
    nb = nseg + degree
    # build multi-coefficient BSpline
    coeffs = np.eye(nb)
    spline = BSpline(knots, coeffs, degree, extrapolate=False)
    # evaluate derivative in one call
    B_deriv_full = spline(x_clipped, nu=deriv_order)
    B_deriv = sp.csr_matrix(B_deriv_full)
    return B_deriv, knots

Usage Examples

Basic B-spline Basis

import numpy as np
from psplines.basis import b_spline_basis

# Create evaluation points
x = np.linspace(0, 1, 50)

# Generate B-spline basis matrix
B, knots = b_spline_basis(x, xl=0, xr=1, nseg=10, degree=3)

print(f"Basis shape: {B.shape}")  # (50, 13) for degree=3, nseg=10
print(f"Number of knots: {len(knots)}")  # 17 knots total

Derivative Basis

from psplines.basis import b_spline_derivative_basis

# Generate first derivative basis
B_deriv, knots = b_spline_derivative_basis(
    x, xl=0, xr=1, nseg=10, degree=3, 
    derivative_order=1, knots=knots
)

print(f"Derivative basis shape: {B_deriv.shape}")

Mathematical Background

B-spline Basis Construction

B-splines of degree \(d\) are defined recursively:

\[B_{i,0}(x) = \begin{cases} 1 & \text{if } t_i \leq x < t_{i+1} \\ 0 & \text{otherwise} \end{cases}\]
\[B_{i,d}(x) = \frac{x - t_i}{t_{i+d} - t_i} B_{i,d-1}(x) + \frac{t_{i+d+1} - x}{t_{i+d+1} - t_{i+1}} B_{i+1,d-1}(x)\]

where \(t_i\) are the knots.

Knot Vector Construction

For \(n\) segments and degree \(d\):

  • Interior knots: \(n+1\) equally spaced points
  • Boundary knots: \(d\) repeated knots at each end
  • Total knots: \(n + 1 + 2d\)
  • Number of basis functions: \(n + d\)

Properties

  • Local support: Each basis function is non-zero over at most \(d+1\) knot spans
  • Partition of unity: \(\sum_i B_{i,d}(x) = 1\) for all \(x\)
  • Non-negativity: \(B_{i,d}(x) \geq 0\) for all \(i, x\)
  • Smoothness: \(B_{i,d}(x)\) is \(C^{d-1}\) continuous

Derivative Computation

The \(k\)-th derivative of a B-spline is computed using the recursive formula:

\[\frac{d^k}{dx^k} B_{i,d}(x) = \frac{d!}{(d-k)!} \sum_{j=0}^k (-1)^{k-j} \binom{k}{j} \frac{B_{i+j,d-k}(x)}{(t_{i+d+1-j} - t_{i+j})^k}\]

Implementation Details

Sparse Matrix Format

  • All basis matrices are returned as scipy.sparse.csr_matrix
  • Efficient storage for large problems
  • Fast matrix-vector operations

Numerical Considerations

  • Uses SciPy's BSpline class for robust evaluation
  • Handles edge cases at boundaries
  • Maintains numerical stability for high-degree splines

Performance Tips

  1. Reuse knots: Pass the same knot vector to derivative functions
  2. Batch evaluation: Evaluate multiple points simultaneously
  3. Appropriate degree: Degree 3 (cubic) is usually sufficient
  4. Segment selection: More segments increase flexibility but computational cost