Data interpolation

High-level interpolation and fitting of gridded data using B-splines.

Functions

BSplineKit.SplineInterpolations.interpolateFunction
interpolate(x, y, BSplineOrder(k), [bc = nothing])

Interpolate values y at locations x using B-splines of order k.

Grid points x must be real-valued and are assumed to be in increasing order.

Returns a SplineInterpolation which can be evaluated at any intermediate point.

Optionally, one may pass one of the boundary conditions listed in the Boundary conditions section. Currently, the Natural and Periodic boundary conditions are available.

See also interpolate!.

Periodic boundary conditions

Periodic boundary conditions should be used if the interpolated data is supposed to represent a periodic signal. In this case, pass bc = Period(L), where L is the period of the x-axis. Note that the endpoint x[begin] + L should not be included in the x vector.

Cubic periodic splines

Cubic periodic splines (BSplineOrder(4)) are particularly well optimised compared to periodic splines of other orders. Just note that interpolations using cubic periodic splines modify their input (including x and y values).

Examples

julia> xs = -1:0.1:1;

julia> ys = cospi.(xs);

julia> S = interpolate(xs, ys, BSplineOrder(4))
SplineInterpolation containing the 21-element Spline{Float64}:
 basis: 21-element BSplineBasis of order 4, domain [-1.0, 1.0]
 order: 4
 knots: [-1.0, -1.0, -1.0, -1.0, -0.8, -0.7, -0.6, -0.5, -0.4, -0.3  …  0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 1.0, 1.0, 1.0, 1.0]
 coefficients: [-1.0, -1.00111, -0.8975, -0.597515, -0.314147, 1.3265e-6, 0.314142, 0.597534, 0.822435, 0.96683  …  0.96683, 0.822435, 0.597534, 0.314142, 1.3265e-6, -0.314147, -0.597515, -0.8975, -1.00111, -1.0]
 interpolation points: -1.0:0.1:1.0

julia> S(-1)
-1.0

julia> (Derivative(1) * S)(-1)
-0.01663433622896893

julia> (Derivative(2) * S)(-1)
10.527273287554928

julia> Snat = interpolate(xs, ys, BSplineOrder(4), Natural())
SplineInterpolation containing the 21-element Spline{Float64}:
 basis: 21-element RecombinedBSplineBasis of order 4, domain [-1.0, 1.0], BCs {left => (D{2},), right => (D{2},)}
 order: 4
 knots: [-1.0, -1.0, -1.0, -1.0, -0.9, -0.8, -0.7, -0.6, -0.5, -0.4  …  0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0, 1.0, 1.0, 1.0]
 coefficients: [-0.833333, -0.647516, -0.821244, -0.597853, -0.314057, -2.29076e-5, 0.314148, 0.597532, 0.822435, 0.96683  …  0.96683, 0.822435, 0.597532, 0.314148, -2.29076e-5, -0.314057, -0.597853, -0.821244, -0.647516, -0.833333]
 interpolation points: -1.0:0.1:1.0

julia> Snat(-1)
-1.0

julia> (Derivative(1) * Snat)(-1)
0.2872618670889516

julia> (Derivative(2) * Snat)(-1)
-3.33066907387547e-14

Periodic boundary conditions

Interpolate $f(x) = \cos(πx)$ for $x ∈ [-1, 1)$. Note that the period is $L = 2$ and that the endpoint ($x = 1$) must not be included in the data points.

julia> xp = -1:0.1:0.9;

julia> yp = cospi.(xp);

julia> Sper = interpolate(xp, yp, BSplineOrder(4), Periodic(2))
SplineInterpolation containing the 20-element Spline{Float64}:
 basis: 20-element PeriodicBSplineBasis of order 4, domain [-1.0, 1.0), period 2.0
 order: 4
 knots: [..., -1.2, -1.1, -1.0, -0.9, -0.8, -0.7, -0.6, -0.5, -0.4, -0.3  …  0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0, 1.1, ...]
 coefficients: [..., -1.01659, -0.96683, -0.822435, -0.597534, -0.314142, 1.10589e-17, 0.314142, 0.597534, 0.822435, 0.96683, 1.01659, 0.96683, 0.822435, 0.597534, 0.314142, 1.51788e-17, -0.314142, -0.597534, -0.822435, -0.96683, ...]
 interpolation points: -1.0:0.1:0.9

As expected, the periodic spline does a better job at approximating the periodic function $f(x) = \cos(πx)$ near the boundaries than the other interpolations:

julia> x = -0.99; cospi(x), Sper(x), Snat(x), S(x)
(-0.9995065603657316, -0.9995032595823043, -0.9971071640321145, -0.9996420091470221)

julia> x = 0.998; cospi(x), Sper(x), Snat(x), S(x)
(-0.9999802608561371, -0.9999801044078943, -0.9994253145274461, -1.0000122303614758)
source
BSplineKit.SplineInterpolations.fitFunction
fit(xs, ys, λ::Real, [BSplineOrder(4)], [bc = Natural()]; [weights = nothing])

Fit a cubic smoothing spline to the given data.

Returns a cubic spline which roughly passes through the data (points (xs[i], ys[i])) given some smoothing parameter $λ$. Note that $λ = 0$ means no smoothing and the results are equivalent to all the data.

One can optionally pass a weights vector if one desires to give different weights to different data points (e.g. if one wants the curve to pass as close as possible to a specific point). By default all weights are $w_i = 1$.

More precisely, the returned spline $S(x)$ minimises:

\[∑_{i = 1}^N w_i |y_i - S(x_i)|^2 + λ ∫_{x_1}^{x_N} \left[ S''(x) \right]^2 \, \mathrm{d}x\]

Only cubic splines (BSplineOrder(4)) are currently supported. Moreover, the boundary condition (bc) must be Periodic for a periodic spline or Natural otherwise (this is the default). (Currently, the periodic case can be much slower than the default natural condition.)

Examples

julia> xs = (0:0.01:1).^2;

julia> ys = @. cospi(2 * xs) + 0.1 * sinpi(200 * xs);  # smooth + highly fluctuating components

julia> λ = 0.001;  # smoothing parameter

julia> S = fit(xs, ys, λ)
101-element Spline{Float64}:
 basis: 101-element RecombinedBSplineBasis of order 4, domain [0.0, 1.0], BCs {left => (D{2},), right => (D{2},)}
 order: 4
 knots: [0.0, 0.0, 0.0, 0.0, 0.0001, 0.0004, 0.0009, 0.0016, 0.0025, 0.0036  …  0.8836, 0.9025, 0.9216, 0.9409, 0.9604, 0.9801, 1.0, 1.0, 1.0, 1.0]
 coefficients: [0.946872, 0.631018, 1.05101, 1.04986, 1.04825, 1.04618, 1.04366, 1.04067, 1.03722, 1.03331  …  0.437844, 0.534546, 0.627651, 0.716043, 0.798813, 0.875733, 0.947428, 1.01524, 0.721199, 0.954231]
source

Types

BSplineKit.SplineInterpolations.SplineInterpolationType
SplineInterpolation

Spline interpolation.

This is the type returned by interpolate.

A SplineInterpolation I can be evaluated at any point x using the I(x) syntax.

It can also be updated with new data on the same data points using interpolate!.


SplineInterpolation(undef, B::AbstractBSplineBasis, x::AbstractVector, [T = eltype(x)])

Initialise a SplineInterpolation from B-spline basis and a set of interpolation (or collocation) points x.

Note that the length of x must be equal to the number of B-splines.

Use interpolate! to actually interpolate data known on the x locations.

source