Reference
The following provides detailled documentation about types and methods provided by the InterpolationKernels
package. This information is also available from the REPL by typing ?
followed by the name of a method or a type.
Kernel supports
InterpolationKernels.Support
— TypeInterpolationKernels.Support{T,S,L,R}
is the abstract type for the support of an interpolation kernel parameterized by the floating-point type T
used for computations, the integer size S
of the support and the types L
and R
of the left and right bounds which can be InterpolationKernels.Open
or InterpolationKernels.Closed
.
InterpolationKernels.SymmetricSupport
— TypeInterpolationKernels.SymmetricSupport{T,S,L,R}() -> sup
yields an instance of a symmetric support parameterized by the floating-point type T
, the integer size S
of the support and the types L
and R
of the left and right bounds which can be InterpolationKernels.Open
or InterpolationKernels.Closed
.
Depending on L
and R
, the support is:
(L,R) = (Open,Open) -------> sup = (-S/2,S/2)
(L,R) = (Closed,Open) -----> sup = [-S/2,S/2)
(L,R) = (Open,Closed) -----> sup = (-S/2,S/2]
(L,R) = (Closed,Closed) ---> sup = [-S/2,S/2]
InterpolationKernels.LeftAnchoredSupport
— TypeInterpolationKernels.LeftAnchoredSupport{T,S,L,R}(a)
yields an instance of a support with lower bound a
and parameterized by the floating-point type T
, the integer size S
of the support and the types L
and R
of the left and right bounds which can be InterpolationKernels.Open
or InterpolationKernels.Closed
.
Depending on L
and R
, the support is:
(L,R) = (Open,Open) -------> sup = (a,a+S)
(L,R) = (Closed,Open) -----> sup = [a,a+S)
(L,R) = (Open,Closed) -----> sup = (a,a+S]
(L,R) = (Closed,Closed) ---> sup = [a,a+S]
InterpolationKernels.RightAnchoredSupport
— TypeInterpolationKernels.RightAnchoredSupport{T,S,L,R}(b)
yields an instance of a support with upper bound b
and parameterized by the floating-point type T
, the integer size S
of the support and the types L
and R
of the left and right bounds which can be InterpolationKernels.Open
or InterpolationKernels.Closed
.
Depending on L
and R
, the support is:
(L,R) = (Open,Open) -------> sup = (b-S,b)
(L,R) = (Closed,Open) -----> sup = [b-S,b)
(L,R) = (Open,Closed) -----> sup = (b-S,b]
(L,R) = (Closed,Closed) ---> sup = [b-S,b]
InterpolationKernels.infimum
— FunctionInterpolationKernels.infimum(sup) -> a
yields the greatest lower bound a
of the kernel support sup
.
InterpolationKernels.supremum
— FunctionInterpolationKernels.supremum(sup) -> b
yields the least upper bound b
of the kernel support sup
.
Interpolation kernels
InterpolationKernels.Kernel
— TypeInterpolation Kernels
An interpolation kernel Kernel{T,S}
is parameterized by the floating-point type T
of its coefficients and by the integer size S
of its support. For efficiency reasons, only kernels with (small) finite size supports are implemented in InterpolationKernels
.
A kernel ker
is a callable object which may be used as a function with a real argument:
ker(x::Real)
yields kernel value at offset x
. Whatever the type of x
, ker(x)
is always of type T = eltype(ker)
the floating-point type associated with ker
. All kernel supports are symmetric; that is ker(x)
is zero if abs(x) > S/2
with S = length(ker)
the size of the kernel support.
Kernel floating-point type conversion
Called as a function with a real argument, a given kernel returns a value of its associated floating-point type. This has been chosen to have fast interpolation methods. Converting a kernel ker
to use floating-point type T
is simply done by one of:
T(ker)
Kernel{T}(ker)
convert(Kernel{T}, ker)
Beware that changing the floating-point type may lead to a loss of precision if the kernel has numerical parameters.
Common methods
A few common methods are specialized for any interpolation kernel ker
:
eltype(ker) -> T
yields the floating-point type for calculations,
length(ker) -> S
yield the size the support of ker
which is also the number of neighbors involved in an interpolation by this kernel,
values(ker)
yields a tuple of the parameters of ker
such that an identical instance can be built by:
typeof(ker)(values(ker)...)
finally:
compute_offset_and_weights(ker, x) -> off, (w1, w2, ..., wS)
yields the offset off
and an S
-tuple of interpolation weights to interpolate an array at coordinate x
(in fractional index units).
InterpolationKernels.iscardinal
— Functioniscardinal(ker)
yields whether the kernel ker
is zero for any non-zero integer argument. Cardinal kernels are directly suitable for interpolation.
InterpolationKernels.isnormalized
— Functionisnormalized(ker)
yields whether the kernel ker
has the partition of unity property. That is, the sum of the values computed by the kernel ker
on a unit spaced grid is equal to one.
InterpolationKernels.support
— FunctionInterpolationKernels.support(ker) -> sup
yields the support of the interpolation kernel ker
.
InterpolationKernels.brief
— FunctionInterpolationKernels.brief(ker)
yields a brief description of the kernel type or instance ker
.
InterpolationKernels.compute_weights
— FunctionInterpolationKernels.compute_weights(ker, t) -> wgt
computes the interpolation weights returned by InterpolationKernels.compute_offset_and_weights
for kernel ker
with symmetric support. Assuming interpolation is performed at position x
, argument t
is given by:
t = x - floor(x) if length(ker) is even, hence t ∈ [0,1)
t = x - round(x) if length(ker) is odd, hence t ∈ [-1/2,+1/2]
The returned weights are then:
wgt = ntuple(i -> ker(t + k - i), length(ker))
where k = (length(ker) + 1) >> 1
(i.e., integer division of length(ker)+1
by 2). These conventions have been adopted so that, by specializing the compute_weights
method, computing the length(ker)
weights at the same time may be done in much fewer operations than calling ker
as a function for each weight.
InterpolationKernels.compute_offset_and_weights
— FunctionInterpolationKernels.compute_offset_and_weights(ker, x) -> off, wgt
yields the index offset off
and the weights wgt
to interpolate with kernel ker
at position x
in fractional index units. The offset is a scalar and the weights are an n
-tuple with n = length(ker)
the size of the support of the kernel, all returned values have the same floating point type eltype(ker)
as the kernel.
Not taking into account boundary conditions, interpolating a vector A
at position x
would then write:
off, wgt = InterpolationKernels.compute_offset_and_weights(ker, x)
k = Int(off)
result = wgt[1]*A[k+1] + ... + wgt[n]*A[k+n]
Note that 1-based indexing is assumed by compute_offset_and_weights
to interpret the position x
and compute the offset off
. If this is not the case, the code should be:
j1 = first(axes(A,1)) # first index in A
off, wgt = InterpolationKernels.compute_offset_and_weights(ker, x - (j1 - 1))
k = Int(off) + (j1 - 1)
result = wgt[1]*A[k+1] + ... + wgt[n]*A[k+n]
where expression x - (j1 - 1)
is assuming that the position x
is in fractional index for A
, that is x = j1
at the first entry of A
.
For fast computations, this method should be specialized for specific kernel types. For kernels with symmetric support, the method InterpolationKernels.compute_weights
is called by compute_offset_and_weights
to calculate the interpolation weights; for such kernels it is sufficient to specialize compute_weights
instead of compute_offset_and_weights
.
B-splines
InterpolationKernels.BSpline
— TypeBSpline{S,T}()
yields a B-spline (short for basis spline) of order S
that is a piecewise polynomial function of degree S - 1
on a support of length S
. The parameter T
is the floating-point type for computations, T = Float64
is assuled if this parameter is not specified.
Fr now, not all B-spline are implemented in InterpolationKernels
, S
must be: 1
(for a rectangular B-spline), 2
(for a linear B-spline), 3
(for a quadratic B-spline), or 4
(for a cubic B-spline).
If ker
is a B-spline, then ker'
is its derivative which can also be directly constructed by calling BSplinePrime
.
The derivative of B-splines of order S ≤ 2
is not defined everywhere. It is allowed to take their derivative but it (arbitrarily) yields zero where not defined. Returning NaN
would have been more correct but it has been considered that it would do more harm than good in practice.
InterpolationKernels.BSplinePrime
— TypeBSplinePrime{S,T}()
yields the derivative of a B-spline of order S
for floating-point T
.
See the caveats in BSpline
about taking the derivative of B-splines of order S ≤ 2
.
Generic cubic spline
InterpolationKernels.CubicSpline
— TypeCubicSpline{T}(a, b) -> ker
yields an instance of a generic cubic spline for floating-point type T
and parameters a = ker'(1)
and b = ker(1)
the slope and the value of the function ker(x)
at x = 1
.
A cubic spline kernel is at least C¹ continuous, the expression ker'
yields a kernel instance implementing the 1st derivative of the generic cubic spline ker
(see CubicSplinePrime
to directly build a derivative).
Depending on the values of the parameters a
and b
, more specific cubic spline kernels can be emulated:
CubicSpline{T}(-1/2,1/6)
yields a cubic B-spline as built byBSpline{4,T}()
.CubicSpline{T}(a,0)
yields a cardinal cubic spline as built byCardinalCubicSpline{T}(a)
.CubicSpline{T}(-1/2,0)
yields a Catmull-Rom kernel as built byCatmullRomSpline{T}()
.CubicSpline{T}(-b/2-c,b/6)
yields Mitchell & Netravali cubic spline as built byMitchellNetravaliSpline{T}(b,c)
.
Instances of CubicSpline
are very well optimized and, in practice, they may be as fast or even faster than these more specialized counterparts.
InterpolationKernels.CubicSplinePrime
— TypeCubicSplinePrime{T}(a, b)
yields a kernel instance that is the 1st derivative of the generic cubic spline of parameters a
and b
(see CubicSpline
) for floating-point type T
(Float64
by default).
Cardinal cubic splines
InterpolationKernels.CardinalCubicSpline
— TypeCardinalCubicSpline{T}(a)
yields an instance of the Keys family of cardinal cubic splines for floating-point type T
and parameter a = ker'(1)
the slope of the function ker(x)
at x = 1
.
These kernels are C¹ continuous piecewise normalized cardinal cubic spline which depend on one parameter a
and defined by:
ker(x) = 1 + ((2 + a)*|x| - (3 + a))*x^2 if |x| ≤ 1
a*(|x| - 1)*(|x| - 2)^2 if 1 ≤ |x| ≤ 2
0 if |x| ≥ 2
The expression ker'
yields a kernel instance which is the 1st derivative of the Keys kernel ker
(also see the constructor CardinalCubicSplinePrime
).
Reference:
- Keys, Robert, G., "Cubic Convolution Interpolation for Digital Image Processing", IEEE Trans. Acoustics, Speech, and Signal Processing, Vol. ASSP-29, No. 6, December 1981, pp. 1153-1160.
InterpolationKernels.CardinalCubicSplinePrime
— TypeCardinalCubicSplinePrime{T}(a)
yields a kernel instance that is the 1st derivative of the Keys cardinal cubic spline (see CardinalCubicSpline
) for floating-point type T
and parameter a
. This derivative is given by:
ker′(x) = (3(a + 2)*|x| - 2(a + 3))*x if |x| ≤ 1
(3a)*(|x| - 2)*(|x| - 4/3)*sign(x) if 1 ≤ |x| ≤ 2
0 if |x| ≥ 2
Catmull-Rom interpolation kernel
InterpolationKernels.CatmullRomSpline
— TypeCatmullRomSpline{T}()
yields an instance of the Catmull-Rom interpolation kernel for floating-point type T
which is assumed to be Float64
if omitted.
Catmull-Rom interpolation kernel is a cardinal piecewise cubic spline defined by:
ker(x) = ((3/2)*|x| - (5/2))*x^2 + 1 if |x| ≤ 1
(((5/2) - (1/2)*|x|)*|x| - 4)*|x| + 2 if 1 ≤ |x| ≤ 2
0 if |x| ≥ 2
The expression ker'
yields a kernel instance which is the 1st derivative of the Catmull-Rom interpolation kernel ker
(also see the constructor CatmullRomSplinePrime
).
InterpolationKernels.CatmullRomSplinePrime
— TypeCatmullRomSplinePrime{T}()
yields a kernel instance that is the 1st derivative of the Catmull-Rom interpolation kernel (see CatmullRomSpline
) for floating-point type T
which is assumed to be Float64
if omitted.
The 1st derivative of the Catmull-Rom interpolation kernel is given by:
ker′(x) = ((9/2)*|x| - 5)*x if a = |x| ≤ 1
(5 - (3/2)*|x|)*x - 4*sign(x) if 1 ≤ |x| ≤ 2
0 if |x| ≥ 2
Mitchell-Netravali kernels
InterpolationKernels.MitchellNetravaliSpline
— TypeMitchellNetravaliSpline{T}(b=1/3, c=1/3)
yields an instance of the Mitchell & Netravali family of kernels for floating-point type T
and parameters (b,c)
.
These kernels are cubic splines which depends on 2 parameters, b
and c
. Whatever the values of (b,c)
, Mitchell & Netravali kernels are normalized, even and C¹ continuous functions (these kernels and their first derivatives are continuous).
Taking b = 0
yields the family of cardinal cubic splines (see CardinalCubicSpline
) and is a sufficient and necessary condition to have Mitchell & Netravali kernels be cardinal functions.
Using the constraint: b + 2c = 1
yields a cubic filter with, at least, quadratic order approximation.
Some specific values of (b,c)
yield other well known kernels:
(b,c) = (1,0) --> cubic B-spline
(b,c) = (0,-a) --> Keys's cardinal cubic spline CardinalCubicSpline(a)
(b,c) = (0,1/2) --> Catmull-Rom kernel CatmullRomSpline()
(b,c) = (b,0) --> Duff's tensioned B-spline
(b,c) = (6β,-α-3β) --> generic cubic spline CubicSpline(α,β)
(b,c) = (1/3,1/3) --> recommended by Mitchell-Netravali
The expression ker'
yields a kernel instance which is the 1st derivative of the Mitchell & Netravali kernel ker
(also see the constructor MitchellNetravaliSplinePrime
).
Mitchell & Netravali family of kernels are currently instances of CubicSpline
.
Reference:
InterpolationKernels.MitchellNetravaliSplinePrime
— TypeMitchellNetravaliSplinePrime([T=Float64,] [b=1/3, c=1/3,] B=Flat)
yields a kernel instance that is the 1st derivative of the Mitchell & Netravali kernel (see MitchellNetravaliSpline
) for floating-point type T
, parameters b
and c
and boundary conditions B
.
Lanczos re-sampling kernels
InterpolationKernels.LanczosKernel
— TypeLanczosKernel{S,T}()
yields an instance of a Lanczos re-sampling kernel of support size S
(which must be even) and for floating-point type T
.
The Lanczos re-sampling kernels are even cardinal functions which tend to be normalized for large support size. They are defined by:
ker(x) = S/(2*(π*x)^2)*sin(π*x)*sin(2*π*x/S) if |x| ≤ S/2
0 if |x| ≥ S/2
The expression ker'
yields the first derivative of a Lanczos re-sampling kernel ker
(also see the constructor LanczosKernelPrime
).
InterpolationKernels.LanczosKernelPrime
— TypeLanczosKernelPrime{S,T}()
yields a kernel instance that is the 1st derivative of the Lanczos re-sampling kernel (see LanczosKernel
) of support size S
and for floating-point type T
.