Type Reference

Abstract Types

Interpolant Hierarchy

FastInterpolations.AbstractInterpolantType
AbstractInterpolant{T<:AbstractFloat}

Abstract supertype for all interpolant objects.

Type Parameter

  • T: Float type (Float32 or Float64)

Subtypes

  • AbstractSeriesInterpolant{T}: Multi-series interpolants
  • LinearInterpolant{T}: Piecewise linear interpolation
  • ConstantInterpolant{T}: Piecewise constant (step) interpolation
  • QuadraticInterpolant{T}: C1 piecewise quadratic spline
  • CubicInterpolant{T}: C2 natural/clamped/periodic cubic spline

Note

This is a pure type hierarchy - no methods are defined on AbstractInterpolant itself. All functionality is implemented in concrete subtypes.

source
FastInterpolations.AbstractSeriesInterpolantType
AbstractSeriesInterpolant{T<:AbstractFloat}

Abstract supertype for multi-series interpolant objects. Series interpolants handle multiple y-series sharing the same x-grid.

Type Parameter

  • T: Float type (Float32 or Float64)

Subtypes

  • LinearSeriesInterpolant{T}: Multiple linear interpolants sharing x-grid
  • ConstantSeriesInterpolant{T}: Multiple constant interpolants sharing x-grid
  • QuadraticSeriesInterpolant{T}: Multiple quadratic interpolants sharing x-grid
  • CubicSeriesInterpolant{T}: Multiple cubic splines sharing x-grid

Key Features

  • Anchor optimization: compute interval once, evaluate all series
  • Matrix storage: unified storage for optimal SIMD vectorization
  • Zero-allocation batch evaluation with pre-built anchors

Usage

x = collect(range(0.0, 1.0, 101))
y1, y2, y3 = sin.(2π .* x), cos.(2π .* x), exp.(-x)

sitp = cubic_interp(x, [y1, y2, y3])  # Creates CubicSeriesInterpolant

vals = sitp(0.5)            # Returns Vector of 3 values
sitp(output, 0.5)           # In-place evaluation

Note

This is a pure type hierarchy - no methods are defined on AbstractSeriesInterpolant itself. All functionality is implemented in concrete subtypes.

source

Evaluation Operations

FastInterpolations.AbstractEvalOpType
AbstractEvalOp

Abstract type for evaluation operations (value, derivatives). Used for compile-time dispatch to select appropriate kernel.

Subtypes:

  • EvalValue: Evaluate function value f(x)
  • EvalDeriv1: Evaluate first derivative f'(x)
  • EvalDeriv2: Evaluate second derivative f''(x)
  • EvalDeriv3: Evaluate third derivative f'''(x)
source
FastInterpolations.EvalDeriv3Type
EvalDeriv3 <: AbstractEvalOp

Singleton type indicating evaluation of third derivative f'''(x).

Note

For cubic splines, S'''(x) is constant within each interval. Linear/quadratic/constant interpolants always return zero.

Mathematical Background

The cubic spline third derivative is: S'''(x) = (zR - zL) / h

where zL and zR are the second derivative values (moments) at the interval endpoints, and h is the interval width.

source

Boundary Conditions

Cubic Splines (BCPair)

FastInterpolations.AbstractBCType
AbstractBC{T<:AbstractFloat}

Abstract base type for all boundary condition specifications.

Subtypes

  • NaturalBC{T}: Natural BC (zero curvature at both ends) - default
  • ClampedBC{T}: Clamped BC (zero slope at both ends)
  • PeriodicBC{T}: Periodic boundary condition
  • PointBC{T}: Single-point derivative conditions (Deriv1, Deriv2)
  • BCPair{T,L,R}: Pair of left/right boundary conditions
  • Left{T,B}: Endpoint wrapper for BC at left (x[1]) - used by quadratic splines
  • Right{T,B}: Endpoint wrapper for BC at right (x[end]) - used by quadratic splines
source
FastInterpolations.PointBCType
PointBC{T<:AbstractFloat} <: AbstractBC{T}

Abstract type for single-point boundary conditions. Represents a derivative condition at one endpoint.

Subtypes

  • Deriv1{T}: First derivative (slope) BC
  • Deriv2{T}: Second derivative (curvature) BC
source
FastInterpolations.Deriv1Type
Deriv1{T<:AbstractFloat} <: PointBC{T}

First derivative (slope) boundary condition: S'(endpoint) = val

Example

Deriv1(0.5)  # Slope of 0.5 at endpoint
Deriv1(0)    # Zero slope (horizontal tangent)
source
FastInterpolations.Deriv2Type
Deriv2{T<:AbstractFloat} <: PointBC{T}

Second derivative (curvature) boundary condition: S''(endpoint) = val

Example

Deriv2(0)    # Natural BC (zero curvature)
Deriv2(1.5)  # Specified curvature at endpoint
source
FastInterpolations.Deriv3Type
Deriv3{T<:AbstractFloat} <: PointBC{T}

Third derivative boundary condition: S'''(endpoint) = val

For cubic splines, the third derivative is constant within each interval: S'''(x) = (z[i+1] - z[i]) / h[i]. This BC specifies the third derivative value at the first (or last) interval.

Example

Deriv3(0)    # Zero third derivative at endpoint
Deriv3(1.0)  # Specified third derivative
source
FastInterpolations.BCPairType
BCPair{T, L<:PointBC{T}, R<:PointBC{T}} <: AbstractBC{T}

Container for left and right boundary conditions with type parameters for zero-overhead dispatch. The BC types are encoded in the type parameters, enabling compile-time specialization.

Type Parameters

  • T: Float type
  • L: Left boundary condition type (Deriv1{T} or Deriv2{T})
  • R: Right boundary condition type (Deriv1{T} or Deriv2{T})

Example

bc = BCPair(Deriv1(0.5), Deriv2(0))  # Left: slope=0.5, Right: natural
source
FastInterpolations.NaturalBCType
NaturalBC{T<:AbstractFloat} <: AbstractBC{T}

Natural boundary condition: S''(endpoints) = 0 (zero curvature at both ends). Equivalent to BCPair(Deriv2(0), Deriv2(0)).

This is the default BC for cubic spline interpolation.

Example

itp = cubic_interp(x, y; bc=NaturalBC())  # Default
itp = cubic_interp(x, y)                   # Same as above
source
FastInterpolations.ClampedBCType
ClampedBC{T<:AbstractFloat} <: AbstractBC{T}

Clamped boundary condition: S'(endpoints) = 0 (zero slope at both ends). Equivalent to BCPair(Deriv1(0), Deriv1(0)).

Also known as "complete" spline with zero derivative.

Example

itp = cubic_interp(x, y; bc=ClampedBC())
source
FastInterpolations.PeriodicBCType
PeriodicBC{T<:AbstractFloat} <: AbstractBC{T}

Periodic boundary condition: S(x0) = S(xn), S'(x0) = S'(xn), S''(x0) = S''(xn)

This is a user-facing type. Internally, periodic BC uses Sherman-Morrison solver with PeriodicData{T} for the cache.

Example

cache = CubicSplineCache(x; bc=PeriodicBC())
source

Quadratic Splines (Single Endpoint)

FastInterpolations.LeftType
Left{T, B<:PointBC{T}} <: AbstractBC{T}

Wrapper indicating BC is applied at left endpoint (x[1]). Used for quadratic splines where only one endpoint BC is specified.

Example

bc = Left(Deriv1(0.5))   # slope = 0.5 at left endpoint
bc = Left(Deriv2(0.0))   # curvature = 0 at left endpoint
source
FastInterpolations.RightType
Right{T, B<:PointBC{T}} <: AbstractBC{T}

Wrapper indicating BC is applied at right endpoint (x[end]). Used for quadratic splines where only one endpoint BC is specified.

Example

bc = Right(Deriv1(2.0))  # slope = 2.0 at right endpoint
bc = Right(Deriv2(0.0))  # curvature = 0 at right endpoint
source

Series Interpolant Types

SeriesInterpolant stores multiple y-series in a unified matrix with point-contiguous layout for optimal SIMD performance on scalar queries (10-120× faster than composition-based approaches).

Constant Interpolation

FastInterpolations.ConstantSeriesInterpolantType
ConstantSeriesInterpolant{T}

Multi-series constant (step) interpolant with unified matrix storage and SIMD optimization. Shares a single x-grid across N y-series for efficient batch evaluation.

Type Parameters

  • T: Float type (Float32 or Float64)

Fields

  • x::Vector{T}: Shared x-grid
  • y::Matrix{T}: Function values (npoints × nseries) series-contiguous
  • _transpose::LazyTranspose{T}: Lazy point-contiguous layout for scalar SIMD
  • extrap::ExtrapVal: Extrapolation mode
  • side::SideVal: Side selection (:nearest, :left, :right)

Memory Layout

Primary storage is series-contiguous (npoints × nseries):

  • y[i, k] = value of series k at grid point i
  • y[:, k] is contiguous → optimal for vector queries

Lazy transpose (nseries × npoints) for scalar queries:

  • y_point[:, i] is contiguous → optimal for SIMD scalar queries

Usage

x = collect(range(0.0, 1.0, 101))
y1, y2, y3 = sin.(2π .* x), cos.(2π .* x), exp.(-x)

sitp = constant_interp(x, [y1, y2, y3])

# Scalar evaluation
vals = sitp(0.5)            # Returns Vector{Float64} of length 3
sitp(output, 0.5)           # In-place

# Vector evaluation
vals = sitp([0.1, 0.5, 0.9])    # Returns Vector of Vectors
sitp([out1, out2, out3], xq)    # In-place (zero allocation)

Implementation Note: mutable struct with const fields

This type uses mutable struct with all const fields (Julia 1.8+) instead of plain struct for performance reasons. See CubicSeriesInterpolant for details.

source

Linear Interpolation

FastInterpolations.LinearSeriesInterpolantType
LinearSeriesInterpolant{T}

Multi-series linear interpolant with unified matrix storage and SIMD optimization. Shares a single x-grid across N y-series for efficient batch evaluation.

Type Parameters

  • T: Float type (Float32 or Float64)

Fields

  • x::Vector{T}: Shared x-grid
  • y::Matrix{T}: Function values (npoints × nseries) series-contiguous
  • _transpose::LazyTranspose{T}: Lazy point-contiguous layout for scalar SIMD
  • extrap::ExtrapVal: Extrapolation mode

Memory Layout

Primary storage is series-contiguous (npoints × nseries):

  • y[i, k] = value of series k at grid point i
  • y[:, k] is contiguous → optimal for vector queries

Lazy transpose (nseries × npoints) for scalar queries:

  • y_point[:, i] is contiguous → optimal for SIMD scalar queries

Usage

x = collect(range(0.0, 1.0, 101))
y1, y2, y3 = sin.(2π .* x), cos.(2π .* x), exp.(-x)

sitp = linear_interp(x, [y1, y2, y3])

# Scalar evaluation
vals = sitp(0.5)            # Returns Vector{Float64} of length 3
sitp(output, 0.5)           # In-place

# Vector evaluation
vals = sitp([0.1, 0.5, 0.9])    # Returns Vector of Vectors
sitp([out1, out2, out3], xq)    # In-place (zero allocation)

Performance

  • Vector queries use series-contiguous layout directly
  • Scalar queries trigger lazy transpose on first call

Implementation Note: mutable struct with const fields

This type uses mutable struct with all const fields (Julia 1.8+) instead of plain struct for performance reasons. See CubicSeriesInterpolant for details.

source

Quadratic Interpolation

FastInterpolations.QuadraticSeriesInterpolantType
QuadraticSeriesInterpolant{T}

Multi-series quadratic spline interpolant with unified matrix storage and SIMD optimization. Shares a single x-grid across N y-series for efficient batch evaluation.

Type Parameters

  • T: Float type (Float32 or Float64)

Fields

  • x::Vector{T}: Grid points (sorted)
  • y::Matrix{T}: Function values (npoints × nseries) series-contiguous
  • a::Matrix{T}: Quadratic coefficients (npoints × nseries) series-contiguous
  • d::Matrix{T}: Slope coefficients (npoints × nseries) series-contiguous
  • h::Vector{T}: Grid spacing (shared across all series)
  • _transpose::LazyTransposeTriple{T}: Lazy point-contiguous layout for SIMD
  • extrap::ExtrapVal: Extrapolation mode

Memory Layout

Primary storage is series-contiguous (npoints × nseries):

  • y[i, k] = value of series k at grid point i
  • y[:, k] is contiguous → optimal for vector queries

Lazy transpose (nseries × npoints) for scalar queries:

  • y_point[:, i] is contiguous → optimal for SIMD scalar queries

Usage

x = collect(range(0.0, 1.0, 101))
y1, y2, y3 = sin.(2π .* x), cos.(2π .* x), exp.(-x)

sitp = quadratic_interp(x, [y1, y2, y3])

# Scalar evaluation
vals = sitp(0.5)            # Returns Vector{Float64} of length 3
sitp(output, 0.5)           # In-place

# Vector evaluation
vals = sitp([0.1, 0.5, 0.9])    # Returns Vector of Vectors
sitp([out1, out2, out3], xq)    # In-place (zero allocation)

# Derivatives
d1 = sitp(0.5; deriv=1)     # First derivatives
d2 = sitp(0.5; deriv=2)     # Second derivatives

Performance

  • Vector queries use series-contiguous layout directly
  • Scalar queries trigger lazy transpose on first call
  • All series share same h[] array (O(1) memory overhead)

Implementation Note: mutable struct with const fields

This type uses mutable struct with all const fields (Julia 1.8+) instead of plain struct for performance reasons. See CubicSeriesInterpolant for details.

source

Cubic Interpolation

FastInterpolations.CubicSeriesInterpolantType
CubicSeriesInterpolant{T, C, B}

Multi-series cubic spline interpolant with unified matrix storage and SIMD optimization. Shares a single x-grid across N y-series for efficient batch evaluation.

Type Parameters

  • T: Float type (Float32 or Float64)
  • C: Cache type (CubicSplineCache{T})
  • B: Boundary condition config type (BCPair or PeriodicData)

Fields

  • cache::C: Shared CubicSplineCache with LU factorization
  • bc_for_solve::B: BC configuration for solving systems
  • y::Matrix{T}: Function values (npoints × nseries) series-contiguous
  • z::Matrix{T}: Second derivatives (npoints × nseries) series-contiguous
  • _point_snapshot: Atomic field for lazy point-contiguous layout
  • extrap::ExtrapVal: Extrapolation mode

Memory Layout

Primary storage is series-contiguous (npoints × nseries):

  • y[i, k] = value of series k at grid point i
  • y[:, k] is contiguous → optimal for vector queries

Lazy transpose (nseries × npoints) for scalar queries:

  • y_point[:, i] is contiguous → optimal for SIMD scalar queries

Usage

x = collect(range(0.0, 1.0, 101))
y1, y2, y3 = sin.(2π .* x), cos.(2π .* x), exp.(-x)

sitp = cubic_interp(x, [y1, y2, y3])

# Scalar evaluation
vals = sitp(0.5)            # Returns Vector{Float64} of length 3
sitp(output, 0.5)           # In-place

# Vector evaluation
vals = sitp([0.1, 0.5, 0.9])    # Returns Vector of Vectors
sitp([out1, out2, out3], xq)    # In-place (zero allocation)

Performance

  • Vector queries use series-contiguous layout directly
  • Scalar queries trigger lazy transpose on first call
  • All series share same cache (O(1) memory overhead)

Implementation Note: mutable struct with const fields

This type uses mutable struct with all const fields (Julia 1.8+) instead of plain struct for performance reasons. The const annotation ensures fields cannot be reassigned while allowing heap allocation. This pattern provides:

  • Stable memory addresses for the compiler to optimize field access
  • Better inlining of field accesses compared to plain immutable structs
  • Compatibility with mutable inner types (LazyTransposePair with @atomic)

Benchmarks show ~15% regression when using plain struct instead.

source
FastInterpolations.precompute_transpose!Function
precompute_transpose!(lt::LazyTranspose, src::Matrix) -> LazyTranspose

Pre-compute transpose before hot loops. Returns the holder for chaining.

source
precompute_transpose!(ltp::LazyTransposePair, y::Matrix, z::Matrix) -> LazyTransposePair

Pre-compute transpose pair before hot loops. Returns the holder for chaining.

source
precompute_transpose!(ltt::LazyTransposeTriple, y::Matrix, a::Matrix, d::Matrix) -> LazyTransposeTriple

Pre-compute transpose triple before hot loops. Returns the holder for chaining.

source
precompute_transpose!(sitp::LinearSeriesInterpolant) -> sitp

Pre-allocate point-contiguous matrix for scalar queries. Call before hot loops to avoid first-call latency.

source
precompute_transpose!(sitp::ConstantSeriesInterpolant) -> sitp

Pre-allocate point-contiguous matrix for scalar queries. Call before hot loops to avoid first-call latency.

source
precompute_transpose!(sitp::CubicSeriesInterpolant) -> sitp

Pre-allocate point-contiguous matrices for scalar queries. Call before hot loops to avoid first-call latency.

source

Index