diff --git a/docs/src/index.md b/docs/src/index.md index 36cb305f..707ed396 100644 --- a/docs/src/index.md +++ b/docs/src/index.md @@ -20,6 +20,7 @@ corresponding to `(u,t)` pairs. - `SmoothedConstantInterpolation(u,t)` - An integral preserving continuously differentiable approximation of constant interpolation. - `LinearInterpolation(u,t)` - A linear interpolation. + - `SmoothedLinearInterpolation(u,t; λ = 0.25)` - A continuously differentiable linear interpolation with an interval around each data point replaced by spline section. - `QuadraticInterpolation(u,t)` - A quadratic interpolation. - `LagrangeInterpolation(u,t,n)` - A Lagrange interpolation of order `n`. - `QuadraticSpline(u,t)` - A quadratic spline interpolation. diff --git a/docs/src/manual.md b/docs/src/manual.md index f11c11d4..26740b40 100644 --- a/docs/src/manual.md +++ b/docs/src/manual.md @@ -2,6 +2,7 @@ ```@docs LinearInterpolation +SmoothedLinearInterpolation QuadraticInterpolation LagrangeInterpolation AkimaInterpolation diff --git a/docs/src/methods.md b/docs/src/methods.md index c5a3aeaa..70fe933a 100644 --- a/docs/src/methods.md +++ b/docs/src/methods.md @@ -25,6 +25,17 @@ A = LinearInterpolation(u, t) plot(A) ``` +## Smoothed Linear Interpolation + +This is a continuously differentiable linear interpolation with an interval around each data point replaced by spline section. + +```@example tutorial +A = LinearInterpolation(u, t) +plot(A) +A = SmoothedLinearInterpolation(u, t) +plot!(A) +``` + ## Quadratic Interpolation This function fits a parabola passing through the two nearest points from the input diff --git a/ext/DataInterpolationsSparseConnectivityTracerExt.jl b/ext/DataInterpolationsSparseConnectivityTracerExt.jl index 7f1df97b..120cca59 100644 --- a/ext/DataInterpolationsSparseConnectivityTracerExt.jl +++ b/ext/DataInterpolationsSparseConnectivityTracerExt.jl @@ -3,26 +3,24 @@ module DataInterpolationsSparseConnectivityTracerExt using SparseConnectivityTracer: AbstractTracer, Dual, primal, tracer using SparseConnectivityTracer: GradientTracer, gradient_tracer_1_to_1 using SparseConnectivityTracer: HessianTracer, hessian_tracer_1_to_1 -using FillArrays: Fill # from FillArrays.jl +using FillArrays: Fill using DataInterpolations: - AbstractInterpolation, - LinearInterpolation, - QuadraticInterpolation, - LagrangeInterpolation, - AkimaInterpolation, - ConstantInterpolation, - QuadraticSpline, - CubicSpline, - BSplineInterpolation, - BSplineApprox, - CubicHermiteSpline, - # PCHIPInterpolation, - QuinticHermiteSpline, - output_size + AbstractInterpolation, + AkimaInterpolation, + BSplineApprox, + BSplineInterpolation, + ConstantInterpolation, + CubicHermiteSpline, + CubicSpline, + LagrangeInterpolation, + LinearInterpolation, + QuadraticInterpolation, + QuadraticSpline, + QuinticHermiteSpline, + SmoothedLinearInterpolation, + output_size -#===========# # Utilities # -#===========# # Limit support to `u` begin an AbstractVector{<:Number} or AbstractMatrix{<:Number}, # to avoid any cases where the output size is dependent on the input value. @@ -33,8 +31,8 @@ function _sct_interpolate( uType::Type{<:AbstractVector{<:Number}}, t::GradientTracer, is_der_1_zero, - is_der_2_zero, - ) + is_der_2_zero +) return gradient_tracer_1_to_1(t, is_der_1_zero) end function _sct_interpolate( @@ -42,8 +40,8 @@ function _sct_interpolate( uType::Type{<:AbstractVector{<:Number}}, t::HessianTracer, is_der_1_zero, - is_der_2_zero, - ) + is_der_2_zero +) return hessian_tracer_1_to_1(t, is_der_1_zero, is_der_2_zero) end function _sct_interpolate( @@ -51,8 +49,8 @@ function _sct_interpolate( uType::Type{<:AbstractMatrix{<:Number}}, t::GradientTracer, is_der_1_zero, - is_der_2_zero, - ) + is_der_2_zero +) t = gradient_tracer_1_to_1(t, is_der_1_zero) N = only(output_size(interp)) return Fill(t, N) @@ -62,48 +60,47 @@ function _sct_interpolate( uType::Type{<:AbstractMatrix{<:Number}}, t::HessianTracer, is_der_1_zero, - is_der_2_zero, - ) + is_der_2_zero +) t = hessian_tracer_1_to_1(t, is_der_1_zero, is_der_2_zero) N = only(output_size(interp)) return Fill(t, N) end -#===========# # Overloads # -#===========# # We assume that with the exception of ConstantInterpolation and LinearInterpolation, # all interpolations have a non-zero second derivative at some point in the input domain. for (I, is_der1_zero, is_der2_zero) in ( - (:ConstantInterpolation, true, true), - (:LinearInterpolation, false, true), - (:QuadraticInterpolation, false, false), - (:LagrangeInterpolation, false, false), - (:AkimaInterpolation, false, false), - (:QuadraticSpline, false, false), - (:CubicSpline, false, false), - (:BSplineInterpolation, false, false), - (:BSplineApprox, false, false), - (:CubicHermiteSpline, false, false), - (:QuinticHermiteSpline, false, false), - ) + (:AkimaInterpolation, false, false), + (:BSplineApprox, false, false), + (:BSplineInterpolation, false, false), + (:ConstantInterpolation, true, true), + (:CubicHermiteSpline, false, false), + (:CubicSpline, false, false), + (:LagrangeInterpolation, false, false), + (:LinearInterpolation, false, true), + (:QuadraticInterpolation, false, false), + (:QuadraticSpline, false, false), + (:QuinticHermiteSpline, false, false), + (:SmoothedLinearInterpolation, false, false) +) @eval function (interp::$(I){uType})( t::AbstractTracer - ) where {uType <: AbstractArray{<:Number}} + ) where {uType <: AbstractArray{<:Number}} return _sct_interpolate(interp, uType, t, $is_der1_zero, $is_der2_zero) end end # Some Interpolations require custom overloads on `Dual` due to mutation of caches. for I in ( - :LagrangeInterpolation, - :BSplineInterpolation, - :BSplineApprox, - :CubicHermiteSpline, - :QuinticHermiteSpline, - ) + :LagrangeInterpolation, + :BSplineInterpolation, + :BSplineApprox, + :CubicHermiteSpline, + :QuinticHermiteSpline +) @eval function (interp::$(I){uType})(d::Dual) where {uType <: AbstractVector} p = interp(primal(d)) t = interp(tracer(d)) diff --git a/src/DataInterpolations.jl b/src/DataInterpolations.jl index 9b9dd8af..5b128cb5 100644 --- a/src/DataInterpolations.jl +++ b/src/DataInterpolations.jl @@ -118,10 +118,10 @@ _output_size(::AbstractVector{<:Number}) = () _output_size(u::AbstractVector) = size(first(u)) _output_size(u::AbstractArray) = Base.front(size(u)) -export LinearInterpolation, QuadraticInterpolation, LagrangeInterpolation, - AkimaInterpolation, ConstantInterpolation, SmoothedConstantInterpolation, - QuadraticSpline, CubicSpline, BSplineInterpolation, BSplineApprox, - CubicHermiteSpline, PCHIPInterpolation, QuinticHermiteSpline, +export LinearInterpolation, SmoothedLinearInterpolation, QuadraticInterpolation, + LagrangeInterpolation, AkimaInterpolation, ConstantInterpolation, + SmoothedConstantInterpolation, QuadraticSpline, CubicSpline, BSplineInterpolation, + BSplineApprox, CubicHermiteSpline, PCHIPInterpolation, QuinticHermiteSpline, SmoothArcLengthInterpolation, LinearInterpolationIntInv, ConstantInterpolationIntInv, ExtrapolationType export output_dim, output_size diff --git a/src/derivatives.jl b/src/derivatives.jl index b713540a..9ad04ae6 100644 --- a/src/derivatives.jl +++ b/src/derivatives.jl @@ -100,6 +100,26 @@ function _derivative(A::LinearInterpolation, t::Number, iguess) slope end +function _derivative(A::SmoothedLinearInterpolation, t::Number, iguess) + idx = get_idx(A, t, iguess) + + # Check against boundary points between linear and spline interpolation + left = (t < A.p.t_tilde[2idx]) + right = (t > A.p.t_tilde[2idx + 1]) + + # Spline interpolation + if left && idx != 1 + return U_deriv(A, t, idx) + end + + if right && idx != length(A.u) - 1 + return U_deriv(A, t, idx + 1) + end + + # Linear interpolation + return A.p.slope[idx + 1] +end + function _derivative(A::SmoothedConstantInterpolation{<:AbstractVector}, t::Number, iguess) idx = get_idx(A, t, iguess) d_lower, d_upper, c_lower, c_upper = get_parameters(A, idx) diff --git a/src/integrals.jl b/src/integrals.jl index f0e4afed..43b229d3 100644 --- a/src/integrals.jl +++ b/src/integrals.jl @@ -296,3 +296,7 @@ function _integral( integrate_quintic_polynomial(t1, t2, tᵢ, A.u[idx], A.du[idx], A.ddu[idx] / 2, c₁ + Δt * (-c₂ + c₃ * Δt), c₂ - 2c₃ * Δt, c₃) end + +function _integral(A::SmoothedLinearInterpolation, idx::Number, t1::Number, t2::Number) + throw(IntegralNotFoundError()) +end diff --git a/src/interpolation_caches.jl b/src/interpolation_caches.jl index 47c0773b..0024d7f4 100644 --- a/src/interpolation_caches.jl +++ b/src/interpolation_caches.jl @@ -64,6 +64,76 @@ function LinearInterpolation( cache_parameters, assume_linear_t) end +""" + SmoothedLinearInterpolation( + u, t; extrapolation::ExtrapolationType.T = ExtrapolationType.None, + extrapolation_left::ExtrapolationType.T = ExtrapolationType.None, + extrapolation_right::ExtrapolationType.T = ExtrapolationType.None, + cache_parameters = true, assume_linear_t = 1e-2, λ = 0.25) + + A linear interpolation method where a section around each data point is replaced by spline + section to make the transition at the datapoint differentiable. The spline section is guaranteed + to be in the convex hull of the boundary points of the spline and the original data point. + + ## Arguments + + - `u`: data points. + - `t`: time points. + +## Keyword Arguments + + - `extrapolation`: The extrapolation type applied left and right of the data. Possible options + are `ExtrapolationType.None` (default), `ExtrapolationType.Constant`, `ExtrapolationType.Linear` + `ExtrapolationType.Extension`, `ExtrapolationType.Periodic` and `ExtrapolationType.Reflective`. + - `extrapolation_left`: The extrapolation type applied left of the data. See `extrapolation` for + the possible options. This keyword is ignored if `extrapolation != Extrapolation.none`. + - `extrapolation_right`: The extrapolation type applied right of the data. See `extrapolation` for + the possible options. This keyword is ignored if `extrapolation != Extrapolation.none`. + - `cache_parameters`: precompute parameters at initialization for faster interpolation + computations. Note: if activated, `u` and `t` should not be modified. Defaults to `false`. + - `assume_linear_t`: boolean value to specify a faster index lookup behavior for + evenly-distributed abscissae. Alternatively, a numerical threshold may be specified + for a test based on the normalized standard deviation of the difference with respect + to the straight line (see [`looks_linear`](@ref)). Defaults to 1e-2. + - `λ`: The fraction of the interval on either side of each data point that is described by a spline section. Choose + in [0, 1.0], defaults to 0.25 +""" +struct SmoothedLinearInterpolation{uType, tType, IType, pType, T} <: + AbstractInterpolation{T} + u::uType + t::tType + I::IType + p::SmoothedLinearParameterCache{pType} + extrapolation_left::ExtrapolationType.T + extrapolation_right::ExtrapolationType.T + iguesser::Guesser{tType} + cache_parameters::Bool + linear_lookup::Bool + function SmoothedLinearInterpolation( + u, t, I, p, extrapolation_left, extrapolation_right, + cache_parameters, assume_linear_t) + linear_lookup = seems_linear(assume_linear_t, t) + new{typeof(u), typeof(t), typeof(I), typeof(p.slope), eltype(u)}( + u, t, I, p, extrapolation_left, extrapolation_right, + Guesser(t), cache_parameters, linear_lookup) + end +end + +function SmoothedLinearInterpolation( + u, t; extrapolation::ExtrapolationType.T = ExtrapolationType.None, + extrapolation_left::ExtrapolationType.T = ExtrapolationType.None, + extrapolation_right::ExtrapolationType.T = ExtrapolationType.None, + cache_parameters = true, assume_linear_t = 1e-2, λ = 0.25) + extrapolation_left, + extrapolation_right = munge_extrapolation( + extrapolation, extrapolation_left, extrapolation_right) + u, t = munge_data(u, t) + p = SmoothedLinearParameterCache(u, t, λ, cache_parameters) + return SmoothedLinearInterpolation( + u, t, nothing, p, extrapolation_left, + extrapolation_right, cache_parameters, assume_linear_t) +end + """ QuadraticInterpolation(u, t, mode = :Forward; extrapolation_left::ExtrapolationType.T = ExtrapolationType.None, extrapolation::ExtrapolationType.T = ExtrapolationType.None, extrapolation_right::ExtrapolationType.T = ExtrapolationType.None, diff --git a/src/interpolation_methods.jl b/src/interpolation_methods.jl index c30b0cff..39ab9174 100644 --- a/src/interpolation_methods.jl +++ b/src/interpolation_methods.jl @@ -129,6 +129,27 @@ function _interpolate(A::LinearInterpolation{<:AbstractArray}, t::Number, iguess return A.u[ax..., idx] + slope * Δt end +# Smoothed linear interpolation +function _interpolate(A::SmoothedLinearInterpolation, t::Number, iguess) + idx = get_idx(A, t, iguess) + + # Check against boundary points between linear and spline interpolation + left = (t < A.p.t_tilde[2idx]) + right = (t > A.p.t_tilde[2idx + 1]) + + # Spline interpolation + if left && idx != 1 + return U(A, t, idx) + end + + if right && idx != length(A.u) - 1 + return U(A, t, idx + 1) + end + + # Linear interpolation + return A.p.u_tilde[2 * idx] + A.p.slope[idx + 1] * (t - A.p.t_tilde[2 * idx]) +end + # Quadratic Interpolation function _interpolate(A::QuadraticInterpolation, t::Number, iguess) idx = get_idx(A, t, iguess) diff --git a/src/interpolation_utils.jl b/src/interpolation_utils.jl index ff16d56f..1499489d 100644 --- a/src/interpolation_utils.jl +++ b/src/interpolation_utils.jl @@ -415,3 +415,73 @@ function get_transition_ts(A::SmoothedConstantInterpolation) end get_transition_ts(A::AbstractInterpolation) = A.t + +# Compute the spline parametrization value s ∈ [0,1] +function S(A::SmoothedLinearInterpolation, t, idx) + (; Δt, ΔΔt, degenerate_ΔΔt, t_tilde, λ) = A.p + Δtᵢ = Δt[idx] + ΔΔtᵢ = ΔΔt[idx] + tdiff = t - t_tilde[2 * idx - 1] + @assert tdiff >= 0 + + s = if degenerate_ΔΔt[idx] + # Degenerate case Δtᵢ₊₁ ≈ Δtᵢ + tdiff / (λ * Δtᵢ) + else + (-Δtᵢ + sqrt(Δtᵢ^2 + 2 * ΔΔtᵢ * tdiff / λ)) / ΔΔtᵢ + end + ε = 1e-5 + @assert -ε<=s<=1+ε "s should be in [0,1], got $s." + return s +end + +# Compute the value of a spline section of SmoothedLinearInterpolation +# as a function of t +function U(A::SmoothedLinearInterpolation, t, idx) + s = S(A, t, idx) + return U_s(A, s, idx) +end + +# Compute the value of a spline section of SmoothedLinearInterpolation +# as a function of the spline parametrization value s +function U_s(A::SmoothedLinearInterpolation, s, idx) + (; Δu, ΔΔu, u_tilde, λ) = A.p + Δuᵢ = Δu[idx] + ΔΔuᵢ = ΔΔu[idx] + return λ / 2 * ΔΔuᵢ * s^2 + λ * Δuᵢ * s + u_tilde[2 * idx - 1] +end + +# Compute the derivative of a spline section of SmoothedLinearInterpolation +# as a function of t +function U_deriv(A::SmoothedLinearInterpolation, t, idx) + s = S(A, t, idx) + s_deriv = S_deriv(A, t, idx) + return U_s_deriv(A, s, idx) * s_deriv +end + +# Compute the derivative of the spline parametrization value s +# with respect to t +function S_deriv(A::SmoothedLinearInterpolation, t, idx) + (; Δt, ΔΔt, degenerate_ΔΔt, t_tilde, λ) = A.p + Δtᵢ = Δt[idx] + ΔΔtᵢ = ΔΔt[idx] + tdiff = t - t_tilde[2 * idx - 1] + @assert tdiff >= 0 + + if degenerate_ΔΔt[idx] + # Degenerate case Δtᵢ₊₁ ≈ Δtᵢ + s_deriv = 1 / (λ * Δtᵢ) + else + s_deriv = 1 / (λ * sqrt(Δtᵢ^2 + 2 * ΔΔtᵢ * tdiff / λ)) + end + return s_deriv +end + +# Compute the derivative of a spline section of SmoothedLinearInterpolation +# as a function of the spline parametrization value s +function U_s_deriv(A::AbstractInterpolation, s, idx) + (; Δu, ΔΔu, λ) = A.p + Δuᵢ = Δu[idx] + ΔΔuᵢ = ΔΔu[idx] + return λ * ΔΔuᵢ * s + λ * Δuᵢ +end diff --git a/src/parameter_caches.jl b/src/parameter_caches.jl index 31206829..bac057d5 100644 --- a/src/parameter_caches.jl +++ b/src/parameter_caches.jl @@ -32,6 +32,48 @@ function linear_interpolation_parameters(u::AbstractArray{T, N}, t, idx) where { return slope end +struct SmoothedLinearParameterCache{uType, tType, λType} + Δu::uType + Δt::tType + ΔΔu::uType + ΔΔt::tType + u_tilde::uType + t_tilde::tType + slope::uType + # Whether ΔΔt is sufficiently close to 0 + degenerate_ΔΔt::Vector{Bool} + λ::λType +end + +function get_spline_ends(u, Δu, λ) + u_tilde = zeros(2 * length(u)) + u_tilde[1] = u[1] + u_tilde[2:2:(end - 1)] = u[1:(end - 1)] .+ (λ / 2) .* Δu[2:(end - 1)] + u_tilde[3:2:end] = u[2:end] .- (λ / 2) .* Δu[2:(end - 1)] + u_tilde[end] = u[end] + return u_tilde +end + +function SmoothedLinearParameterCache(u::AbstractVector, t, λ, cache_parameters) + @assert cache_parameters "Parameter caching is mandatory for SmoothedLinearInterpolation" + + Δu = diff(u) + Δt = diff(t) + pushfirst!(Δt, last(Δt)) + push!(Δt, first(Δt)) + pushfirst!(Δu, last(Δu)) + push!(Δu, first(Δu)) + ΔΔu = diff(Δu) + ΔΔt = diff(Δt) + u_tilde = get_spline_ends(u, Δu, λ) + t_tilde = get_spline_ends(t, Δt, λ) + slope = Δu ./ Δt + degenerate_ΔΔt = collect(isapprox.(ΔΔt, 0, atol = 1e-5)) + + return SmoothedLinearParameterCache( + Δu, Δt, ΔΔu, ΔΔt, u_tilde, t_tilde, slope, degenerate_ΔΔt, λ) +end + struct SmoothedConstantParameterCache{dType, cType} d::dType c::cType diff --git a/test/derivative_tests.jl b/test/derivative_tests.jl index 081de086..bd27de21 100644 --- a/test/derivative_tests.jl +++ b/test/derivative_tests.jl @@ -35,7 +35,8 @@ function test_derivatives(method; args = [], kwargs = [], name::String) # Interpolation transition points for _t in t[2:(end - 1)] - if func isa Union{SmoothedConstantInterpolation, BSplineInterpolation, BSplineApprox} + if func isa + Union{SmoothedConstantInterpolation, BSplineInterpolation, BSplineApprox} # TODO fix interpolations continue else @@ -106,6 +107,16 @@ end LinearInterpolation; args = [u, t], name = "Linear Interpolation with two points") end +@testset "Smoothed Linear Interpolation" begin + u = [3.8, 6.7, 1.8, 2.3] + t = [1.05, 2.6, 6.7, 8.9] + + for λ in [0.0, 0.25, 0.5] + test_derivatives( + SmoothedLinearInterpolation, args = [u, t], name = "Smoothed Linear Interpolation (λ = $λ)") + end +end + @testset "Quadratic Interpolation" begin u = [1.0, 4.0, 9.0, 16.0] t = [1.0, 2.0, 3.0, 4.0] diff --git a/test/interpolation_tests.jl b/test/interpolation_tests.jl index 3677da4f..7d05393c 100644 --- a/test/interpolation_tests.jl +++ b/test/interpolation_tests.jl @@ -239,6 +239,20 @@ end @test_throws DataInterpolations.LeftExtrapolationError A([-1.0, 11.0]) end +@testset "Smoothed Linear Interpolation" begin + test_interpolation_type(SmoothedLinearInterpolation) + + u = [3.8, 6.7, 1.8, 2.3] + t = [1.0, 2.6, 6.7, 8.9] + A = SmoothedLinearInterpolation(u, t) + + @test A.p.t_tilde ≈ [1.0, 1.2, 2.4, 3.1125, 6.1875, 6.975, 8.625, 8.9] + @test A.p.u_tilde ≈ [3.8, 4.1625, 6.3375, 6.0875, 2.4125, 1.8625, 2.2375, 2.3] + + t_eval = A.p.t_tilde[1:(end - 1)] + diff(A.p.t_tilde) / 2 + @test A(t_eval)≈[3.98125, 5.25, 6.41932, 4.25, 2.01298, 2.05, 2.26875] rtol=1e-4 +end + @testset "Quadratic Interpolation" begin test_interpolation_type(QuadraticInterpolation) diff --git a/test/sparseconnectivitytracer_tests.jl b/test/sparseconnectivitytracer_tests.jl index fc6a7307..a35155bd 100644 --- a/test/sparseconnectivitytracer_tests.jl +++ b/test/sparseconnectivitytracer_tests.jl @@ -10,9 +10,7 @@ using Test myprimal(x) = x myprimal(d::Dual) = primal(d) -#===========# -# Test data # -#===========# +# Test data t = [0.0, 1.0, 2.5, 4.0, 6.0]; t_scalar = 2.0 @@ -24,9 +22,7 @@ u = sin.(t) # vector du = cos.(t) ddu = -sin.(t) -#==================# # Test definitions # -#==================# struct InterpolationTest{N, I <: AbstractInterpolation} # N = output dim. of interpolation interp::I @@ -34,18 +30,16 @@ struct InterpolationTest{N, I <: AbstractInterpolation} # N = output dim. of int is_der2_zero::Bool end function InterpolationTest( - interp::I; is_der1_zero = false, is_der2_zero = false - ) where {T, I <: AbstractInterpolation{T}} + interp::I; is_der1_zero = false, + is_der2_zero = false +) where {T, I <: AbstractInterpolation{T}} sz = output_size(interp) N = isempty(sz) ? 1 : only(sz) return InterpolationTest{N, I}(interp, is_der1_zero, is_der2_zero) end testname(t::InterpolationTest{N}) where {N} = "$N-dim $(typeof(t.interp))" -#================# -# Jacobian Tests # -#================# - +# Jacobian Tests function test_jacobian(t::InterpolationTest) return @testset "Jacobian" begin for input in test_inputs @@ -116,9 +110,7 @@ function test_jacobian(t::InterpolationTest{N}, input::AbstractVector) where {N} end end -#===============# -# Hessian Tests # -#===============# +# Hessian Tests function test_hessian(t::InterpolationTest) return @testset "Hessian" begin @@ -213,9 +205,9 @@ function test_output(t::InterpolationTest) @test s_tracer == s_ref end @testset "$T" for T in ( - Dual{eltype(input), DEFAULT_GRADIENT_TRACER}, - Dual{eltype(input), DEFAULT_HESSIAN_TRACER}, - ) + Dual{eltype(input), DEFAULT_GRADIENT_TRACER}, + Dual{eltype(input), DEFAULT_HESSIAN_TRACER} + ) t_dual = trace_input(T, input) out_dual = t.interp(t_dual) s_dual = size(out_dual) @@ -226,27 +218,25 @@ function test_output(t::InterpolationTest) end end -#===========# -# Run tests # -#===========# +# Run tests @testset "1D Interpolations" begin - @testset "$(testname(t))" for t in ( - InterpolationTest( - ConstantInterpolation(u, t); is_der1_zero = true, is_der2_zero = true - ), - InterpolationTest(LinearInterpolation(u, t); is_der2_zero = true), - InterpolationTest(QuadraticInterpolation(u, t)), - InterpolationTest(LagrangeInterpolation(u, t)), - InterpolationTest(AkimaInterpolation(u, t)), - InterpolationTest(QuadraticSpline(u, t)), - InterpolationTest(CubicSpline(u, t)), - InterpolationTest(BSplineInterpolation(u, t, 3, :ArcLen, :Average)), - InterpolationTest(BSplineApprox(u, t, 3, 4, :ArcLen, :Average)), - InterpolationTest(PCHIPInterpolation(u, t)), - InterpolationTest(CubicHermiteSpline(du, u, t)), - InterpolationTest(QuinticHermiteSpline(ddu, du, u, t)), - ) + @testset "$(testname(t))" for t in (InterpolationTest(AkimaInterpolation(u, t)), + InterpolationTest(BSplineApprox(u, t, 3, 4, :ArcLen, :Average)), + InterpolationTest(BSplineInterpolation(u, t, 3, :ArcLen, :Average)), + InterpolationTest( + ConstantInterpolation(u, t); is_der1_zero = true, is_der2_zero = true + ), + InterpolationTest(CubicHermiteSpline(du, u, t)), + InterpolationTest(CubicSpline(u, t)), + InterpolationTest(LagrangeInterpolation(u, t)), + InterpolationTest(LinearInterpolation(u, t); is_der2_zero = true), + InterpolationTest(PCHIPInterpolation(u, t)), + InterpolationTest(QuadraticInterpolation(u, t)), + InterpolationTest(QuadraticSpline(u, t)), + InterpolationTest(QuinticHermiteSpline(ddu, du, u, t)), + InterpolationTest(SmoothedLinearInterpolation(u, t)) + ) test_jacobian(t) test_hessian(t) test_output(t) @@ -259,20 +249,13 @@ for N in (2, 5) @testset "$(N)D Interpolations" begin @testset "$(testname(t))" for t in ( - InterpolationTest( - ConstantInterpolation(um, t); is_der1_zero = true, is_der2_zero = true - ), - InterpolationTest(LinearInterpolation(um, t); is_der2_zero = true), - InterpolationTest(QuadraticInterpolation(um, t)), - InterpolationTest(LagrangeInterpolation(um, t)), - ## The following interpolations appear to not be supported on N dimensions as of DataInterpolations v6.2.0: - # InterpolationTest(AkimaInterpolation(um, t)), - # InterpolationTest(BSplineApprox(um, t, 3, 4, :ArcLen, :Average)), - # InterpolationTest(QuadraticSpline(um, t)), - # InterpolationTest(CubicSpline(um, t)), - # InterpolationTest(BSplineInterpolation(um, t, 3, :ArcLen, :Average)), - # InterpolationTest(PCHIPInterpolation(um, t)), - ) + InterpolationTest( + ConstantInterpolation(um, t); is_der1_zero = true, is_der2_zero = true + ), + InterpolationTest(LinearInterpolation(um, t); is_der2_zero = true), + InterpolationTest(QuadraticInterpolation(um, t)), + InterpolationTest(LagrangeInterpolation(um, t)) ## The following interpolations appear to not be supported on N dimensions as of DataInterpolations v6.2.0: # InterpolationTest(AkimaInterpolation(um, t)), # InterpolationTest(BSplineApprox(um, t, 3, 4, :ArcLen, :Average)), # InterpolationTest(QuadraticSpline(um, t)), # InterpolationTest(CubicSpline(um, t)), # InterpolationTest(BSplineInterpolation(um, t, 3, :ArcLen, :Average)), # InterpolationTest(PCHIPInterpolation(um, t)), + ) test_jacobian(t) test_hessian(t) test_output(t)