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
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
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:
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:
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¶
- Reuse knots: Pass the same knot vector to derivative functions
- Batch evaluation: Evaluate multiple points simultaneously
- Appropriate degree: Degree 3 (cubic) is usually sufficient
- Segment selection: More segments increase flexibility but computational cost