Local Cubic Hermite API

Four C^1 interpolation methods sharing the same Hermite basis kernel, differing only in how slopes are determined.

Overview

Hermite (user-supplied slopes)

FunctionDescription
cubic_interp(x, Hermite(y, dy), xq)Hermite interpolation at point(s) xq
cubic_interp!(out, x, Hermite(y, dy), xq)In-place Hermite interpolation
itp = cubic_interp(x, Hermite(y, dy))Create callable interpolant

PCHIP (monotone-preserving)

FunctionDescription
pchip_interp(x, y, xq)PCHIP interpolation at point(s) xq
pchip_interp!(out, x, y, xq)In-place PCHIP interpolation
itp = pchip_interp(x, y)Create callable interpolant

Cardinal / Catmull-Rom

FunctionDescription
cardinal_interp(x, y, xq; tension=0.0)Cardinal spline at point(s) xq
cardinal_interp!(out, x, y, xq; tension=0.0)In-place cardinal interpolation
itp = cardinal_interp(x, y; tension=0.0)Create callable interpolant

Akima (outlier-robust)

FunctionDescription
akima_interp(x, y, xq)Akima interpolation at point(s) xq
akima_interp!(out, x, y, xq)In-place Akima interpolation
itp = akima_interp(x, y)Create callable interpolant

Data Wrapper

FastInterpolations.HermiteType
Hermite(y, dy)

Wrap user-supplied function values and first derivatives for Hermite interpolation. Use with cubic_interp to skip the global spline solve and use the provided slopes directly.

1D Usage

y  = sin.(x)
dy = cos.(x)   # exact first derivative
val = cubic_interp(x, Hermite(y, dy), 0.5)
itp = cubic_interp(x, Hermite(y, dy))
itp(0.5)

Notes

  • No boundary conditions (bc) — user slopes replace them.
  • No auto-caching (autocache) — no global solve to cache.
  • C$^1$ continuity (not C$^2$ like cubic spline).
  • deriv, extrap, and search kwargs are supported as usual.
source

Functions — Hermite

Hermite interpolation uses cubic_interp with a Hermite(y, dy) wrapper. See the Cubic API for the cubic_interp / cubic_interp! reference.

Functions — PCHIP

FastInterpolations.pchip_interpFunction
pchip_interp(x, y, xq; extrap=NoExtrap(), deriv=EvalValue(), search=AutoSearch(), hint=nothing)

PCHIP (monotone-preserving) interpolation at a single query point.

Computes Fritsch-Carlson slopes internally, then evaluates via cubic Hermite kernel. C$^1$ continuous, monotonicity guaranteed for monotone input data.

source
pchip_interp(x, y, x_query; extrap=NoExtrap(), deriv=EvalValue(), search=AutoSearch(), hint=nothing)

PCHIP interpolation at multiple query points. Returns Vector{Tv}.

source
pchip_interp(x, y; extrap=NoExtrap(), search=AutoSearch()) -> PchipInterpolant1D

Create a callable PCHIP interpolant with monotone-preserving slopes.

Slopes are computed once at construction via the Fritsch-Carlson algorithm. Subsequent itp(xq) calls just evaluate — zero slope recomputation.

Arguments

  • x::AbstractVector: x-coordinates (must be sorted, ≥2 points)
  • y::AbstractVector: y-values
  • extrap::AbstractExtrap: Extrapolation mode (default: NoExtrap())
  • search::AbstractSearchPolicy: Default search policy (default: AutoSearch())

Returns

PchipInterpolant1D callable object. Supports:

  • Scalar: itp(0.5)
  • Vector: itp(xq_vec)
  • In-place: itp(out, xq_vec)
  • Broadcast: itp.(xq) or @. coef * itp(rho)
  • Derivative: itp(0.5; deriv=DerivOp(1))

Example

x = [0.0, 1.0, 2.0, 3.0, 4.0]
y = [0.0, 1.0, 0.5, 0.8, 0.2]  # non-monotone data
itp = pchip_interp(x, y)
itp(1.5)                          # monotone within each monotone segment
itp(1.5; deriv=DerivOp(1))       # first derivative
source
FastInterpolations.pchip_interp!Function
pchip_interp!(output, x, y, x_query; extrap=NoExtrap(), deriv=EvalValue(), search=AutoSearch(), hint=nothing)

In-place PCHIP interpolation with monotone-preserving slopes.

source

Functions — Cardinal

FastInterpolations.cardinal_interpFunction
cardinal_interp(x, y, xq; tension=0.0, extrap=NoExtrap(), deriv=EvalValue(), search=AutoSearch(), hint=nothing)

Cardinal spline interpolation at a single query point. Default tension=0 is Catmull-Rom. C$^1$ continuous.

source
cardinal_interp(x, y, x_query; tension=0.0, extrap=NoExtrap(), deriv=EvalValue(), search=AutoSearch(), hint=nothing)

Cardinal spline interpolation at multiple query points. Returns Vector{Tv}.

source
cardinal_interp(x, y; tension=0.0, extrap=NoExtrap(), search=AutoSearch()) -> CardinalInterpolant1D

Create a callable cardinal spline interpolant.

Default tension=0 gives Catmull-Rom (central finite difference slopes). Increasing tension reduces overshoot; tension=1 gives zero slopes at knots (smooth S-curves between knots).

Example

itp = cardinal_interp(x, y)                # CatmullRom (tension=0)
itp = cardinal_interp(x, y; tension=0.5)   # tighter curve
itp(0.5)
source
FastInterpolations.cardinal_interp!Function
cardinal_interp!(output, x, y, x_query; tension=0.0, extrap=NoExtrap(), deriv=EvalValue(), search=AutoSearch(), hint=nothing)

In-place cardinal spline interpolation.

source

Functions — Akima

FastInterpolations.akima_interpFunction
akima_interp(x, y, xq; extrap=NoExtrap(), deriv=EvalValue(), search=AutoSearch(), hint=nothing)

Akima interpolation at a single query point. Outlier-robust, C$^1$ continuous.

source
akima_interp(x, y, x_query; extrap=NoExtrap(), deriv=EvalValue(), search=AutoSearch(), hint=nothing)

Akima interpolation at multiple query points. Returns Vector{Tv}.

source
akima_interp(x, y; extrap=NoExtrap(), search=AutoSearch()) -> AkimaInterpolant1D

Create a callable Akima interpolant with outlier-robust slopes.

Slopes are computed once at construction. Requires ≥ 2 grid points.

Example

itp = akima_interp(x, y)
itp(0.5)
itp(0.5; deriv=DerivOp(1))
source

Interpolant Types

FastInterpolations.CubicHermiteInterpolant1DType
CubicHermiteInterpolant1D{Tg, Tv, X, Y, DY, S, E, P}

Callable interpolant for cubic Hermite interpolation with user-supplied slopes. Returned by cubic_interp(x, Hermite(y, dy)) (2-argument form).

Stores precomputed grid spacing for O(1) h/inv_h lookup on uniform grids. Evaluation uses _hermite_kernel_1d (derivative-based Hermite basis functions), NOT _cubic_kernel (moment-based spline formulation).

Type Parameters

  • Tg<:AbstractFloat: Grid coordinate type
  • Tv: Value type (unconstrained — supports Complex, etc.)
  • X<:AbstractVector{Tg}: Grid vector type (preserves Range for O(1) lookup)
  • Y<:AbstractVector{Tv}: Values vector type
  • DY<:AbstractVector{Tv}: Slopes vector type
  • S<:AbstractGridSpacing{Tg}: Grid spacing type
  • E<:AbstractExtrap: Extrapolation mode type
  • P<:AbstractSearchPolicy: Search policy type

Usage

itp = cubic_interp(x, Hermite(y, dy))
itp(0.5)                          # scalar query
itp(xq_vec)                       # vector query (allocating)
itp(out, xq_vec)                  # vector query (in-place)
itp(0.5; deriv=DerivOp(1))       # first derivative
itp(0.5; search=BinarySearch())  # override search
source
FastInterpolations.PchipInterpolant1DType
PchipInterpolant1D{Tg, Tv, X, Y, DY, S, E, P}

Callable interpolant for PCHIP (Piecewise Cubic Hermite Interpolating Polynomial). Returned by pchip_interp(x, y) (2-argument form).

Slopes are computed once at construction via the Fritsch-Carlson monotone-preserving algorithm. Evaluation uses the same cubic Hermite kernel as CubicHermiteInterpolant1D.

Properties

  • C$^1$ continuous (continuous first derivative, discontinuous second derivative at knots)
  • Monotonicity preserving: monotone input data → monotone interpolant
  • Local: each slope depends only on neighboring data — no global solve

Usage

itp = pchip_interp(x, y)
itp(0.5)                          # scalar query
itp(xq_vec)                       # vector query (allocating)
itp(out, xq_vec)                  # vector query (in-place)
itp(0.5; deriv=DerivOp(1))       # first derivative
source
FastInterpolations.CardinalInterpolant1DType
CardinalInterpolant1D{Tg, Tv, X, Y, DY, S, E, P}

Callable interpolant for cardinal spline interpolation. Returned by cardinal_interp(x, y) (2-argument form).

Slopes are computed once at construction via: d_k = (1 - tension) * (y[k+1] - y[k-1]) / (x[k+1] - x[k-1])

When tension=0, this is the Catmull-Rom spline (simple central FD).

Properties

  • C$^1$ continuous (continuous first derivative)
  • Local: each slope depends only on immediate neighbors
  • Tension parameter: controls overshoot (0 = CatmullRom, 1 = zero slopes at knots)

Usage

itp = cardinal_interp(x, y)                     # default: CatmullRom (tension=0)
itp = cardinal_interp(x, y; tension=0.5)        # tighter curve
itp(0.5)
itp(0.5; deriv=DerivOp(1))
source
FastInterpolations.AkimaInterpolant1DType
AkimaInterpolant1D{Tg, Tv, X, Y, DY, S, E, P}

Callable interpolant for Akima interpolation. Returned by akima_interp(x, y) (2-argument form).

Slopes are computed once at construction via the Akima (1970) weighted-average algorithm. Evaluation uses the same cubic Hermite kernel.

Properties

  • C$^1$ continuous (continuous first derivative)
  • Local: each slope depends on 4 adjacent secants (5-point stencil)
  • Outlier-robust: deviant secants receive less weight

Usage

itp = akima_interp(x, y)
itp(0.5)
itp(0.5; deriv=DerivOp(1))
source