Skip to content

Commit fa05ca6

Browse files
authored
Merge pull request #228 from tpapp/tp/do-not-require-Function
Do not require that function arguments are ::Function.
2 parents 083aab5 + 65f5333 commit fa05ca6

File tree

2 files changed

+30
-16
lines changed

2 files changed

+30
-16
lines changed

src/methods.jl

Lines changed: 16 additions & 16 deletions
Original file line numberDiff line numberDiff line change
@@ -153,13 +153,13 @@ function FiniteDifferenceMethod(grid::AbstractVector{Int}, q::Int; kw_args...)
153153
end
154154

155155
"""
156-
(m::FiniteDifferenceMethod)(f::Function, x::T) where T<:AbstractFloat
156+
(m::FiniteDifferenceMethod)(f, x::T) where T<:AbstractFloat
157157
158158
Estimate the derivative of `f` at `x` using the finite differencing method `m` and an
159159
automatically determined step size.
160160
161161
# Arguments
162-
- `f::Function`: Function to estimate derivative of.
162+
- `f`: Function to estimate derivative of.
163163
- `x::T`: Input to estimate derivative at.
164164
165165
# Returns
@@ -188,12 +188,12 @@ julia> FiniteDifferences.estimate_step(fdm, sin, 1.0) # Computes step size and
188188
# We loop over all concrete subtypes of `FiniteDifferenceMethod` for Julia v1.0 compatibility.
189189
for T in (UnadaptedFiniteDifferenceMethod, AdaptedFiniteDifferenceMethod)
190190
@eval begin
191-
function (m::$T)(f::TF, x::Real) where TF<:Function
191+
function (m::$T)(f::TF, x::Real) where TF
192192
x = float(x) # Assume that converting to float is desired, if it isn't already.
193193
step = first(estimate_step(m, f, x))
194194
return m(f, x, step)
195195
end
196-
function (m::$T{P,0})(f::TF, x::Real) where {P,TF<:Function}
196+
function (m::$T{P,0})(f::TF, x::Real) where {P,TF}
197197
# The automatic step size calculation fails if `Q == 0`, so handle that edge
198198
# case.
199199
return f(x)
@@ -202,13 +202,13 @@ for T in (UnadaptedFiniteDifferenceMethod, AdaptedFiniteDifferenceMethod)
202202
end
203203

204204
"""
205-
(m::FiniteDifferenceMethod)(f::Function, x::T, step::Real) where T<:AbstractFloat
205+
(m::FiniteDifferenceMethod)(f, x::T, step::Real) where T<:AbstractFloat
206206
207207
Estimate the derivative of `f` at `x` using the finite differencing method `m` and a given
208208
step size.
209209
210210
# Arguments
211-
- `f::Function`: Function to estimate derivative of.
211+
- `f`: Function to estimate derivative of.
212212
- `x::T`: Input to estimate derivative at.
213213
- `step::Real`: Step size.
214214
@@ -235,7 +235,7 @@ julia> fdm(sin, 1, 1e-3) - cos(1) # Check the error.
235235
# We loop over all concrete subtypes of `FiniteDifferenceMethod` for 1.0 compatibility.
236236
for T in (UnadaptedFiniteDifferenceMethod, AdaptedFiniteDifferenceMethod)
237237
@eval begin
238-
function (m::$T{P,Q})(f::TF, x::Real, step::Real) where {P,Q,TF<:Function}
238+
function (m::$T{P,Q})(f::TF, x::Real, step::Real) where {P,Q,TF}
239239
x = float(x) # Assume that converting to float is desired, if it isn't already.
240240
fs = _eval_function(m, f, x, step)
241241
return _compute_estimate(m, fs, x, step, m.coefs)
@@ -245,7 +245,7 @@ end
245245

246246
function _eval_function(
247247
m::FiniteDifferenceMethod, f::TF, x::T, step::Real,
248-
) where {TF<:Function,T<:AbstractFloat}
248+
) where {TF,T<:AbstractFloat}
249249
return f.(x .+ T(step) .* m.grid)
250250
end
251251

@@ -336,7 +336,7 @@ end
336336
"""
337337
function estimate_step(
338338
m::FiniteDifferenceMethod,
339-
f::Function,
339+
f,
340340
x::T
341341
) where T<:AbstractFloat
342342
@@ -345,7 +345,7 @@ estimate of the derivative.
345345
346346
# Arguments
347347
- `m::FiniteDifferenceMethod`: Finite difference method to estimate the step size for.
348-
- `f::Function`: Function to evaluate the derivative of.
348+
- `f`: Function to evaluate the derivative of.
349349
- `x::T`: Point to estimate the derivative at.
350350
351351
# Returns
@@ -355,13 +355,13 @@ estimate of the derivative.
355355
"""
356356
function estimate_step(
357357
m::UnadaptedFiniteDifferenceMethod, f::TF, x::T,
358-
) where {TF<:Function,T<:AbstractFloat}
358+
) where {TF,T<:AbstractFloat}
359359
step, acc = _compute_step_acc_default(m, x)
360360
return _limit_step(m, x, step, acc)
361361
end
362362
function estimate_step(
363363
m::AdaptedFiniteDifferenceMethod{P,Q}, f::TF, x::T,
364-
) where {P,Q,TF<:Function,T<:AbstractFloat}
364+
) where {P,Q,TF,T<:AbstractFloat}
365365
∇f_magnitude, f_magnitude = _estimate_magnitudes(m.bound_estimator, f, x)
366366
if ∇f_magnitude == 0.0 || f_magnitude == 0.0
367367
step, acc = _compute_step_acc_default(m, x)
@@ -373,7 +373,7 @@ end
373373

374374
function _estimate_magnitudes(
375375
m::FiniteDifferenceMethod{P,Q}, f::TF, x::T,
376-
) where {P,Q,TF<:Function,T<:AbstractFloat}
376+
) where {P,Q,TF,T<:AbstractFloat}
377377
step = first(estimate_step(m, f, x))
378378
fs = _eval_function(m, f, x, step)
379379
# Estimate magnitude of `∇f` in a neighbourhood of `x`.
@@ -551,7 +551,7 @@ end
551551
"""
552552
extrapolate_fdm(
553553
m::FiniteDifferenceMethod,
554-
f::Function,
554+
f,
555555
x::Real,
556556
initial_step::Real=10,
557557
power::Int=1,
@@ -567,7 +567,7 @@ automatically sets `power = 2` if `m` is symmetric and `power = 1`. Moreover, it
567567
568568
# Arguments
569569
- `m::FiniteDifferenceMethod`: Finite difference method to estimate the step size for.
570-
- `f::Function`: Function to evaluate the derivative of.
570+
- `f`: Function to evaluate the derivative of.
571571
- `x::Real`: Point to estimate the derivative at.
572572
- `initial_step::Real=10`: Initial step size.
573573
@@ -576,7 +576,7 @@ automatically sets `power = 2` if `m` is symmetric and `power = 1`. Moreover, it
576576
"""
577577
function extrapolate_fdm(
578578
m::FiniteDifferenceMethod,
579-
f::Function,
579+
f,
580580
x::Real,
581581
initial_step::Real=10;
582582
power::Int=1,

test/methods.jl

Lines changed: 14 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1,3 +1,10 @@
1+
"""
2+
Wrapper used in test set “do not require f::Function” below, moved outside so that it
3+
works on Julia 1.0.
4+
"""
5+
struct NotAFunction end # not <: Function on purpose, cf #224
6+
(::NotAFunction)(x) = abs2(x)
7+
18
@testset "Methods" begin
29
@testset "Correctness" begin
310
# Finite difference methods to test.
@@ -162,4 +169,11 @@
162169
end
163170
end
164171
end
172+
173+
@testset "do not require f::Function" begin
174+
x = 0.7
175+
for f in [forward_fdm, central_fdm, backward_fdm]
176+
f(5, 1)(NotAFunction(), x) 2 * x
177+
end
178+
end
165179
end

0 commit comments

Comments
 (0)