API Reference

Public API

interpolation — standalone interpolation & spline utilities.

All public symbols are re-exported from interpolation.core so callers can simply write:

from interpolation import cubic_spline_not_a_knot, bspline_interpolate

This package is a self-contained extraction of the interpolation utilities originally developed as part of the ITDS gas-turbine design platform.

class interpolation.BSplineRepresentation(knots, coefficients, order)

Bases: object

Univariate B-spline representation.

Parameters:
  • knots (ndarray)

  • coefficients (ndarray)

  • order (int)

knots

Knot vector of shape (n_coeff + order,).

Type:

np.ndarray

coefficients

B-spline coefficients, shape (n_coeff,).

Type:

np.ndarray

order

Spline order (degree + 1).

Type:

int

class interpolation.PiecewisePolynomial(breakpoints, coefficients, order)

Bases: object

Piecewise polynomial representation on contiguous intervals.

Parameters:
  • breakpoints (ndarray)

  • coefficients (ndarray)

  • order (int)

breakpoints

Sorted interval boundaries, shape (n+1,).

Type:

np.ndarray

coefficients

Per-interval coefficients in ascending power order, shape (n, order).

Type:

np.ndarray

order

Polynomial order (degree + 1).

Type:

int

class interpolation.SplineConstraint(derivative, point, type, value)

Bases: object

1D spline linear constraint metadata.

Parameters:
  • derivative (int)

  • point (float)

  • type (str)

  • value (float)

derivative

Derivative order constrained at point.

Type:

int

point

Constraint location.

Type:

float

type

Constraint type: '==', '>=', or '<='.

Type:

str

value

Constraint target value.

Type:

float

class interpolation.SplineKnots(degree, breakpoints)

Bases: object

Spline knot metadata for high-level fitting APIs.

Parameters:
  • degree (int)

  • breakpoints (ndarray)

degree

Polynomial degree.

Type:

int

breakpoints

Knot breakpoints.

Type:

np.ndarray

class interpolation.SurfaceConstraint(dx_order, dy_order, point, type, value)

Bases: object

2D spline linear constraint metadata.

Parameters:
  • dx_order (int)

  • dy_order (int)

  • point (tuple)

  • type (str)

  • value (float)

dx_order

Partial derivative order in x.

Type:

int

dy_order

Partial derivative order in y.

Type:

int

point

Constraint point (x, y).

Type:

tuple

type

Constraint type: '==', '>=', or '<='.

Type:

str

value

Constraint target value.

Type:

float

interpolation.piecewise_poly_derivative(pp, x, k=1)

Evaluate derivative of piecewise polynomial at one point.

imsl_name: PPDER

Parameters:
  • pp (PiecewisePolynomial) – Piecewise polynomial representation.

  • x (float) – Query coordinate.

  • k (int) – Derivative order.

Returns:

Derivative value.

Return type:

float

Raises:

ValueError – If derivative order is negative.

interpolation.piecewise_poly_evaluate(pp, x)

Evaluate piecewise polynomial at one point.

imsl_name: PPVAL

Parameters:
  • pp (PiecewisePolynomial) – Piecewise polynomial representation.

  • x (float) – Query coordinate.

Returns:

Polynomial value.

Return type:

float

interpolation.piecewise_poly_evaluate_grid(pp, points)

Evaluate piecewise polynomial on a vector of points.

imsl_name: PP1GD

Parameters:
  • pp (PiecewisePolynomial) – Piecewise polynomial representation.

  • points (Sequence[float]) – Query coordinates.

Returns:

Evaluated values.

Return type:

np.ndarray

interpolation.piecewise_poly_integrate(pp, a, b)

Integrate piecewise polynomial over [a, b].

imsl_name: PPITG

Parameters:
  • pp (PiecewisePolynomial) – Piecewise polynomial representation.

  • a (float) – Lower bound.

  • b (float) – Upper bound.

Returns:

Definite integral value.

Return type:

float

interpolation.cubic_spline_akima(x, y)

Construct Akima (CSAKM-style) C1 cubic spline interpolation.

imsl_name: CSAKM

Computes a shape-aware piecewise-cubic interpolant using Akima’s weighted local slopes to suppress oscillatory wiggles near steep gradient changes. The spline is C1 continuous and uses the input abscissas as breakpoints. As in the IMSL documentation, this method is nonlinear: linear data are reproduced exactly, while generic cubic polynomials are not guaranteed to be reproduced exactly.

Parameters:
  • x (Sequence[float]) – Strictly increasing x coordinates.

  • y (Sequence[float]) – Data values.

Returns:

Akima cubic spline representation with order=4 and breakpoints == x.

Return type:

PiecewisePolynomial

Raises:

ValueError – If x/y inputs are invalid or x is not strictly increasing.

interpolation.cubic_spline_derivative(pp, x, k=1)

Evaluate k-th derivative of cubic spline.

imsl_name: CSDER

Parameters:
  • pp (PiecewisePolynomial) – Piecewise polynomial spline.

  • x (float) – Query coordinate.

  • k (int) – Derivative order.

Returns:

Derivative value.

Return type:

float

interpolation.cubic_spline_derivative_bc(x, y, dydx_start, dydx_end)

Construct clamped cubic spline with endpoint derivative conditions.

imsl_name: CSDEC

Parameters:
  • x (Sequence[float]) – Strictly increasing x coordinates.

  • y (Sequence[float]) – Data values.

  • dydx_start (float) – First derivative at x[0].

  • dydx_end (float) – First derivative at x[-1].

Returns:

Clamped cubic spline.

Return type:

PiecewisePolynomial

interpolation.cubic_spline_easy(x, y, points)

Interpolate with cubic not-a-knot and evaluate in one call.

imsl_name: CSIEZ

Parameters:
  • x (Sequence[float]) – Strictly increasing x coordinates.

  • y (Sequence[float]) – Data values.

  • points (Sequence[float]) – Query coordinates.

Returns:

Interpolated values.

Return type:

np.ndarray

interpolation.cubic_spline_evaluate(pp, x)

Evaluate cubic spline at one point.

imsl_name: CSVAL

Parameters:
Returns:

Spline value.

Return type:

float

interpolation.cubic_spline_evaluate_grid(pp, points)

Evaluate cubic spline at a vector of query points.

imsl_name: CS1GD

Parameters:
  • pp (PiecewisePolynomial) – Piecewise polynomial spline.

  • points (Sequence[float]) – Query coordinates.

Returns:

Evaluated values.

Return type:

np.ndarray

interpolation.cubic_spline_hermite(x, y, dydx)

Construct cubic Hermite spline from values and derivatives.

imsl_name: CSHER

Parameters:
  • x (Sequence[float]) – Strictly increasing knot coordinates.

  • y (Sequence[float]) – Knot values.

  • dydx (Sequence[float]) – Knot derivatives.

Returns:

Cubic Hermite piecewise polynomial.

Return type:

PiecewisePolynomial

Raises:

ValueError – If derivative vector length mismatches data length.

interpolation.cubic_spline_integrate(pp, a, b)

Integrate cubic spline over [a, b].

imsl_name: CSITG

Parameters:
  • pp (PiecewisePolynomial) – Piecewise polynomial spline.

  • a (float) – Lower integration bound.

  • b (float) – Upper integration bound.

Returns:

Definite integral value.

Return type:

float

interpolation.cubic_spline_not_a_knot(x, y)

Construct not-a-knot-like cubic spline interpolation.

imsl_name: CSINT

Parameters:
  • x (Sequence[float]) – Strictly increasing x coordinates.

  • y (Sequence[float]) – Data values.

Returns:

Cubic piecewise polynomial interpolation.

Return type:

PiecewisePolynomial

interpolation.cubic_spline_periodic(x, y)

Construct periodic cubic spline approximation.

imsl_name: CSPER

Parameters:
  • x (Sequence[float]) – Strictly increasing x coordinates.

  • y (Sequence[float]) – Data values where y[0] ~= y[-1].

Returns:

Periodic cubic spline approximation.

Return type:

PiecewisePolynomial

Raises:

ValueError – If periodic endpoint values do not match closely.

interpolation.cubic_spline_shape_preserving(x, y)

Construct monotone shape-preserving cubic spline (Fritsch-Carlson).

imsl_name: CSCON

Parameters:
  • x (Sequence[float]) – Strictly increasing x coordinates.

  • y (Sequence[float]) – Data values.

Returns:

Monotone cubic spline.

Return type:

PiecewisePolynomial

interpolation.cubic_spline_smooth(x, y, lam, weights=None)

Compute cubic smoothing spline using Tikhonov-style regularization.

imsl_name: CSSMH

Parameters:
  • x (Sequence[float]) – Strictly increasing x coordinates.

  • y (Sequence[float]) – Data values.

  • lam (float) – Smoothing parameter lambda (>= 0).

  • weights (Sequence[float] | None) – Optional non-negative weights.

Returns:

Smoothed cubic spline representation.

Return type:

PiecewisePolynomial

Raises:

ValueError – If lambda is negative or inputs are invalid.

interpolation.cubic_spline_smooth_cv(x, y, lambdas=None, weights=None)

Compute smoothing spline with GCV-selected lambda.

imsl_name: CSSCV

Parameters:
  • x (Sequence[float]) – Strictly increasing x coordinates.

  • y (Sequence[float]) – Data values.

  • lambdas (Sequence[float] | None) – Candidate lambda values.

  • weights (Sequence[float] | None) – Optional sample weights.

Returns:

Keys lambda, spline, gcv_scores.

Return type:

dict

interpolation.cubic_spline_smooth_error(x, y, error_scale=1.0, weights=None)

Compute smoothing spline with automatic lambda from residual scale.

imsl_name: CSSED

Parameters:
  • x (Sequence[float]) – Strictly increasing x coordinates.

  • y (Sequence[float]) – Data values.

  • error_scale (float) – Multiplicative lambda scale.

  • weights (Sequence[float] | None) – Optional sample weights.

Returns:

Smoothed cubic spline representation.

Return type:

PiecewisePolynomial

interpolation.bspline_derivative(rep, x, k=1, eps=1e-06)

Evaluate derivative of B-spline via finite differences.

imsl_name: BSDER

Parameters:
  • rep (BSplineRepresentation) – B-spline representation.

  • x (float) – Query coordinate.

  • k (int) – Derivative order.

  • eps (float) – Finite-difference step.

Returns:

Derivative value.

Return type:

float

interpolation.bspline_easy(x, y, points, order=4)

Interpolate and evaluate B-spline in one call.

imsl_name: SPLEZ

Parameters:
  • x (Sequence[float]) – Sample coordinates.

  • y (Sequence[float]) – Sample values.

  • points (Sequence[float]) – Query points.

  • order (int) – Spline order.

Returns:

Evaluated values at query points.

Return type:

np.ndarray

interpolation.bspline_evaluate(rep, x)

Evaluate a univariate B-spline.

imsl_name: BSVAL

Parameters:
Returns:

Spline value.

Return type:

float

interpolation.bspline_evaluate_grid(rep, points)

Evaluate a univariate B-spline on a vector of points.

imsl_name: BS1GD

Parameters:
  • rep (BSplineRepresentation) – B-spline representation.

  • points (Sequence[float]) – Query coordinates.

Returns:

Evaluated values.

Return type:

np.ndarray

interpolation.bspline_integrate(rep, a, b, n=512)

Integrate a univariate B-spline over [a, b] using trapezoidal rule.

imsl_name: BSITG

Parameters:
  • rep (BSplineRepresentation) – B-spline representation.

  • a (float) – Lower bound.

  • b (float) – Upper bound.

  • n (int) – Number of quadrature points.

Returns:

Definite integral value.

Return type:

float

interpolation.bspline_interpolate(x, y, order=4)

Interpolate 1D data with a B-spline representation.

imsl_name: BSINT

Parameters:
  • x (Sequence[float]) – Sample coordinates.

  • y (Sequence[float]) – Sample values.

  • order (int) – Spline order.

Returns:

Interpolating B-spline representation.

Return type:

BSplineRepresentation

interpolation.bspline_knots(x, order=4)

Generate an open knot sequence for interpolation.

imsl_name: BSNAK

Parameters:
  • x (Sequence[float]) – Sample coordinates.

  • order (int) – Spline order.

Returns:

Knot vector.

Return type:

np.ndarray

interpolation.bspline_optimal_knots(x, order=4, n_coeff=None)

Generate a heuristic near-optimal knot vector via quantiles.

imsl_name: BSOPK

Parameters:
  • x (Sequence[float]) – Sample coordinates.

  • order (int) – Spline order.

  • n_coeff (int | None) – Number of coefficients.

Returns:

Knot vector.

Return type:

np.ndarray

interpolation.bspline_to_piecewise_poly(rep, n_segments=None)

Convert a B-spline representation to cubic piecewise polynomial form.

imsl_name: BSCPP

Parameters:
  • rep (BSplineRepresentation) – B-spline representation.

  • n_segments (int | None) – Number of output segments. Defaults to unique knot spans.

Returns:

Approximate cubic piecewise polynomial.

Return type:

PiecewisePolynomial

interpolation.bspline_2d_derivative(rep, x, y, dx_order=1, dy_order=0, eps=1e-06)

Evaluate partial derivative of a 2D B-spline via finite differences.

imsl_name: BS2DR

Parameters:
  • rep (dict) – Tensor B-spline representation.

  • x (float) – Query x.

  • y (float) – Query y.

  • dx_order (int) – Derivative order in x.

  • dy_order (int) – Derivative order in y.

  • eps (float) – Finite-difference step.

Returns:

Partial derivative value.

Return type:

float

interpolation.bspline_2d_evaluate(rep, x, y)

Evaluate a 2D tensor-product B-spline.

imsl_name: BS2VL

Parameters:
  • rep (dict) – Tensor B-spline representation.

  • x (float) – Query x.

  • y (float) – Query y.

Returns:

Spline value.

Return type:

float

interpolation.bspline_2d_evaluate_grid(rep, xq, yq)

Evaluate a 2D B-spline on a tensor query grid.

imsl_name: BS2GD

Parameters:
  • rep (dict) – Tensor B-spline representation.

  • xq (Sequence[float]) – Query x points.

  • yq (Sequence[float]) – Query y points.

Returns:

Evaluated grid values with shape (len(xq), len(yq)).

Return type:

np.ndarray

interpolation.bspline_2d_integrate(rep, ax, bx, ay, by, nx=128, ny=128)

Integrate 2D B-spline over rectangle [ax, bx] x [ay, by].

imsl_name: BS2IG

Parameters:
  • rep (dict) – Tensor B-spline representation.

  • ax (float) – x-lower bound.

  • bx (float) – x-upper bound.

  • ay (float) – y-lower bound.

  • by (float) – y-upper bound.

  • nx (int) – Number of x quadrature points.

  • ny (int) – Number of y quadrature points.

Returns:

Double integral value.

Return type:

float

interpolation.bspline_2d_interpolate(x, y, z, order_x=4, order_y=4)

Build 2D tensor-product B-spline interpolation.

imsl_name: BS2IN

Parameters:
  • x (Sequence[float]) – x coordinates.

  • y (Sequence[float]) – y coordinates.

  • z (np.ndarray) – Values of shape (len(x), len(y)).

  • order_x (int) – Spline order in x.

  • order_y (int) – Spline order in y.

Returns:

Tensor B-spline representation.

Return type:

dict

interpolation.bspline_3d_derivative(rep, x, y, z, dx_order=1, dy_order=0, dz_order=0, eps=0.0001)

Evaluate partial derivative of a 3D B-spline via finite differences.

imsl_name: BS3DR

For a mixed derivative of total order n (e.g. dx_order + dy_order + dz_order = 3), nested central differences divide by (2*eps)^n. With the default eps=1e-4 the round-off error for a 3rd-order mixed derivative is O(eps_machine / eps^3) ~ 1e-3, which is acceptable. Using eps=1e-6 would give ~1e3 (catastrophic cancellation).

Parameters:
  • rep (dict) – 3D tensor B-spline representation.

  • x (float) – Query x.

  • y (float) – Query y.

  • z (float) – Query z.

  • dx_order (int) – Derivative order in x.

  • dy_order (int) – Derivative order in y.

  • dz_order (int) – Derivative order in z.

  • eps (float) – Finite-difference step size. Defaults to 1e-4, which balances round-off and truncation for mixed derivatives up to order 3.

Returns:

Partial derivative value.

Return type:

float

interpolation.bspline_3d_evaluate(rep, x, y, z)

Evaluate a 3D tensor-product B-spline.

imsl_name: BS3VL

Parameters:
  • rep (dict) – 3D tensor B-spline representation.

  • x (float) – Query x.

  • y (float) – Query y.

  • z (float) – Query z.

Returns:

Spline value.

Return type:

float

interpolation.bspline_3d_evaluate_grid(rep, xq, yq, zq)

Evaluate a 3D B-spline on a tensor query grid.

imsl_name: BS3GD

Parameters:
  • rep (dict) – 3D tensor B-spline representation.

  • xq (Sequence[float]) – Query x values.

  • yq (Sequence[float]) – Query y values.

  • zq (Sequence[float]) – Query z values.

Returns:

Evaluated values with shape (len(xq), len(yq), len(zq)).

Return type:

np.ndarray

interpolation.bspline_3d_integrate(rep, ax, bx, ay, by, az, bz, nx=48, ny=48, nz=48)

Integrate 3D B-spline over box [ax, bx] x [ay, by] x [az, bz].

imsl_name: BS3IG

Parameters:
  • rep (dict) – 3D tensor B-spline representation.

  • ax (float) – x-lower bound.

  • bx (float) – x-upper bound.

  • ay (float) – y-lower bound.

  • by (float) – y-upper bound.

  • az (float) – z-lower bound.

  • bz (float) – z-upper bound.

  • nx (int) – x quadrature count.

  • ny (int) – y quadrature count.

  • nz (int) – z quadrature count.

Returns:

Triple integral value.

Return type:

float

interpolation.bspline_3d_interpolate(x, y, z, values, order_x=4, order_y=4, order_z=4)

Build 3D tensor-product B-spline interpolation.

imsl_name: BS3IN

Parameters:
  • x (Sequence[float]) – x coordinates.

  • y (Sequence[float]) – y coordinates.

  • z (Sequence[float]) – z coordinates.

  • values (np.ndarray) – Values of shape (len(x), len(y), len(z)).

  • order_x (int) – Spline order in x.

  • order_y (int) – Spline order in y.

  • order_z (int) – Spline order in z.

Returns:

Tensor B-spline representation.

Return type:

dict

interpolation.quad_grid_derivative_1d(values, x0, dx, xq)

Evaluate 1D derivative of quadratic interpolation on a uniform grid.

imsl_name: QDDER

Parameters:
  • values (Sequence[float]) – Sample values on uniform grid.

  • x0 (float) – Grid origin.

  • dx (float) – Grid spacing.

  • xq (float) – Query x coordinate.

Returns:

Derivative value.

Return type:

float

interpolation.quad_grid_derivative_2d(values, x0, dx, y0, dy, xq, yq, dx_order=1, dy_order=0, eps=1e-06)

Evaluate 2D partial derivative for quadratic gridded interpolation.

imsl_name: QD2DR

Parameters:
  • values (np.ndarray) – Grid values.

  • x0 (float) – x-origin.

  • dx (float) – x-spacing.

  • y0 (float) – y-origin.

  • dy (float) – y-spacing.

  • xq (float) – Query x coordinate.

  • yq (float) – Query y coordinate.

  • dx_order (int) – x derivative order.

  • dy_order (int) – y derivative order.

  • eps (float) – Finite-difference step.

Returns:

Partial derivative value.

Return type:

float

interpolation.quad_grid_derivative_3d(values, x0, dx, y0, dy, z0, dz, xq, yq, zq, dx_order=1, dy_order=0, dz_order=0, eps=1e-06)

Evaluate partial derivative for 3D quadratic gridded interpolation.

imsl_name: QD3DR

Parameters:
  • values (np.ndarray) – Grid values.

  • x0 (float) – x-origin.

  • dx (float) – x-spacing.

  • y0 (float) – y-origin.

  • dy (float) – y-spacing.

  • z0 (float) – z-origin.

  • dz (float) – z-spacing.

  • xq (float) – Query x coordinate.

  • yq (float) – Query y coordinate.

  • zq (float) – Query z coordinate.

  • dx_order (int) – x derivative order.

  • dy_order (int) – y derivative order.

  • dz_order (int) – z derivative order.

  • eps (float) – Finite-difference step.

Returns:

Partial derivative value.

Return type:

float

interpolation.quad_grid_evaluate_1d(values, x0, dx, xq)

Evaluate 1D quadratic interpolation on a uniform grid.

imsl_name: QDVAL

Parameters:
  • values (Sequence[float]) – Sample values on uniform grid.

  • x0 (float) – Grid origin.

  • dx (float) – Grid spacing.

  • xq (float) – Query x coordinate.

Returns:

Interpolated value.

Return type:

float

interpolation.quad_grid_evaluate_2d(values, x0, dx, y0, dy, xq, yq)

Evaluate 2D quadratic interpolation on a uniform tensor grid.

imsl_name: QD2VL

Parameters:
  • values (np.ndarray) – Grid values of shape (nx, ny).

  • x0 (float) – x-origin.

  • dx (float) – x-spacing.

  • y0 (float) – y-origin.

  • dy (float) – y-spacing.

  • xq (float) – Query x coordinate.

  • yq (float) – Query y coordinate.

Returns:

Interpolated value.

Return type:

float

interpolation.quad_grid_evaluate_3d(values, x0, dx, y0, dy, z0, dz, xq, yq, zq)

Evaluate 3D quadratic interpolation on a uniform tensor grid.

imsl_name: QD3VL

Parameters:
  • values (np.ndarray) – Grid values of shape (nx, ny, nz).

  • x0 (float) – x-origin.

  • dx (float) – x-spacing.

  • y0 (float) – y-origin.

  • dy (float) – y-spacing.

  • z0 (float) – z-origin.

  • dz (float) – z-spacing.

  • xq (float) – Query x coordinate.

  • yq (float) – Query y coordinate.

  • zq (float) – Query z coordinate.

Returns:

Interpolated value.

Return type:

float

interpolation.least_squares_bspline_2d(x, y, z, order_x=4, order_y=4)

Fit 2D tensor-product B-spline least squares surface.

imsl_name: BSLS2

Parameters:
  • x (Sequence[float]) – x coordinates.

  • y (Sequence[float]) – y coordinates.

  • z (np.ndarray) – Grid values.

  • order_x (int) – x spline order.

  • order_y (int) – y spline order.

Returns:

2D tensor B-spline representation.

Return type:

dict

interpolation.least_squares_bspline_3d(x, y, z, values, order_x=4, order_y=4, order_z=4)

Fit 3D tensor-product B-spline least-squares volume.

imsl_name: BSLS3

Parameters:
  • x (Sequence[float]) – x coordinates.

  • y (Sequence[float]) – y coordinates.

  • z (Sequence[float]) – z coordinates.

  • values (np.ndarray) – 3D grid values.

  • order_x (int) – x spline order.

  • order_y (int) – y spline order.

  • order_z (int) – z spline order.

Returns:

3D tensor B-spline representation.

Return type:

dict

interpolation.least_squares_bspline_constrained(x, y, knots, order, constraints_matrix, constraints_rhs, penalty=1000000.0)

Fit constrained B-spline least squares via quadratic penalty method.

imsl_name: CONFT

Parameters:
  • x (Sequence[float]) – Sample x coordinates.

  • y (Sequence[float]) – Sample y values.

  • knots (Sequence[float]) – Knot vector.

  • order (int) – Spline order.

  • constraints_matrix (np.ndarray) – Linear constraints matrix C.

  • constraints_rhs (Sequence[float]) – Constraint target vector d.

  • penalty (float) – Constraint penalty weight.

Returns:

Fitted constrained B-spline representation.

Return type:

BSplineRepresentation

interpolation.least_squares_bspline_fixed_knots(x, y, knots, order)

Fit B-spline least-squares coefficients with fixed knots.

imsl_name: BSLSQ

Parameters:
  • x (Sequence[float]) – Sample x coordinates.

  • y (Sequence[float]) – Sample y values.

  • knots (Sequence[float]) – Knot vector.

  • order (int) – Spline order.

Returns:

Fitted B-spline representation.

Return type:

BSplineRepresentation

interpolation.least_squares_bspline_variable_knots(x, y, order=4, n_coeff=None)

Fit B-spline least squares using quantile-based variable knots.

imsl_name: BSVLS

Parameters:
  • x (Sequence[float]) – Sample x coordinates.

  • y (Sequence[float]) – Sample y values.

  • order (int) – Spline order.

  • n_coeff (int | None) – Number of coefficients.

Returns:

Fitted B-spline representation.

Return type:

BSplineRepresentation

interpolation.least_squares_general(x, y, basis_functions)

Fit linear combination of arbitrary basis functions.

imsl_name: FNLSQ

Parameters:
  • x (Sequence[float]) – Sample x values.

  • y (Sequence[float]) – Sample y values.

  • basis_functions (Sequence[Callable]) – Basis functions evaluated at x.

Returns:

Fit result with coefficients and residual norm.

Return type:

dict

interpolation.least_squares_linear(x, y)

Fit linear model y = a + b*x by least squares.

imsl_name: RLINE

Parameters:
  • x (Sequence[float]) – Sample x values.

  • y (Sequence[float]) – Sample y values.

Returns:

Fit result with keys intercept, slope, residual_norm.

Return type:

dict

interpolation.least_squares_polynomial(x, y, degree)

Fit polynomial model by least squares.

imsl_name: RCURV

Parameters:
  • x (Sequence[float]) – Sample x values.

  • y (Sequence[float]) – Sample y values.

  • degree (int) – Polynomial degree.

Returns:

Fit result with descending polynomial coefficients and residual norm.

Return type:

dict

interpolation.akima_surface_fit(points, values, query_points, power=2.0)

Interpolate scattered 2D data using inverse-distance weighting.

imsl_name: SURF

Parameters:
  • points (np.ndarray) – Input sample points of shape (n, 2).

  • values (Sequence[float]) – Sample values of shape (n,).

  • query_points (np.ndarray) – Query points of shape (m, 2).

  • power (float) – Distance weighting exponent.

Returns:

Interpolated values of shape (m,).

Return type:

np.ndarray

Raises:

ValueError – If shapes are invalid.

interpolation.multidim_interpolate(points, values, query_points, derivative_axis=None, eps=1e-05)

Interpolate scattered N-D data and optionally estimate derivative.

imsl_name: SURFND

Parameters:
  • points (np.ndarray) – Sample points of shape (n, d).

  • values (Sequence[float]) – Sample values of shape (n,).

  • query_points (np.ndarray) – Query points of shape (m, d).

  • derivative_axis (int | None) – Axis along which to return derivative. If None, return interpolated values.

  • eps (float) – Finite-difference step for derivative.

Returns:

Interpolated values or derivative estimates.

Return type:

np.ndarray

interpolation.rational_chebyshev(x, y, numerator_degree, denominator_degree)

Compute rational approximation via linearized least-squares fit.

imsl_name: RATCH

Parameters:
  • x (Sequence[float]) – Sample x values.

  • y (Sequence[float]) – Sample y values.

  • numerator_degree (int) – Numerator polynomial degree.

  • denominator_degree (int) – Denominator polynomial degree.

Returns:

Rational fit with keys numerator, denominator,

max_abs_error.

Return type:

dict

Raises:

ValueError – If polynomial degrees are negative.

interpolation.spline_fitting(x, y, weights=None, order=4, n_coeff=None)

Fit a weighted least-squares univariate B-spline.

imsl_name: SPLINE_FITTING

Parameters:
  • x (Sequence[float]) – Sample x coordinates.

  • y (Sequence[float]) – Sample y values.

  • weights (Sequence[float] | None) – Optional non-negative weights.

  • order (int) – Spline order (degree + 1).

  • n_coeff (int | None) – Number of spline coefficients. Defaults to len(x).

Returns:

Fitted spline representation.

Return type:

BSplineRepresentation

Raises:

ValueError – If inputs are invalid.

interpolation.spline_values(rep, points)

Evaluate fitted spline values at query points.

imsl_name: SPLINE_VALUES

Parameters:
  • rep (BSplineRepresentation) – B-spline representation.

  • points (Sequence[float]) – Query coordinates.

Returns:

Evaluated values.

Return type:

np.ndarray

interpolation.surface_fitting(x, y, z, order_x=4, order_y=4, n_coeff_x=None, n_coeff_y=None)

Fit a tensor-product B-spline surface by least squares.

imsl_name: SURFACE_FITTING

Parameters:
  • x (Sequence[float]) – x grid coordinates of shape (nx,).

  • y (Sequence[float]) – y grid coordinates of shape (ny,).

  • z (np.ndarray) – Grid values of shape (nx, ny).

  • order_x (int) – Spline order in x direction.

  • order_y (int) – Spline order in y direction.

  • n_coeff_x (int | None) – Number of x-direction coefficients.

  • n_coeff_y (int | None) – Number of y-direction coefficients.

Returns:

Surface representation with keys knots_x, knots_y, coefficients, order_x, order_y.

Return type:

dict

Raises:

ValueError – If dimensions are invalid.

interpolation.surface_values(surface_rep, xq, yq)

Evaluate a tensor-product B-spline surface on query grid.

imsl_name: SURFACE_VALUES

Parameters:
  • surface_rep (dict) – Surface representation from surface_fitting.

  • xq (Sequence[float]) – Query x coordinates.

  • yq (Sequence[float]) – Query y coordinates.

Returns:

Values on query grid with shape (len(xq), len(yq)).

Return type:

np.ndarray

interpolation.launch_documentation(opener=<function open>, package_dir=None)

Launch package documentation in a web browser.

The launcher opens local HTML documentation only.

Parameters:
  • opener (Callable[[str], bool]) – Browser opener function. Defaults to webbrowser.open.

  • package_dir (Path | None) – Optional package directory override used for tests. Defaults to the directory of this module.

Returns:

The file://... URI opened.

Return type:

str

Raises:

FileNotFoundError – If no local documentation index can be located.