Skip to content

Commit ccd0b77

Browse files
committed
Co-authored-by: Songchen Tan <i@tansongchen.com>
1 parent e451f40 commit ccd0b77

File tree

6 files changed

+130
-2
lines changed

6 files changed

+130
-2
lines changed

lib/OrdinaryDiffEqTaylorSeries/Project.toml

+6
Original file line numberDiff line numberDiff line change
@@ -4,11 +4,14 @@ authors = ["ParamThakkar123 <paramthakkar864@gmail.com>"]
44
version = "1.1.0"
55

66
[deps]
7+
ADTypes = "47edcb42-4c32-4615-8424-f2b9edc5f35b"
78
DiffEqBase = "2b5f629d-d688-5b77-993f-72d75c75574e"
89
FastBroadcast = "7034ab61-46d4-4ed7-9d0f-46aef9175898"
910
LinearAlgebra = "37e2e46d-f89d-539d-b4ee-838fcccc9c8e"
1011
MuladdMacro = "46d2c3a1-f734-5fdb-9937-b9b9aeba4221"
1112
OrdinaryDiffEqCore = "bbf590c4-e513-4bbe-9b18-05decba2e5d8"
13+
OrdinaryDiffEqDifferentiation = "4302a76b-040a-498a-8c04-15b101fed76b"
14+
OrdinaryDiffEqNonlinearSolve = "127b3ac7-2247-4354-8eb6-78cf4e7c58e8"
1215
PrecompileTools = "aea7be01-6a6a-4083-8856-8a6e6704d82a"
1316
Preferences = "21216c6a-2e73-6563-6e65-726566657250"
1417
RecursiveArrayTools = "731186ca-8d62-57ce-b412-fbd966d074cd"
@@ -20,12 +23,15 @@ TaylorDiff = "b36ab563-344f-407b-a36a-4f200bebf99c"
2023
TruncatedStacktraces = "781d530d-4396-4725-bb49-402e4bee1e77"
2124

2225
[compat]
26+
ADTypes = "1.14.0"
2327
DiffEqBase = "6.152.2"
2428
DiffEqDevTools = "2.44.4"
2529
FastBroadcast = "0.3.5"
2630
LinearAlgebra = "<0.0.1, 1"
2731
MuladdMacro = "0.2.4"
2832
OrdinaryDiffEqCore = "1.1"
33+
OrdinaryDiffEqDifferentiation = "1.6.0"
34+
OrdinaryDiffEqNonlinearSolve = "1.6.0"
2935
PrecompileTools = "1.2.1"
3036
Preferences = "1.4.3"
3137
Random = "<0.0.1, 1"

lib/OrdinaryDiffEqTaylorSeries/src/OrdinaryDiffEqTaylorSeries.jl

+9-2
Original file line numberDiff line numberDiff line change
@@ -2,7 +2,9 @@ module OrdinaryDiffEqTaylorSeries
22

33
import OrdinaryDiffEqCore: alg_order, get_current_adaptive_order, get_current_alg_order,
44
alg_stability_size, explicit_rk_docstring,
5-
OrdinaryDiffEqAdaptiveAlgorithm, OrdinaryDiffEqMutableCache,
5+
OrdinaryDiffEqAdaptiveAlgorithm,
6+
OrdinaryDiffEqNewtonAdaptiveAlgorithm,
7+
OrdinaryDiffEqNewtonAlgorithm,OrdinaryDiffEqMutableCache,
68
alg_cache,
79
OrdinaryDiffEqConstantCache, @fold, trivial_limiter!,
810
step_accept_controller!,
@@ -15,6 +17,11 @@ import OrdinaryDiffEqCore: alg_order, get_current_adaptive_order, get_current_al
1517
CompositeAlgorithm, _ode_addsteps!, copyat_or_push!,
1618
AutoAlgSwitch, get_fsalfirstlast,
1719
full_cache, DerivativeOrderNotPossibleError
20+
using OrdinaryDiffEqDifferentiation: UJacobianWrapper, dolinsolve
21+
using OrdinaryDiffEqNonlinearSolve: du_alias_or_new, markfirststage!, build_nlsolver,
22+
nlsolve!, nlsolvefail, isnewton, get_W, set_new_W!,
23+
NLNewton, COEFFICIENT_MULTISTEP
24+
import ADTypes: AutoForwardDiff
1825
import Static: False
1926
import MuladdMacro: @muladd
2027
import FastBroadcast: @..
@@ -59,6 +66,6 @@ PrecompileTools.@compile_workload begin
5966
solver_list = nothing
6067
end
6168

62-
export ExplicitTaylor2, ExplicitTaylor, ExplicitTaylorAdaptiveOrder
69+
export ExplicitTaylor2, ExplicitTaylor, ExplicitTaylorAdaptiveOrder, ImplicitTaylor2
6370

6471
end

lib/OrdinaryDiffEqTaylorSeries/src/TaylorSeries_caches.jl

+14
Original file line numberDiff line numberDiff line change
@@ -150,3 +150,17 @@ function alg_cache(
150150
ExplicitTaylorAdaptiveOrderConstantCache(
151151
alg.min_order, alg.max_order, current_order, jets)
152152
end
153+
154+
struct ImplicitTaylor2ConstantCache{N} <: OrdinaryDiffEqConstantCache
155+
nlsolver::N
156+
end
157+
158+
function alg_cache(alg::ImplicitTaylor2, u, rate_prototype, ::Type{uEltypeNoUnits},
159+
::Type{uBottomEltypeNoUnits},
160+
::Type{tTypeNoUnits}, uprev, uprev2, f, t, dt, reltol, p, calck,
161+
::Val{false}) where {uEltypeNoUnits, uBottomEltypeNoUnits, tTypeNoUnits}
162+
γ, c = 1, 1
163+
nlsolver = build_nlsolver(alg, u, uprev, p, t, dt, f, rate_prototype, uEltypeNoUnits,
164+
uBottomEltypeNoUnits, tTypeNoUnits, γ, c, Val(false))
165+
ImplicitTaylor2ConstantCache(nlsolver)
166+
end

lib/OrdinaryDiffEqTaylorSeries/src/TaylorSeries_perform_step.jl

+68
Original file line numberDiff line numberDiff line change
@@ -207,3 +207,71 @@ end
207207
end
208208
return nothing
209209
end
210+
211+
function initialize!(integrator, cache::ImplicitTaylor2ConstantCache)
212+
integrator.kshortsize = 2
213+
integrator.k = typeof(integrator.k)(undef, integrator.kshortsize)
214+
integrator.fsalfirst = integrator.f(integrator.uprev, integrator.p, integrator.t) # Pre-start fsal
215+
OrdinaryDiffEqCore.increment_nf!(integrator.stats, 1)
216+
217+
# Avoid undefined entries if k is an array of arrays
218+
integrator.fsallast = zero(integrator.fsalfirst)
219+
integrator.k[1] = integrator.fsalfirst
220+
integrator.k[2] = integrator.fsallast
221+
end
222+
223+
@muladd function perform_step!(integrator, cache::ImplicitTaylor2ConstantCache,
224+
repeat_step = false)
225+
@unpack t, dt, uprev, u, f, p = integrator
226+
nlsolver = cache.nlsolver
227+
alg = unwrap_alg(integrator, true)
228+
markfirststage!(nlsolver)
229+
230+
# initial guess
231+
if alg.extrapolant == :linear
232+
nlsolver.z = dt * integrator.fsalfirst
233+
else # :constant
234+
nlsolver.z = zero(u)
235+
end
236+
237+
nlsolver.tmp = uprev
238+
nlsolver.γ = 1
239+
z = nlsolve!(nlsolver, integrator, cache, repeat_step)
240+
nlsolvefail(nlsolver) && return
241+
u = nlsolver.tmp + z
242+
243+
if integrator.opts.adaptive && integrator.success_iter > 0
244+
# local truncation error (LTE) bound by dt^2/2*max|y''(t)|
245+
# use 2nd divided differences (DD) a la SPICE and Shampine
246+
247+
# TODO: check numerical stability
248+
uprev2 = integrator.uprev2
249+
tprev = integrator.tprev
250+
251+
dt1 = dt * (t + dt - tprev)
252+
dt2 = (t - tprev) * (t + dt - tprev)
253+
c = 7 / 12 # default correction factor in SPICE (LTE overestimated by DD)
254+
r = c * dt^2 # by mean value theorem 2nd DD equals y''(s)/2 for some s
255+
256+
tmp = r *
257+
integrator.opts.internalnorm.((u - uprev) / dt1 - (uprev - uprev2) / dt2, t)
258+
atmp = calculate_residuals(tmp, uprev, u, integrator.opts.abstol,
259+
integrator.opts.reltol, integrator.opts.internalnorm, t)
260+
integrator.EEst = integrator.opts.internalnorm(atmp, t)
261+
else
262+
integrator.EEst = 1
263+
end
264+
265+
integrator.fsallast = f(u, p, t + dt)
266+
267+
if integrator.opts.adaptive && integrator.differential_vars !== nothing
268+
atmp = @. ifelse(!integrator.differential_vars, integrator.fsallast, false) ./
269+
integrator.opts.abstol
270+
integrator.EEst += integrator.opts.internalnorm(atmp, t)
271+
end
272+
273+
OrdinaryDiffEqCore.increment_nf!(integrator.stats, 1)
274+
integrator.k[1] = integrator.fsalfirst
275+
integrator.k[2] = integrator.fsallast
276+
integrator.u = u
277+
end

lib/OrdinaryDiffEqTaylorSeries/src/algorithms.jl

+28
Original file line numberDiff line numberDiff line change
@@ -35,3 +35,31 @@ Base.@kwdef struct ExplicitTaylorAdaptiveOrder{StageLimiter, StepLimiter, Thread
3535
step_limiter!::StepLimiter = trivial_limiter!
3636
thread::Thread = False()
3737
end
38+
39+
@doc explicit_rk_docstring(
40+
"An implicit Taylor series method.",
41+
"ImplicitTaylor")
42+
struct ImplicitTaylor2{CS, AD, F, F2, P, FDT, ST, CJ, StepLimiter} <:
43+
OrdinaryDiffEqNewtonAdaptiveAlgorithm{CS, AD, FDT, ST, CJ}
44+
linsolve::F
45+
nlsolve::F2
46+
precs::P
47+
extrapolant::Symbol
48+
controller::Symbol
49+
step_limiter!::StepLimiter
50+
autodiff::AD
51+
end
52+
53+
function ImplicitTaylor2(; chunk_size = Val{0}(), autodiff = AutoForwardDiff(),
54+
standardtag = Val{true}(), concrete_jac = nothing,
55+
diff_type = Val{:forward}(),
56+
linsolve = nothing, precs = DEFAULT_PRECS, nlsolve = NLNewton(),
57+
extrapolant = :constant,
58+
controller = :PI, step_limiter! = trivial_limiter!)
59+
AD_choice, chunk_size, diff_type = _process_AD_choice(autodiff, chunk_size, diff_type)
60+
61+
ImplicitTaylor2{_unwrap_val(chunk_size), typeof(AD_choice), typeof(linsolve),
62+
typeof(nlsolve), typeof(precs), diff_type, _unwrap_val(standardtag),
63+
_unwrap_val(concrete_jac), typeof(step_limiter!)}(linsolve,
64+
nlsolve, precs, extrapolant, controller, step_limiter!, AD_choice)
65+
end

lib/OrdinaryDiffEqTaylorSeries/test/runtests.jl

+5
Original file line numberDiff line numberDiff line change
@@ -33,3 +33,8 @@ end
3333
sol = solve(prob_ode_linear, ExplicitTaylorAdaptiveOrder())
3434
@test length(sol) < 20
3535
end
36+
37+
@testset "Implicit Taylor Tests" begin
38+
sol = solve(prob_ode_linear, ImplicitTaylor2())
39+
@test length(sol) < 100
40+
end

0 commit comments

Comments
 (0)