-
Notifications
You must be signed in to change notification settings - Fork 112
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
Monotonic interpolations #238
Changes from 6 commits
c3d3d0b
5854513
28ad194
1bfd8d2
9d65102
9c687b4
c3763e3
89dec3b
a19eff7
514185a
a2144cf
e6a38cb
File filter
Filter by extension
Conversations
Jump to
Diff view
Diff view
There are no files selected for viewing
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,355 @@ | ||
|
||
#= | ||
Interpolate using cubic Hermite splines. The breakpoints in arrays xbp and ybp are assumed to be sorted. | ||
Evaluate the function in all points of the array xeval. | ||
Methods: | ||
"Linear" yuck | ||
"FiniteDifference" classic cubic interpolation, no tension parameter | ||
Finite difference can overshoot for non-monotonic data | ||
"Cardinal" cubic cardinal splines, uses tension parameter which must be between [0,1] | ||
cubin cardinal splines can overshoot for non-monotonic data | ||
(increasing tension decreases overshoot) | ||
"FritschCarlson" monotonic - tangents are first initialized, then adjusted if they are not monotonic | ||
can overshoot for non-monotonic data | ||
"FritschButland" monotonic - faster algorithm (only requires one pass) but somewhat higher apparent "tension" | ||
"Steffen" monotonic - also only one pass, results usually between FritschCarlson and FritschButland | ||
Sources: | ||
Fritsch & Carlson (1980), "Monotone Piecewise Cubic Interpolation", doi:10.1137/0717021. | ||
Fritsch & Butland (1984), "A Method for Constructing Local Monotone Piecewise Cubic Interpolants", doi:10.1137/0905021. | ||
Steffen (1990), "A Simple Method for Monotonic Interpolation in One Dimension", http://adsabs.harvard.edu/abs/1990A%26A...239..443S | ||
|
||
Implementation based on http://bl.ocks.org/niclasmattsson/7bceb05fba6c71c78d507adae3d29417 | ||
=# | ||
|
||
export | ||
LinearMonotonicInterpolation, | ||
FiniteDifferenceMonotonicInterpolation, | ||
CardinalMonotonicInterpolation, | ||
FritschCarlsonMonotonicInterpolation, | ||
FritschButlandMonotonicInterpolation, | ||
SteffenMonotonicInterpolation | ||
|
||
""" | ||
MonotonicInterpolationType | ||
|
||
Abstract class for all types of monotonic interpolation. | ||
""" | ||
abstract type MonotonicInterpolationType <: InterpolationType end | ||
|
||
""" | ||
LinearMonotonicInterpolation | ||
|
||
Simple linear interpolation. | ||
""" | ||
struct LinearMonotonicInterpolation <: MonotonicInterpolationType | ||
end | ||
|
||
""" | ||
FiniteDifferenceMonotonicInterpolation | ||
|
||
Classic cubic interpolation, no tension parameter. | ||
Finite difference can overshoot for non-monotonic data. | ||
""" | ||
struct FiniteDifferenceMonotonicInterpolation <: MonotonicInterpolationType | ||
end | ||
|
||
""" | ||
CardinalMonotonicInterpolation(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). | ||
""" | ||
struct CardinalMonotonicInterpolation{T<:Number} <: MonotonicInterpolationType | ||
tension :: T # must be in [0, 1] | ||
end | ||
|
||
""" | ||
FritschCarlsonMonotonicInterpolation | ||
|
||
Monotonic 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 | ||
""" | ||
struct FritschCarlsonMonotonicInterpolation <: MonotonicInterpolationType | ||
end | ||
|
||
""" | ||
FritschButlandMonotonicInterpolation | ||
|
||
Monotonic 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". | ||
""" | ||
struct FritschButlandMonotonicInterpolation <: MonotonicInterpolationType | ||
end | ||
|
||
""" | ||
SteffenMonotonicInterpolation | ||
|
||
Monotonic 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. | ||
""" | ||
struct SteffenMonotonicInterpolation <: MonotonicInterpolationType | ||
end | ||
|
||
""" | ||
MonotonicInterpolation | ||
|
||
Monotonic interpolation up to third order represented by type, knots and | ||
coefficients. | ||
""" | ||
struct MonotonicInterpolation{T, TCoeffs, Tel, Type<:MonotonicInterpolationType, | ||
TKnots<:AbstractVector{<:Number}, AType <: AbstractArray{Tel,1}} <: AbstractInterpolation{T,1,DimSpec{Type}} | ||
|
||
it::Type | ||
knots::TKnots | ||
A::AType # constant parts of piecewise polynomials | ||
m::Vector{TCoeffs} # coefficients of linear parts of piecewise polynomials | ||
c::Vector{TCoeffs} # coefficients of quadratic parts of piecewise polynomials | ||
d::Vector{TCoeffs} # coefficients of cubic parts of piecewise polynomials | ||
end | ||
|
||
|
||
size(A::MonotonicInterpolation) = size(A.knots) | ||
axes(A::MonotonicInterpolation) = axes(A.knots) | ||
|
||
function MonotonicInterpolation(::Type{TWeights}, it::IType, knots::TKnots, A::AbstractArray{Tel,1}, | ||
m::Vector{TCoeffs}, c::Vector{TCoeffs}, d::Vector{TCoeffs}) where {TWeights, TCoeffs, Tel, IType<:MonotonicInterpolationType, TKnots<:AbstractVector{<:Number}} | ||
|
||
isconcretetype(IType) || error("The spline type must be a leaf type (was $IType)") | ||
isconcretetype(TCoeffs) || warn("For performance reasons, consider using an array of a concrete type (eltype(A) == $(eltype(A)))") | ||
|
||
check_monotonic(knots, A) | ||
|
||
cZero = zero(TWeights) | ||
if isempty(A) | ||
T = Base.promote_op(*, typeof(cZero), eltype(A)) | ||
else | ||
T = typeof(cZero * first(A)) | ||
end | ||
|
||
MonotonicInterpolation{T, TCoeffs, Tel, IType, TKnots, typeof(A)}(it, knots, A, m, c, d) | ||
end | ||
|
||
function interpolate(::Type{TWeights}, ::Type{TCoeffs}, knots::TKnots, | ||
A::AbstractArray{Tel,1}, it::IT) where {TWeights,TCoeffs,Tel,TKnots<:AbstractVector{<:Number},IT<:MonotonicInterpolationType} | ||
|
||
check_monotonic(knots, A) | ||
|
||
# first we need to determine tangents (m) | ||
n = length(knots) | ||
m, Δ = calcTangents(TCoeffs, knots, A, it) | ||
c = Vector{TCoeffs}(undef, n-1) | ||
d = Vector{TCoeffs}(undef, n-1) | ||
for k ∈ 1:n-1 | ||
if IT == LinearMonotonicInterpolation | ||
c[k] = d[k] = zero(TCoeffs) | ||
else | ||
xdiff = knots[k+1] - knots[k] | ||
c[k] = (3*Δ[k] - 2*m[k] - m[k+1]) / xdiff | ||
d[k] = (m[k] + m[k+1] - 2*Δ[k]) / (xdiff * xdiff) | ||
end | ||
end | ||
|
||
MonotonicInterpolation(TWeights, it, knots, A, m, c, d) | ||
end | ||
|
||
function interpolate(knots::AbstractVector{<:Number}, A::AbstractArray{Tel,1}, | ||
it::IT) where {Tel,IT<:MonotonicInterpolationType} | ||
|
||
interpolate(tweight(A), tcoef(A), knots, A, it) | ||
end | ||
|
||
function (itp::MonotonicInterpolation)(x::Number) | ||
x >= itp.knots[1] || error("Given number $x is outside of interpolated range [$(itp.knots[1]), $(itp.knots[end])]") | ||
x <= itp.knots[end] || error("Given number $x is outside of interpolated range [$(itp.knots[1]), $(itp.knots[end])]") | ||
k = searchsortedfirst(itp.knots, x) | ||
if k > 1 | ||
k -= 1 | ||
end | ||
xdiff = x - itp.knots[k] | ||
return itp.A[k] + itp.m[k]*xdiff + itp.c[k]*xdiff*xdiff + itp.d[k]*xdiff*xdiff*xdiff | ||
end | ||
|
||
function derivative(itp::MonotonicInterpolation, x::Number) | ||
k = searchsortedfirst(itp.knots, x) | ||
if k > 1 | ||
k -= 1 | ||
end | ||
xdiff = x - itp.knots[k] | ||
return itp.m[k] + 2*itp.c[k]*xdiff + 3*itp.d[k]*xdiff*xdiff | ||
end | ||
|
||
@inline function check_monotonic(knots, A) | ||
axes(knots) == axes(A) || throw(DimensionMismatch("knot vector must have the same axes as the corresponding array")) | ||
issorted(knots) || error("knot-vector must be sorted in increasing order") | ||
end | ||
|
||
function calcTangents(::Type{TCoeffs}, x::AbstractVector{<:Number}, | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. These methods could do with some more documentation, that explains the mathematical reasoning behind the algorithm (or at least refers to resources that explain it). The Fritsch & Carlson version has good docs, but the other ones have none... There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. This code was taken almost straight from the implementation provided by niclasmattsson in discussion for issue #105 . I haven't really looked into the details except for the interpretation of coefficients |
||
y::AbstractVector{Tel}, method::LinearMonotonicInterpolation) where {TCoeffs, Tel} | ||
|
||
n = length(x) | ||
Δ = Vector{TCoeffs}(undef, n-1) | ||
m = Vector{TCoeffs}(undef, n) | ||
for k in 1:n-1 | ||
Δk = (y[k+1] - y[k]) / (x[k+1] - x[k]) | ||
Δ[k] = Δk | ||
m[k] = Δk | ||
end | ||
m[n] = Δ[n-1] | ||
return (m, Δ) | ||
end | ||
|
||
function calcTangents(::Type{TCoeffs}, x::AbstractVector{<:Number}, | ||
y::AbstractVector{Tel}, method::FiniteDifferenceMonotonicInterpolation) where {TCoeffs, Tel} | ||
|
||
n = length(x) | ||
Δ = Vector{TCoeffs}(undef, n-1) | ||
m = Vector{TCoeffs}(undef, n) | ||
for k in 1:n-1 | ||
Δk = (y[k+1] - y[k]) / (x[k+1] - x[k]) | ||
Δ[k] = Δk | ||
if k == 1 # left endpoint | ||
m[k] = Δk | ||
else | ||
m[k] = (Δ[k-1] + Δk) / 2 | ||
end | ||
end | ||
m[n] = Δ[n-1] | ||
return (m, Δ) | ||
end | ||
|
||
function calcTangents(::Type{TCoeffs}, x::AbstractVector{<:Number}, | ||
y::AbstractVector{Tel}, method::CardinalMonotonicInterpolation{T}) where {T, TCoeffs, Tel} | ||
|
||
n = length(x) | ||
Δ = Vector{TCoeffs}(undef, n-1) | ||
m = Vector{TCoeffs}(undef, n) | ||
for k in 1:n-1 | ||
Δk = (y[k+1] - y[k]) / (x[k+1] - x[k]) | ||
Δ[k] = Δk | ||
if k == 1 # left endpoint | ||
m[k] = Δk | ||
else | ||
m[k] = (oneunit(T) - method.tension) * (y[k+1] - y[k-1]) / (x[k+1] - x[k-1]) | ||
end | ||
end | ||
m[n] = Δ[n-1] | ||
return (m, Δ) | ||
end | ||
|
||
function calcTangents(::Type{TCoeffs}, x::AbstractVector{<:Number}, | ||
y::AbstractVector{Tel}, method::FritschCarlsonMonotonicInterpolation) where {TCoeffs, Tel} | ||
|
||
n = length(x) | ||
Δ = Vector{TCoeffs}(undef, n-1) | ||
m = Vector{TCoeffs}(undef, n) | ||
for k in 1:n-1 | ||
Δk = (y[k+1] - y[k]) / (x[k+1] - x[k]) | ||
Δ[k] = Δk | ||
if k == 1 # left endpoint | ||
m[k] = Δk | ||
else | ||
# If any consecutive secant lines change sign (i.e. curve changes direction), initialize the tangent to zero. | ||
# This is needed to make the interpolation monotonic. Otherwise set tangent to the average of the secants. | ||
if Δk <= zero(Δk) | ||
m[k] = zero(TCoeffs) | ||
else | ||
m[k] = Δ[k-1] * (Δ[k-1] + Δk) / 2.0 | ||
end | ||
end | ||
end | ||
m[n] = Δ[n-1] | ||
#= | ||
Fritsch & Carlson derived necessary and sufficient conditions for monotonicity in their 1980 paper. Splines will be | ||
monotonic if all tangents are in a certain region of the alpha-beta plane, with alpha and beta as defined below. | ||
A robust choice is to put alpha & beta within a circle around origo with radius 3. The FritschCarlson algorithm | ||
makes simple initial estimates of tangents and then does another pass over data points to move any outlier tangents | ||
into the monotonic region. FritschButland & Steffen algorithms make more elaborate first estimates of tangents that | ||
are guaranteed to lie in the monotonic region, so no second pass is necessary. | ||
=# | ||
|
||
# Second pass of FritschCarlson: adjust any non-monotonic tangents. | ||
for k in 1:n-1 | ||
Δk = Δ[k] | ||
if Δk == zero(TCoeffs) | ||
m[k] = zero(TCoeffs) | ||
m[k+1] = zero(TCoeffs) | ||
continue | ||
end | ||
α = m[k] / Δk | ||
β = m[k+1] / Δk | ||
τ = 3.0 / sqrt(α^2 + β^2) | ||
if τ < 1.0 # if we're outside the circle with radius 3 then move onto the circle | ||
m[k] = τ * α * Δk | ||
m[k+1] = τ * β * Δk | ||
end | ||
end | ||
return (m, Δ) | ||
end | ||
|
||
function calcTangents(::Type{TCoeffs}, x::AbstractVector{<:Number}, | ||
y :: AbstractVector{Tel}, method :: FritschButlandMonotonicInterpolation) where {TCoeffs, Tel} | ||
|
||
# based on Fritsch & Butland (1984), | ||
# "A Method for Constructing Local Monotone Piecewise Cubic Interpolants", | ||
# doi:10.1137/0905021. | ||
|
||
n = length(x) | ||
Δ = Vector{TCoeffs}(undef, n-1) | ||
m = Vector{TCoeffs}(undef, n) | ||
for k in 1:n-1 | ||
Δk = (y[k+1] - y[k]) / (x[k+1] - x[k]) | ||
Δ[k] = Δk | ||
if k == 1 # left endpoint | ||
m[k] = Δk | ||
elseif Δ[k-1] * Δk <= zero(TCoeffs) | ||
m[k] = zero(TCoeffs) | ||
else | ||
α = (1.0 + (x[k+1] - x[k]) / (x[k+1] - x[k-1])) / 3.0 | ||
m[k] = Δ[k-1] * Δk / (α*Δk + (1.0 - α)*Δ[k-1]) | ||
end | ||
end | ||
m[n] = Δ[n-1] | ||
return (m, Δ) | ||
end | ||
|
||
function calcTangents(::Type{TCoeffs}, x::AbstractVector{<:Number}, | ||
y::AbstractVector{Tel}, method::SteffenMonotonicInterpolation) where {TCoeffs, Tel} | ||
|
||
# Steffen (1990), | ||
# "A Simple Method for Monotonic Interpolation in One Dimension", | ||
# http://adsabs.harvard.edu/abs/1990A%26A...239..443S | ||
|
||
n = length(x) | ||
Δ = Vector{TCoeffs}(undef, n-1) | ||
m = Vector{TCoeffs}(undef, n) | ||
for k in 1:n-1 | ||
Δk = (y[k+1] - y[k]) / (x[k+1] - x[k]) | ||
Δ[k] = Δk | ||
if k == 1 # left endpoint | ||
m[k] = Δk | ||
else | ||
p = ((x[k+1] - x[k]) * Δ[k-1] + (x[k] - x[k-1]) * Δk) / (x[k+1] - x[k-1]) | ||
m[k] = (sign(Δ[k-1]) + sign(Δk)) * | ||
min(abs(Δ[k-1]), abs(Δk), 0.5*abs(p)) | ||
end | ||
end | ||
m[n] = Δ[n-1] | ||
return (m, Δ) | ||
end | ||
|
||
# How many non-NoInterp dimensions are there? | ||
count_interp_dims(::Type{<:MonotonicInterpolationType}) = 1 | ||
|
||
lbounds(itp::MonotonicInterpolation) = first(itp.knots) | ||
ubounds(itp::MonotonicInterpolation) = last(itp.knots) |
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Why is this function not named
gradient
(returning an array of length 1) orgradient1
(returning a scalar)?There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
I've used the convention from ForwardDiff.jl but
gradient1
seems like a more fitting choice. Should I write a new method for thegradient
function as well (receiving and returning one element arrays)?There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
The behavior of
gradient
should, in order to match how other interpolation schemes work, beI think (correct me if I'm wrong) that the monotonic implementation here only supports 1D interpolations, which means this would effectively be
gradient(itp, x) = [derivative(itp, x)]
.gradient1
is just for convenience when the user knows it's a 1D interpolation object, but it can also avoid the array allocation:gradient1(itp, x) = derivative(itp, x)
. In other words, if you renamederivative
togradient1
and then addgradient(itp, x) = [gradient1(itp, x)]
(with typed arguments) you'll be good.As of #240 there's a fallback defined for
gradient!
that just callsgradient
.There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
But
gradient
should return anSVector
rather thanVector
---there's a ton of overhead forVector
but essentially none forSVector
.There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
I've renamed
derivative
togradient1
and introduced new functions:gradient
,hessian
andhessian1
.hessian
andgradient
returnSVector
s.And yes, only 1D -> 1D monotonic interpolation is supported by this code.