-
Notifications
You must be signed in to change notification settings - Fork 26
/
Copy pathmethods.jl
594 lines (528 loc) · 20.2 KB
/
methods.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
export FiniteDifferenceMethod, fdm, backward_fdm, forward_fdm, central_fdm, extrapolate_fdm
"""
FiniteDifferences.DEFAULT_CONDITION
The default value for the [condition number](https://en.wikipedia.org/wiki/Condition_number)
of finite difference method. The condition number specifies the amplification of the ∞-norm
when passed to the function's derivative.
"""
const DEFAULT_CONDITION = 10
"""
FiniteDifferences.DEFAULT_FACTOR
The default factor number. The factor number specifies the multiple to amplify the estimated
round-off errors by.
"""
const DEFAULT_FACTOR = 1
abstract type FiniteDifferenceMethod{P,Q} end
"""
UnadaptedFiniteDifferenceMethod{P,Q} <: FiniteDifferenceMethod{P,Q}
A finite difference method that estimates a `Q`th order derivative from `P` function
evaluations. This method does not dynamically adapt its step size.
# Fields
- `grid::SVector{P,Int}`: Multiples of the step size that the function will be evaluated at.
- `coefs::SVector{P,Float64}`: Coefficients corresponding to the grid functions that the
function evaluations will be weighted by.
- `coefs_neighbourhood::NTuple{3,SVector{P,Float64}}`: Sets of coefficients used for
estimating the magnitude of the derivative in a neighbourhood.
- `bound_estimator::FiniteDifferenceMethod`: A finite difference method that is tuned to
perform adaptation for this finite difference method.
- `condition::Float64`: Condition number. See See [`DEFAULT_CONDITION`](@ref).
- `factor::Float64`: Factor number. See See [`DEFAULT_FACTOR`](@ref).
- `max_range::Float64`: Maximum distance that a function is evaluated from the input at
which the derivative is estimated.
- `∇f_magnitude_mult::Float64`: Internally computed quantity.
- `f_error_mult::Float64`: Internally computed quantity.
"""
struct UnadaptedFiniteDifferenceMethod{P,Q} <: FiniteDifferenceMethod{P,Q}
grid::SVector{P,Int}
coefs::SVector{P,Float64}
coefs_neighbourhood::NTuple{3,SVector{P,Float64}}
condition::Float64
factor::Float64
max_range::Float64
∇f_magnitude_mult::Float64
f_error_mult::Float64
end
"""
AdaptedFiniteDifferenceMethod{
P, Q, E<:FiniteDifferenceMethod
} <: FiniteDifferenceMethod{P,Q}
A finite difference method that estimates a `Q`th order derivative from `P` function
evaluations.
This method dynamically adapts its step size. The adaptation works by explicitly estimating
the truncation error and round-off error, and choosing the step size to optimally balance
those. The truncation error is given by the magnitude of the `P`th order derivative, which
will be estimated with another finite difference method (`bound_estimator`). This finite
difference method, `bound_estimator`, will be tasked with estimating the `P`th order
derivative in a _neighbourhood_, not just at some `x`. To do this, it will use a careful
reweighting of the function evaluations to estimate the `P`th order derivative at, in the
case of a central method, `x - h`, `x`, and `x + h`, where `h` is the step size. The
coefficients for this estimate, the _neighbourhood estimate_, are given by the three sets of
coefficients in `bound_estimator.coefs_neighbourhood`. The round-off error is estimated by the
round-off error of the function evaluations performed by `bound_estimator`. The truncation
error is amplified by `condition`, and the round-off error is amplified by `factor`. The
quantities `∇f_magnitude_mult` and `f_error_mult` are precomputed quantities that facilitate
the step size adaptation procedure.
# Fields
- `grid::SVector{P,Int}`: Multiples of the step size that the function will be evaluated at.
- `coefs::SVector{P,Float64}`: Coefficients corresponding to the grid functions that the
function evaluations will be weighted by.
- `coefs_neighbourhood::NTuple{3,SVector{P,Float64}}`: Sets of coefficients used for
estimating the magnitude of the derivative in a neighbourhood.
- `condition::Float64`: Condition number. See See [`DEFAULT_CONDITION`](@ref).
- `factor::Float64`: Factor number. See See [`DEFAULT_FACTOR`](@ref).
- `max_range::Float64`: Maximum distance that a function is evaluated from the input at
which the derivative is estimated.
- `∇f_magnitude_mult::Float64`: Internally computed quantity.
- `f_error_mult::Float64`: Internally computed quantity.
- `bound_estimator::FiniteDifferenceMethod`: A finite difference method that is tuned to
perform adaptation for this finite difference method.
"""
struct AdaptedFiniteDifferenceMethod{
P, Q, E<:FiniteDifferenceMethod
} <: FiniteDifferenceMethod{P,Q}
grid::SVector{P,Int}
coefs::SVector{P,Float64}
coefs_neighbourhood::NTuple{3,SVector{P,Float64}}
condition::Float64
factor::Float64
max_range::Float64
∇f_magnitude_mult::Float64
f_error_mult::Float64
bound_estimator::E
end
"""
FiniteDifferenceMethod(
grid::AbstractVector{Int},
q::Int;
condition::Real=DEFAULT_CONDITION,
factor::Real=DEFAULT_FACTOR,
max_range::Real=Inf
)
Construct a finite difference method.
# Arguments
- `grid::Vector{Int}`: The grid. See [`AdaptedFiniteDifferenceMethod`](@ref) or
[`UnadaptedFiniteDifferenceMethod`](@ref).
- `q::Int`: Order of the derivative to estimate.
# Keywords
- `condition::Real`: Condition number. See [`DEFAULT_CONDITION`](@ref).
- `factor::Real`: Factor number. See [`DEFAULT_FACTOR`](@ref).
- `max_range::Real=Inf`: Maximum distance that a function is evaluated from the input at
which the derivative is estimated.
# Returns
- `FiniteDifferenceMethod`: Specified finite difference method.
"""
function FiniteDifferenceMethod(
grid::SVector{P,Int},
q::Int;
condition::Real=DEFAULT_CONDITION,
factor::Real=DEFAULT_FACTOR,
max_range::Real=Inf,
) where P
_check_p_q(P, q)
coefs, coefs_neighbourhood, ∇f_magnitude_mult, f_error_mult = _coefs_mults(grid, q)
return UnadaptedFiniteDifferenceMethod{P,q}(
grid,
coefs,
coefs_neighbourhood,
condition,
factor,
max_range,
∇f_magnitude_mult,
f_error_mult
)
end
function FiniteDifferenceMethod(grid::AbstractVector{Int}, q::Int; kw_args...)
return FiniteDifferenceMethod(SVector{length(grid)}(grid), q; kw_args...)
end
"""
(m::FiniteDifferenceMethod)(f, x::T) where T<:AbstractFloat
Estimate the derivative of `f` at `x` using the finite differencing method `m` and an
automatically determined step size.
# Arguments
- `f`: Function to estimate derivative of.
- `x::T`: Input to estimate derivative at.
# Returns
- Estimate of the derivative.
# Examples
```julia-repl
julia> fdm = central_fdm(5, 1)
FiniteDifferenceMethod:
order of method: 5
order of derivative: 1
grid: [-2, -1, 0, 1, 2]
coefficients: [0.08333333333333333, -0.6666666666666666, 0.0, 0.6666666666666666, -0.08333333333333333]
julia> fdm(sin, 1)
0.5403023058681607
julia> fdm(sin, 1) - cos(1) # Check the error.
2.098321516541546e-14
julia> FiniteDifferences.estimate_step(fdm, sin, 1.0) # Computes step size and estimates the error.
(0.001065235154086019, 1.9541865128909085e-13)
```
"""
# We loop over all concrete subtypes of `FiniteDifferenceMethod` for Julia v1.0 compatibility.
for T in (UnadaptedFiniteDifferenceMethod, AdaptedFiniteDifferenceMethod)
@eval begin
function (m::$T)(f::TF, x::Real) where TF
x = float(x) # Assume that converting to float is desired, if it isn't already.
step = first(estimate_step(m, f, x))
return m(f, x, step)
end
function (m::$T{P,0})(f::TF, x::Real) where {P,TF}
# The automatic step size calculation fails if `Q == 0`, so handle that edge
# case.
return f(x)
end
end
end
"""
(m::FiniteDifferenceMethod)(f, x::T, step::Real) where T<:AbstractFloat
Estimate the derivative of `f` at `x` using the finite differencing method `m` and a given
step size.
# Arguments
- `f`: Function to estimate derivative of.
- `x::T`: Input to estimate derivative at.
- `step::Real`: Step size.
# Returns
- Estimate of the derivative.
# Examples
```julia-repl
julia> fdm = central_fdm(5, 1)
FiniteDifferenceMethod:
order of method: 5
order of derivative: 1
grid: [-2, -1, 0, 1, 2]
coefficients: [0.08333333333333333, -0.6666666666666666, 0.0, 0.6666666666666666, -0.08333333333333333]
julia> fdm(sin, 1, 1e-3)
0.5403023058679624
julia> fdm(sin, 1, 1e-3) - cos(1) # Check the error.
-1.7741363933510002e-13
```
"""
# We loop over all concrete subtypes of `FiniteDifferenceMethod` for 1.0 compatibility.
for T in (UnadaptedFiniteDifferenceMethod, AdaptedFiniteDifferenceMethod)
@eval begin
function (m::$T{P,Q})(f::TF, x::Real, step::Real) where {P,Q,TF}
x = float(x) # Assume that converting to float is desired, if it isn't already.
fs = _eval_function(m, f, x, step)
return _compute_estimate(m, fs, x, step, m.coefs)
end
end
end
function _eval_function(
m::FiniteDifferenceMethod, f::TF, x::T, step::Real,
) where {TF,T<:AbstractFloat}
return f.(x .+ T(step) .* m.grid)
end
function _compute_estimate(
m::FiniteDifferenceMethod{P,Q},
fs::SVector{P,TF},
x::T,
step::Real,
coefs::SVector{P,Float64},
) where {P,Q,TF,T<:AbstractFloat}
# If we substitute `T.(coefs)` in the expression below, then allocations occur. We
# therefore perform the broadcasting first. See
# https://github.com/JuliaLang/julia/issues/39151.
_coefs = T.(coefs)
return sum(fs .* _coefs) ./ T(step)^Q
end
# Check the method and derivative orders for consistency.
function _check_p_q(p::Integer, q::Integer)
q >= 0 || throw(DomainError(q, "order of derivative (`q`) must be non-negative"))
q < p || throw(DomainError(
(q, p),
"order of the method (`p`) must be strictly greater than that of the derivative " *
"(`q`)",
))
end
function _coefs(grid, p, q)
# For high precision on the `\`, we use `Rational`, and to prevent overflows we use
# `BigInt`. At the end we go to `Float64` for fast floating point math, rather than
# rational math.
C = [Rational{BigInt}(g)^i for i in 0:(p - 1), g in grid]
x = zeros(Rational{BigInt}, p)
x[q + 1] = factorial(big(q))
return SVector{p}(Float64.(C \ x))
end
const _COEFFS_MULTS_CACHE = Dict{
Tuple{SVector,Integer}, # Keys: (grid, q)
# Values: (coefs, coefs_neighbourhood, ∇f_magnitude_mult, f_error_mult)
Tuple{SVector,Tuple{Vararg{SVector}},Float64,Float64}
}()
# Compute coefficients for the method and cache the result.
function _coefs_mults(grid::SVector{P, Int}, q::Integer) where P
return get!(_COEFFS_MULTS_CACHE, (grid, q)) do
# Compute coefficients for derivative estimate.
coefs = _coefs(grid, P, q)
# Compute coefficients for a neighbourhood estimate.
if all(grid .>= 0)
coefs_neighbourhood = (
_coefs(grid .- 2, P, q),
_coefs(grid .- 1, P, q),
_coefs(grid, P, q)
)
elseif all(grid .<= 0)
coefs_neighbourhood = (
_coefs(grid, P, q),
_coefs(grid .+ 1, P, q),
_coefs(grid .+ 2, P, q)
)
else
coefs_neighbourhood = (
_coefs(grid .- 1, P, q),
_coefs(grid, P, q),
_coefs(grid .+ 1, P, q)
)
end
# Compute multipliers.
∇f_magnitude_mult = sum(abs.(coefs .* grid .^ P)) / factorial(big(P))
f_error_mult = sum(abs.(coefs))
return coefs, coefs_neighbourhood, ∇f_magnitude_mult, f_error_mult
end
end
function Base.show(
io::IO,
m::MIME"text/plain",
x::FiniteDifferenceMethod{P, Q},
) where {P, Q}
@printf io "FiniteDifferenceMethod:\n"
@printf io " order of method: %d\n" P
@printf io " order of derivative: %d\n" Q
@printf io " grid: %s\n" x.grid
@printf io " coefficients: %s\n" x.coefs
end
"""
function estimate_step(
m::FiniteDifferenceMethod,
f,
x::T
) where T<:AbstractFloat
Estimate the step size for a finite difference method `m`. Also estimates the error of the
estimate of the derivative.
# Arguments
- `m::FiniteDifferenceMethod`: Finite difference method to estimate the step size for.
- `f`: Function to evaluate the derivative of.
- `x::T`: Point to estimate the derivative at.
# Returns
- `Tuple{<:AbstractFloat, <:AbstractFloat}`: Estimated step size and an estimate of the
error of the finite difference estimate. The error will be `NaN` if the method failed
to estimate the error.
"""
function estimate_step(
m::UnadaptedFiniteDifferenceMethod, f::TF, x::T,
) where {TF,T<:AbstractFloat}
step, acc = _compute_step_acc_default(m, x)
return _limit_step(m, x, step, acc)
end
function estimate_step(
m::AdaptedFiniteDifferenceMethod{P,Q}, f::TF, x::T,
) where {P,Q,TF,T<:AbstractFloat}
∇f_magnitude, f_magnitude = _estimate_magnitudes(m.bound_estimator, f, x)
if ∇f_magnitude == 0.0 || f_magnitude == 0.0
step, acc = _compute_step_acc_default(m, x)
else
step, acc = _compute_step_acc(m, ∇f_magnitude, eps(f_magnitude))
end
return _limit_step(m, x, step, acc)
end
function _estimate_magnitudes(
m::FiniteDifferenceMethod{P,Q}, f::TF, x::T,
) where {P,Q,TF,T<:AbstractFloat}
step = first(estimate_step(m, f, x))
fs = _eval_function(m, f, x, step)
# Estimate magnitude of `∇f` in a neighbourhood of `x`.
∇fs = SVector{3}(
_compute_estimate(m, fs, x, step, m.coefs_neighbourhood[1]),
_compute_estimate(m, fs, x, step, m.coefs_neighbourhood[2]),
_compute_estimate(m, fs, x, step, m.coefs_neighbourhood[3])
)
∇f_magnitude = maximum(maximum.(abs, ∇fs))
# Estimate magnitude of `f` in a neighbourhood of `x`.
f_magnitude = maximum(maximum.(abs, fs))
return ∇f_magnitude, f_magnitude
end
function _compute_step_acc_default(m::FiniteDifferenceMethod, x::T) where {T<:AbstractFloat}
# Compute a default step size using a heuristic and [`DEFAULT_CONDITION`](@ref).
return _compute_step_acc(m, m.condition, eps(T))
end
function _compute_step_acc(
m::FiniteDifferenceMethod{P,Q}, ∇f_magnitude::Real, f_error::Real,
) where {P,Q}
# Set the step size by minimising an upper bound on the error of the estimate.
C₁ = f_error * m.f_error_mult * m.factor
C₂ = ∇f_magnitude * m.∇f_magnitude_mult
step = (Q / (P - Q) * (C₁ / C₂))^(1 / P)
# Estimate the accuracy of the method.
acc = C₁ * step^(-Q) + C₂ * step^(P - Q)
return step, acc
end
function _limit_step(
m::FiniteDifferenceMethod, x::T, step::Real, acc::Real,
) where {T<:AbstractFloat}
# First, limit the step size based on the maximum range.
step_max = m.max_range / maximum(abs.(m.grid))
if step > step_max
step = step_max
acc = NaN
end
# Second, prevent very large step sizes, which can occur for high-order methods or
# slowly-varying functions.
step_default, _ = _compute_step_acc_default(m, x)
step_max_default = 1000step_default
if step > step_max_default
step = step_max_default
acc = NaN
end
return step, acc
end
for direction in [:forward, :central, :backward]
fdm_fun = Symbol(direction, "_fdm")
grid_fun = Symbol("_", direction, "_grid")
@eval begin
function $fdm_fun(
p::Int,
q::Int;
adapt::Int=1,
condition::Real=DEFAULT_CONDITION,
factor::Real=DEFAULT_FACTOR,
max_range::Real=Inf,
geom::Bool=false
)
_check_p_q(p, q)
grid = $grid_fun(p)
geom && (grid = _exponentiate_grid(grid))
coefs, coefs_nbhd, ∇f_magnitude_mult, f_error_mult = _coefs_mults(grid, q)
if adapt >= 1
bound_estimator = $fdm_fun(
# We need to increase the order by two to be able to estimate the
# magnitude of the derivative of `f` in a neighbourhood of `x`.
p + 2,
p;
adapt=adapt - 1,
condition=condition,
factor=factor,
max_range=max_range,
geom=geom
)
return AdaptedFiniteDifferenceMethod{p, q, typeof(bound_estimator)}(
grid,
coefs,
coefs_nbhd,
condition,
factor,
max_range,
∇f_magnitude_mult,
f_error_mult,
bound_estimator
)
else
return UnadaptedFiniteDifferenceMethod{p, q}(
grid,
coefs,
coefs_nbhd,
condition,
factor,
max_range,
∇f_magnitude_mult,
f_error_mult
)
end
end
@doc """
$($(Meta.quot(fdm_fun)))(
p::Int,
q::Int;
adapt::Int=1,
condition::Real=DEFAULT_CONDITION,
factor::Real=DEFAULT_FACTOR,
max_range::Real=Inf,
geom::Bool=false
)
Construct a finite difference method at a $($(Meta.quot(direction))) grid of `p` points.
# Arguments
- `p::Int`: Number of grid points.
- `q::Int`: Order of the derivative to estimate.
# Keywords
- `adapt::Int=1`: Use another finite difference method to estimate the magnitude of the
`p`th order derivative, which is important for the step size computation. Recurse
this procedure `adapt` times.
- `condition::Real`: Condition number. See [`DEFAULT_CONDITION`](@ref).
- `factor::Real`: Factor number. See [`DEFAULT_FACTOR`](@ref).
- `max_range::Real=Inf`: Maximum distance that a function is evaluated from the input at
which the derivative is estimated.
- `geom::Bool`: Use geometrically spaced points instead of linearly spaced points.
# Returns
- `FiniteDifferenceMethod`: The specified finite difference method.
""" $fdm_fun
end
end
_forward_grid(p::Int) = SVector{p}(0:(p - 1))
_backward_grid(p::Int) = SVector{p}((1 - p):0)
function _central_grid(p::Int)
if isodd(p)
return SVector{p}(div(1 - p, 2):div(p - 1, 2))
else
return SVector{p}((div(-p, 2):-1)..., (1:div(p, 2))...)
end
end
function _exponentiate_grid(grid::SVector{P,Int}, base::Int=3) where P
return sign.(grid) .* div.(base .^ abs.(grid), base)
end
function _is_symmetric(
vec::SVector{P};
centre_zero::Bool=false,
negate_half::Bool=false
) where P
half_sign = negate_half ? -1 : 1
if isodd(length(vec))
centre_zero && vec[end ÷ 2 + 1] != 0 && return false
return vec[1:end ÷ 2] == half_sign .* reverse(vec[(end ÷ 2 + 2):end])
else
return vec[1:end ÷ 2] == half_sign .* reverse(vec[(end ÷ 2 + 1):end])
end
end
function _is_symmetric(m::FiniteDifferenceMethod)
grid_symmetric = _is_symmetric(m.grid, centre_zero=true, negate_half=true)
coefs_symmetric =_is_symmetric(m.coefs, negate_half=true)
return grid_symmetric && coefs_symmetric
end
"""
extrapolate_fdm(
m::FiniteDifferenceMethod,
f,
x::Real,
initial_step::Real=10,
power::Int=1,
breaktol::Real=Inf,
kw_args...
)
Use Richardson extrapolation to extrapolate a finite difference method.
Takes further in keyword arguments for `Richardson.extrapolate`. This method
automatically sets `power = 2` if `m` is symmetric and `power = 1`. Moreover, it defaults
`breaktol = Inf`.
# Arguments
- `m::FiniteDifferenceMethod`: Finite difference method to estimate the step size for.
- `f`: Function to evaluate the derivative of.
- `x::Real`: Point to estimate the derivative at.
- `initial_step::Real=10`: Initial step size.
# Returns
- `Tuple{<:AbstractFloat, <:AbstractFloat}`: Estimate of the derivative and error.
"""
function extrapolate_fdm(
m::FiniteDifferenceMethod,
f,
x::Real,
initial_step::Real=10;
power::Int=1,
breaktol::Real=Inf,
kw_args...
)
(power == 1 && _is_symmetric(m)) && (power = 2)
return extrapolate(
step -> m(f, x, step),
float(initial_step);
power=power,
breaktol=breaktol,
kw_args...
)
end