Extrapolation

Extrapolation controls behavior when query points fall outside the data domain [x[1], x[end]].

Overview

All extrapolation modes are concrete subtypes of AbstractExtrap. Pass them via the extrap keyword argument:

# Using Extrap() factory (recommended)
cubic_interp(x, y, xq; extrap=Extrap(:constant))
linear_interp(x, y, xq; extrap=Extrap(:extend))

# Using direct types (also supported)
cubic_interp(x, y, xq; extrap=ConstExtrap())
linear_interp(x, y, xq; extrap=ExtendExtrap())

# Interpolant: extrap is fixed at creation
itp = cubic_interp(x, y; extrap=Extrap(:extend))
itp(xq)  # uses ExtendExtrap
Factory Functions

The Extrap() factory provides a single discoverable entry point. For ND per-axis configuration, use the multi-arg form: Extrap(:extend, :constant, :none). See Factory Functions for details.

All interpolation methods (linear_interp, quadratic_interp, cubic_interp, constant_interp) support the same extrapolation modes.

TypeBehavior
NoExtrap()Throws DomainError (default)
ConstExtrap()Returns boundary values
ExtendExtrap()Extends boundary polynomial
WrapExtrap()Wraps coordinates periodically (no smoothness enforced)
AbstractExtrap
├── NoExtrap       # DomainError on out-of-domain (default)
├── ConstExtrap    # clamp to y[1] / y[end]
├── ExtendExtrap   # continue boundary polynomial
└── WrapExtrap     # modular coordinate mapping

Examples

using FastInterpolations

# Sample data
x = [0.0, 0.7, 1.5, 2.3, 3.0, 4.2, 5.0, 6.0]
y = [0.2, 1.1, 0.6, 1.8, 1.2, 0.4, 1.5, 0.8]

# Query points (including extrapolation region)
xq = range(x[1] - 1.5, x[end] + 1.5, 300)

NoExtrap() (Default)

Throws DomainError for out-of-domain queries. Use when extrapolation is unexpected.

julia> cubic_interp(x, y, -1.0; extrap=NoExtrap())  # scalar query outside domain
ERROR: DomainError with -1.0:
query point outside interpolation domain [0.0, 6.0]

julia> cubic_interp(x, y, xq; extrap=NoExtrap())  # vector query (xq includes out-of-domain points)
ERROR: DomainError with -1.5:
query point outside interpolation domain [0.0, 6.0]

Only interior queries succeed:

yq = cubic_interp(x, y, range(x[1], x[end], 200); extrap=NoExtrap())
Example block output

ConstExtrap()

Returns boundary values: y[1] for left, y[end] for right.

yq = cubic_interp(x, y, xq; extrap=ConstExtrap())
Example block output

ExtendExtrap()

Extends the boundary polynomial beyond the domain.

yq = cubic_interp(x, y, xq; extrap=ExtendExtrap())
Example block output

WrapExtrap()

Wraps coordinates periodically:

\[S(x + \tau) = S(x), \quad \tau = x_{\text{end}} - x_1\]

This is purely coordinate mapping—it does not enforce any physical conditions at the boundary. The spline may have discontinuities in value, slope, or curvature at the wrap point.

For Smooth Periodicity

If you need C² continuity at the periodic boundary, use bc=PeriodicBC() with cubic_interp. This enforces $S(x_1) = S(x_{\text{end}})$, $S'(x_1) = S'(x_{\text{end}})$, and $S''(x_1) = S''(x_{\text{end}})$.

yq = cubic_interp(x, y, xq; extrap=WrapExtrap())
Example block output

Comparison

# All three modes on same plot
y_const = cubic_interp(x, y, xq; extrap=ConstExtrap())
y_ext   = cubic_interp(x, y, xq; extrap=ExtendExtrap())
y_wrap  = cubic_interp(x, y, xq; extrap=WrapExtrap())
Example block output

Summary

ModeBehaviorUse Case
NoExtrap()DomainErrorStrict domain enforcement (default)
ConstExtrap()Returns boundary valuesPhysical constraints
ExtendExtrap()Continues boundary polynomialSmooth continuation
WrapExtrap()Wraps coordinates (no smoothness)Cyclic data (see PeriodicBC for C² continuity)

See Also