Interpolation methods
Method groups
B-Splines
Interpolate with 3D tensor-product B-Spline basis.
- nitransforms.interp.bspline.grid_bspline_weights(target_grid, ctrl_grid)[source]
Evaluate tensor-product B-Spline weights on a grid.
For each of the \(N\) input locations \(\mathbf{x} = (x_i, x_j, x_k)\) and \(K\) control points or knots \(\mathbf{c} =(c_i, c_j, c_k)\), the tensor-product cubic B-Spline kernel weights are calculated:
\[\Psi^3(\mathbf{x}, \mathbf{c}) = \beta^3(x_i - c_i) \cdot \beta^3(x_j - c_j) \cdot \beta^3(x_k - c_k), \label{eq:bspline_weights}\tag{1}\]where each \(\beta^3\) represents the cubic B-Spline for one dimension. The 1D B-Spline kernel implementation uses
numpy.piecewise
, and is based on the closed-form given by Eq. (6) of [Unser1999].By iterating over dimensions, the data samples that fall outside of the compact support of the tensor-product kernel associated to each control point can be filtered out and dismissed to lighten computation.
Finally, the resulting weights matrix \(\Psi^3(\mathbf{k}, \mathbf{s})\) can be easily identified in Eq. \(\eqref{eq:1}\) and used as the design matrix for approximation of data.
- Parameters:
target_grid (
ImageGrid
ornibabel.spatialimages
) – Regular grid of \(N\) locations at which tensor B-Spline basis will be evaluated.ctrl_grid (
ImageGrid
ornibabel.spatialimages
) – Regular grid of \(K\) control points (knot) where B-Spline basis are defined.
- Returns:
weights – A sparse matrix of interpolating weights \(\Psi^3(\mathbf{k}, \mathbf{s})\) for the N voxels of the target EPI, for each of the total K knots. This sparse matrix can be directly used as design matrix for the fitting step of approximation/extrapolation.
- Return type:
numpy.ndarray
(\(K \times N\))