Public API
Interpolations.CubicSplineInterpolation — Functionetp = CubicSplineInterpolation(knots, A; bc=Line(OnGrid()), extrapolation_bc=Throw())A shorthand for extrapolate(interpolate(knots, A, BSpline(Cubic(bc))), extrapolation_bc).
Interpolations.LinearInterpolation — Functionetp = LinearInterpolation(knots, A; extrapolation_bc=Throw())A shorthand for extrapolate(interpolate(knots, A, scheme), extrapolation_bc), where scheme is either BSpline(Linear()) or Gridded(Linear()) depending on whether knots are ranges or vectors.
Interpolations.bounds — Methodbounds(itp::AbstractInterpolation)Return the bounds of the domain of itp as a tuple of (min, max) pairs for each coordinate. This is best explained by example:
julia> itp = interpolate([1 2 3; 4 5 6], BSpline(Linear()));
julia> bounds(itp)
((1, 2), (1, 3))
julia> data = 1:3;
julia> knots = ([10, 11, 13.5],);
julia> itp = interpolate(knots, data, Gridded(Linear()));
julia> bounds(itp)
((10.0, 13.5),)Interpolations.extrapolate — Methodextrapolate(itp, scheme) adds extrapolation behavior to an interpolation object, according to the provided scheme.
The scheme can take any of these values:
Throw- throws a BoundsError for out-of-bounds indicesFlat- for constant extrapolation, taking the closest in-bounds valueLine- linear extrapolation (the wrapped interpolation object must support gradient)Reflect- reflecting extrapolation (indices must supportmod)Periodic- periodic extrapolation (indices must supportmod)
You can also combine schemes in tuples. For example, the scheme (Line(), Flat()) will use linear extrapolation in the first dimension, and constant in the second.
Finally, you can specify different extrapolation behavior in different direction. ((Line(),Flat()), Flat()) will extrapolate linearly in the first dimension if the index is too small, but use constant etrapolation if it is too large, and always use constant extrapolation in the second dimension.
Interpolations.extrapolate — Methodextrapolate(itp, fillvalue) creates an extrapolation object that returns the fillvalue any time the indexes in itp(x1,x2,...) are out-of-bounds.
Interpolations.interpolate — Methoditp = interpolate(A, interpmode)Interpolate an array A in the mode determined by interpmode. interpmode may be one of
NoInterp()BSpline(Constant())BSpline(Linear())BSpline(Quadratic(bc))(seeBoundaryCondition)BSpline(Cubic(bc))
It may also be a tuple of such values, if you want to use different interpolation schemes along each axis.
Interpolations.interpolate — Methoditp = interpolate((nodes1, nodes2, ...), A, interpmode)Interpolate an array A on a non-uniform but rectangular grid specified by the given nodes, in the mode determined by interpmode.
interpmode may be one of
NoInterp()Gridded(Constant())Gridded(Linear())
It may also be a tuple of such values, if you want to use different interpolation schemes along each axis.
Interpolations.scale — Methodscale(itp, xs, ys, ...) scales an existing interpolation object to allow for indexing using other coordinate axes than unit ranges, by wrapping the interpolation object and transforming the indices from the provided axes onto unit ranges upon indexing.
The parameters xs etc must be either ranges or linspaces, and there must be one coordinate range/linspace for each dimension of the interpolation object.
For every NoInterp dimension of the interpolation object, the range must be exactly 1:size(itp, d).
Interpolations.BSpline — TypeBSpline(degree)A flag signaling BSpline (integer-grid b-spline) interpolation along the corresponding axis. degree is one of Constant, Linear, Quadratic, or Cubic.
Interpolations.CardinalMonotonicInterpolation — TypeCardinalMonotonicInterpolation(tension)Cubic cardinal splines, uses tension parameter which must be between [0,1] Cubin cardinal splines can overshoot for non-monotonic data (increasing tension reduces overshoot).
Interpolations.Constant — TypeConstant b-splines are nearest-neighbor interpolations, and effectively return A[round(Int,x)] when interpolating.
Interpolations.Cubic — TypeCubic(bc::BoundaryCondition)Indicate that the corresponding axis should use quadratic interpolation.
Extended help
Assuming uniform knots with spacing 1, the ith piece of cubic spline implemented here is defined as follows.
y_i(x) = cm p(x-i) + c q(x-i) + cp q(1- (x-i)) + cpp p(1 - (x-i))where
p(δx) = 1/6 * (1-δx)^3
q(δx) = 2/3 - δx^2 + 1/2 δx^3and the values cX for X ∈ {m, _, p, pp} are the pre-filtered coefficients.
For future reference, this expands out to the following polynomial:
y_i(x) = 1/6 cm (1+i-x)^3 + c (2/3 - (x-i)^2 + 1/2 (x-i)^3) +
cp (2/3 - (1+i-x)^2 + 1/2 (1+i-x)^3) + 1/6 cpp (x-i)^3When we derive boundary conditions we will use derivatives y_0'(x) and y_0''(x)
Interpolations.FiniteDifferenceMonotonicInterpolation — TypeFiniteDifferenceMonotonicInterpolationClassic cubic interpolation, no tension parameter. Finite difference can overshoot for non-monotonic data.
Interpolations.Flat — TypeFlat(gt) sets the extrapolation slope to zero
Interpolations.FritschButlandMonotonicInterpolation — TypeFritschButlandMonotonicInterpolationMonotonic interpolation based on Fritsch & Butland (1984), "A Method for Constructing Local Monotone Piecewise Cubic Interpolants", doi:10.1137/0905021.
Faster than FritschCarlsonMonotonicInterpolation (only requires one pass) but somewhat higher apparent "tension".
Interpolations.FritschCarlsonMonotonicInterpolation — TypeFritschCarlsonMonotonicInterpolationMonotonic interpolation based on Fritsch & Carlson (1980), "Monotone Piecewise Cubic Interpolation", doi:10.1137/0717021.
Tangents are first initialized, then adjusted if they are not monotonic can overshoot for non-monotonic data
Interpolations.InPlace — TypeInPlace(gt) is a boundary condition that allows prefiltering to occur in-place (it typically requires padding)
Interpolations.InPlaceQ — TypeInPlaceQ(gt) is similar to InPlace(gt), but is exact when the values being interpolated arise from an underlying quadratic. It is primarily useful for testing purposes, allowing near-exact (to machine precision) comparisons against ground truth.
Interpolations.Line — TypeLine(gt) uses a constant slope for extrapolation
Interpolations.Linear — TypeLinear()Indicate that the corresponding axis should use linear interpolation.
Extended help
Assuming uniform knots with spacing 1, the ith piece of linear b-spline implemented here is defined as follows.
y_i(x) = c p(x) + cp p(1-x)where
p(δx) = xand the values cX for X ∈ {_, p} are the coefficients.
Linear b-splines are naturally interpolating, and require no prefiltering; there is therefore no need for boundary conditions to be provided.
Also, although the implementation is slightly different in order to re-use the framework built for general b-splines, the resulting interpolant is just a piecewise linear function connecting each pair of neighboring data points.
Interpolations.LinearMonotonicInterpolation — TypeLinearMonotonicInterpolationSimple linear interpolation.
Interpolations.NoInterp — TypeNoInterp() indicates that the corresponding axis must use integer indexing (no interpolation is to be performed)
Interpolations.OnCell — TypeOnCell() indicates that the boundary condition applies a half-gridspacing beyond the first & last nodes
Interpolations.OnGrid — TypeOnGrid() indicates that the boundary condition applies at the first & last nodes
Interpolations.Periodic — TypePeriodic(gt) applies periodic boundary conditions
Interpolations.Quadratic — TypeQuadratic(bc::BoundaryCondition)Indicate that the corresponding axis should use quadratic interpolation.
Extended help
Assuming uniform knots with spacing 1, the ith piece of quadratic spline implemented here is defined as follows:
y_i(x) = cm p(x-i) + c q(x) + cp p(1-(x-i))where
p(δx) = (δx - 1)^2 / 2
q(δx) = 3/4 - δx^2and the values for cX for X ∈ {m,_,p} are the pre-filtered coefficients.
For future reference, this expands to the following polynomial:
y_i(x) = cm * 1/2 * (x-i-1)^2 + c * (3/4 - x + i)^2 + cp * 1/2 * (x-i)^2When we derive boundary conditions we will use derivatives y_1'(x-1) and y_1''(x-1)
Interpolations.Reflect — TypeReflect(gt) applies reflective boundary conditions
Interpolations.SteffenMonotonicInterpolation — TypeSteffenMonotonicInterpolationMonotonic interpolation based on Steffen (1990), "A Simple Method for Monotonic Interpolation in One Dimension", http://adsabs.harvard.edu/abs/1990A%26A...239..443S
Only one pass, results usually between FritschCarlson and FritschButland.
Interpolations.Throw — TypeThrow(gt) causes beyond-the-edge extrapolation to throw an error
Internal API
Interpolations.boundstep — FunctionReturns half the width of one step of the range.
This function is used to calculate the upper and lower bounds of OnCell interpolation objects.
Interpolations.count_interp_dims — Methodn = count_interp_dims(ITP)Count the number of dimensions along which type ITP is interpolating. NoInterp dimensions do not contribute to the sum.
Interpolations.gradient_weights — Functionw = gradient_weights(degree, δx)Compute the weights for interpolation of the gradient at an offset δx from the "base" position. degree describes the interpolation scheme.
Example
julia> Interpolations.gradient_weights(Linear(), 0.2)
(-1.0, 1.0)This defines the gradient of a linear interpolation at 3.2 as y[4] - y[3].
Interpolations.hessian_weights — Functionw = hessian_weights(degree, δx)Compute the weights for interpolation of the hessian at an offset δx from the "base" position. degree describes the interpolation scheme.
Example
julia> Interpolations.hessian_weights(Linear(), 0.2)
(0.0, 0.0)Linear interpolation uses straight line segments, so the second derivative is zero.
Interpolations.inner_system_diags — Functiondl, d, du = inner_system_diags{T,IT}(::Type{T}, n::Int, ::Type{IT})Helper function to generate the prefiltering equation system: generates the diagonals for a n-by-n tridiagonal matrix with eltype T corresponding to the interpolation type IT.
dl, d, and du are intended to be used e.g. as in M = Tridiagonal(dl, d, du)
Interpolations.inner_system_diags — MethodCubic: continuity in function value, first and second derivatives yields
2/3 1/6
1/6 2/3 1/6
1/6 2/3 1/6
⋱ ⋱ ⋱Interpolations.prefiltering_system — FunctionM, b = prefiltering_system{T,TC,GT<:GridType,D<:Degree}m(::T, ::Type{TC}, n::Int, ::Type{D}, ::Type{GT})Given element types (T, TC) and interpolation scheme (GT, D) as well the number of rows in the data input (n), compute the system used to prefilter spline coefficients. Boundary conditions determine the values on the first and last rows.
Some of these boundary conditions require that these rows have off-tridiagonal elements (e.g the [1,3] element of the matrix). To maintain the efficiency of solving tridiagonal systems, the Woodbury matrix identity is used to add additional elements off the main 3 diagonals.
The filtered coefficients are given by solving the equation system
M * c = v + bwhere c are the sought coefficients, and v are the data points.
Interpolations.prefiltering_system — MethodQuadratic{Flat} OnCell and Quadratic{Reflect} OnCell amounts to setting y_1'(x) = 0 at x=1/2. Applying this condition yields
-cm + c = 0Interpolations.prefiltering_system — MethodQuadratic{Flat} OnGrid and Quadratic{Reflect} OnGrid amount to setting y_1'(x) = 0 at x=1. Applying this condition yields
-cm + cp = 0Interpolations.prefiltering_system — MethodCubic{Free} OnGrid and Cubic{Free} OnCell amount to requiring an extra continuous derivative at the second-to-last cell boundary; this means y_1'''(2) = y_2'''(2), yielding
1 cm -3 c + 3 cp -1 cpp = 0Interpolations.prefiltering_system — MethodCubic{Periodic} OnGrid closes the system by looking at the coefficients themselves as periodic, yielding
c0 = c(N+1)where N is the number of data points.
Interpolations.prefiltering_system — MethodCubic{Flat}, OnCell amounts to setting y_1'(x) = 0 at x = 1/2. Applying this condition yields
-9/8 cm + 11/8 c - 3/8 cp + 1/8 cpp = 0or, equivalently,
-9 cm + 11 c -3 cp + 1 cpp = 0(Note that we use y_1'(x) although it is strictly not valid in this domain; if we were to use y_0'(x) we would have to introduce new coefficients, so that would not close the system. Instead, we extend the outermost polynomial for an extra half-cell.)
Interpolations.prefiltering_system — MethodCubic{Flat} OnGrid amounts to setting y_1'(x) = 0 at x = 1. Applying this condition yields
-cm + cp = 0Interpolations.prefiltering_system — MethodCubic{Line} OnCell amounts to setting y_1''(x) = 0 at x = 1/2. Applying this condition yields
3 cm -7 c + 5 cp -1 cpp = 0(Note that we use y_1'(x) although it is strictly not valid in this domain; if we were to use y_0'(x) we would have to introduce new coefficients, so that would not close the system. Instead, we extend the outermost polynomial for an extra half-cell.)
Interpolations.prefiltering_system — MethodCubic{Line} OnGrid amounts to setting y_1''(x) = 0 at x = 1. Applying this condition gives:
1 cm -2 c + 1 cp = 0Interpolations.prefiltering_system — MethodQuadratic{Free} OnGrid and Quadratic{Free} OnCell amount to requiring an extra continuous derivative at the second-to-last cell boundary; this means that y_1''(3/2) = y_2''(3/2), yielding
1 cm -3 c + 3 cp - cpp = 0Interpolations.prefiltering_system — MethodQuadratic{Line} OnGrid and Quadratic{Line} OnCell amount to setting y_1''(x) = 0 at x=1 and x=1/2 respectively. Since y_i''(x) is independent of x for a quadratic b-spline, these both yield
1 cm -2 c + 1 cp = 0Interpolations.prefiltering_system — MethodQuadratic{Periodic} OnGrid and Quadratic{Periodic} OnCell close the system by looking at the coefficients themselves as periodic, yielding
c0 = c(N+1)where N is the number of data points.
Interpolations.rescale_gradient — Functionrescale_gradient(r::AbstractRange)
Implements the chain rule dy/dx = dy/du * du/dx for use when calculating gradients with scaled interpolation objects.
Interpolations.show_ranged — Methodshow_ranged(io, X, knots)A replacement for the default array-show for types that may not have the canonical evaluation points. rngs is the tuple of knots along each axis.
Interpolations.value_weights — Functionw = value_weights(degree, δx)Compute the weights for interpolation of the value at an offset δx from the "base" position. degree describes the interpolation scheme.
Example
julia> Interpolations.value_weights(Linear(), 0.2)
(0.8, 0.2)This corresponds to the fact that linear interpolation at x + 0.2 is 0.8*y[x] + 0.2*y[x+1].
Interpolations.weightedindexes — Methodweightedindexes(fs, itpflags, nodes, xs)Compute WeightedIndex values for evaluation at the position xs.... fs is a function or tuple of functions indicating the types of index required, typically value_weights, gradient_weights, and/or hessian_weights. itpflags and nodes can be obtained from itpinfo(itp)....
See the "developer documentation" for further information.
Interpolations.BSplineInterpolation — TypeBSplineInterpolation{T,N,TCoefs,IT,Axs}An interpolant-type for b-spline interpolation on a uniform grid with integer nodes. T indicates the element type for operations like collect(itp), and may also agree with the values obtained from itp(x, y, ...) at least for certain types of x and y. N is the dimensionality of the interpolant. The remaining type-parameters describe the types of fields:
- the
coefsfield holds the interpolation coefficients. Depending on prefiltering, these may or may not be the same as the supplied array of interpolant values. parentaxesholds the axes of the parent. Depending on prefiltering this may be "narrower" than the axes ofcoefs.itholds the interpolation type, e.g.,BSpline(Linear())or(BSpline(Quadratic(OnCell()),BSpline(Linear())).
BSplineInterpolation objects are typically created with interpolate. However, for customized control you may also construct them with
BSplineInterpolation(TWeights, coefs, it, axs)where T gets computed from the product of TWeights and eltype(coefs). (This is equivalent to indicating that you'll be evaluating at locations itp(x::TWeights, y::TWeights, ...).)
Interpolations.BoundaryCondition — TypeBoundaryConditionAn abstract type with one of the following values (see the help for each for details):
Throw(gt)Flat(gt)Line(gt)Free(gt)Periodic(gt)Reflect(gt)InPlace(gt)InPlaceQ(gt)
Interpolations.BoundsCheckStyle — TypeBoundsCheckStyle(itp)A trait to determine dispatch of bounds-checking for itp. Can return NeedsCheck(), in which case bounds-checking is performed, or CheckWillPass() in which case the check will return true.
Interpolations.MonotonicInterpolation — TypeMonotonicInterpolationMonotonic interpolation up to third order represented by type, knots and coefficients.
Interpolations.MonotonicInterpolationType — TypeMonotonicInterpolationTypeAbstract class for all types of monotonic interpolation.
Interpolations.WeightedIndex — Typewi = WeightedIndex(indexes, weights)Construct a weighted index wi, which can be thought of as a generalization of an ordinary array index to the context of interpolation. For an ordinary vector a, a[i] extracts the element at index i. When interpolating, one is typically interested in a range of indexes and the output is some weighted combination of array values at these indexes. For example, for linear interpolation at a location x between integers i and i+1, we have
ret = (1-f)*a[i] + f*a[i+1]where f = x-i lies between 0 and 1. This can be represented as a[wi], where
wi = WeightedIndex(i:i+1, (1-f, f))i.e.,
ret = sum(a[indexes] .* weights)Linear interpolation thus constructs weighted indices using a 2-tuple for weights and a length-2 indexes range. Higher-order interpolation would involve more positions and weights (e.g., 3-tuples for quadratic interpolation, 4-tuples for cubic).
In multiple dimensions, separable interpolation schemes are implemented in terms of multiple weighted indices, accessing A[wi1, wi2, ...] where each wi is the WeightedIndex along the corresponding dimension.
For value interpolation, weights will typically sum to 1. However, for gradient and Hessian computation this will not necessarily be true. For example, the gradient of one-dimensional linear interpolation can be represented as
gwi = WeightedIndex(i:i+1, (-1, 1))
g1 = a[gwi]For a three-dimensional array A, one might compute ∂A/∂x₂ (the second component of the gradient) as A[wi1, gwi2, wi3], where wi1 and wi3 are "value" weights and gwi2 "gradient" weights.
indexes may be supplied as a range or as a tuple of the same length as weights. The latter is applicable, e.g., for periodic boundary conditions.