Skip to content

Commit 3d84dcb

Browse files
Merge pull request #1973 from ranocha/hr/FineRKN5
fix out-of-place version of FineRKN5
2 parents f569818 + 36c5fa6 commit 3d84dcb

File tree

2 files changed

+180
-9
lines changed

2 files changed

+180
-9
lines changed

src/perform_step/rkn_perform_step.jl

Lines changed: 8 additions & 9 deletions
Original file line numberDiff line numberDiff line change
@@ -118,7 +118,7 @@ end
118118

119119
f.f1(k.x[1], du, u, p, t + dt)
120120
f.f2(k.x[2], du, u, p, t + dt)
121-
integrator.stats.nf += 1
121+
integrator.stats.nf += 4
122122
integrator.stats.nf2 += 1
123123
end
124124

@@ -159,7 +159,7 @@ end
159159
abar76 * k6) # abar72 = 0
160160

161161
k7 = f.f1(kdu, ku, p, t + dt * c7)
162-
u = uprev + dt * (duprev + dt * (b1 * k1 + b3 * k3 + b4 * k4 + b5 * k5)) # no b6, b7
162+
u = uprev + dt * (duprev + dt * (b1 * k1 + b3 * k3 + b4 * k4 + b5 * k5)) # no b6, b7
163163
du = duprev + dt * (bbar1 * k1 + bbar3 * k3 + bbar4 * k4 + bbar5 * k5 + bbar6 * k6) # no b2, b7
164164

165165
integrator.u = ArrayPartition((du, u))
@@ -223,7 +223,6 @@ end
223223
@.. broadcast=false ku=uprev +
224224
dt * (c7 * duprev +
225225
dt * (a71 * k1 + a73 * k3 + a74 * k4 + a75 * k5)) # a72 = a76 = 0
226-
227226
@.. broadcast=false kdu=duprev +
228227
dt * (abar71 * k1 + abar73 * k3 + abar74 * k4 +
229228
abar75 * k5 + abar76 * k6) # abar72 = 0
@@ -242,12 +241,12 @@ end
242241
if integrator.opts.adaptive
243242
duhat, uhat = utilde.x
244243
dtsq = dt^2
245-
@.. broadcast=false uhat = dtsq *
246-
(btilde1 * k1 + btilde3 * k3 + btilde4 * k4 +
247-
btilde5 * k5)
248-
@.. broadcast=false duhat = dt *
249-
(bptilde1 * k1 + bptilde3 * k3 + bptilde4 * k4 +
250-
bptilde5 * k5 + bptilde6 * k6 + bptilde7 * k7)
244+
@.. broadcast=false uhat=dtsq *
245+
(btilde1 * k1 + btilde3 * k3 + btilde4 * k4 +
246+
btilde5 * k5)
247+
@.. broadcast=false duhat=dt *
248+
(bptilde1 * k1 + bptilde3 * k3 + bptilde4 * k4 +
249+
bptilde5 * k5 + bptilde6 * k6 + bptilde7 * k7)
251250

252251
calculate_residuals!(atmp, utilde, integrator.uprev, integrator.u,
253252
integrator.opts.abstol, integrator.opts.reltol,

test/algconvergence/partitioned_methods_tests.jl

Lines changed: 172 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -367,3 +367,175 @@ sol = solve(prob, ERKN5(), reltol = 1e-8)
367367
@test length(sol.u) < 34
368368
sol = solve(prob, ERKN7(), reltol = 1e-8)
369369
@test length(sol.u) < 38
370+
371+
# Compare in-place and out-of-place versions
372+
function damped_oscillator(du, u, p, t)
373+
return -u - 0.5 * du
374+
end
375+
function damped_oscillator!(ddu, du, u, p, t)
376+
@. ddu = -u - 0.5 * du
377+
return nothing
378+
end
379+
@testset "in-place vs. out-of-place" begin
380+
ode_i = SecondOrderODEProblem(damped_oscillator!,
381+
[0.0], [1.0],
382+
(0.0, 10.0))
383+
ode_o = SecondOrderODEProblem(damped_oscillator,
384+
[0.0], [1.0],
385+
(0.0, 10.0))
386+
387+
@testset "Nystrom4" begin
388+
alg = Nystrom4()
389+
dt = 0.5
390+
# fixed time step
391+
sol_i = solve(ode_i, alg, dt = dt)
392+
sol_o = solve(ode_o, alg, dt = dt)
393+
@test sol_i.t sol_o.t
394+
@test sol_i.u sol_o.u
395+
@test sol_i.destats.nf == sol_o.destats.nf
396+
@test sol_i.destats.nf2 == sol_o.destats.nf2
397+
@test sol_i.destats.naccept == sol_o.destats.naccept
398+
@test 19 <= sol_i.destats.naccept <= 21
399+
@test abs(sol_i.destats.nf - 4 * sol_i.destats.naccept) < 4
400+
end
401+
402+
@testset "FineRKN5" begin
403+
alg = FineRKN5()
404+
dt = 0.5
405+
# fixed time step
406+
sol_i = solve(ode_i, alg, adaptive = false, dt = dt)
407+
sol_o = solve(ode_o, alg, adaptive = false, dt = dt)
408+
@test sol_i.t sol_o.t
409+
@test sol_i.u sol_o.u
410+
@test sol_i.destats.nf == sol_o.destats.nf
411+
@test sol_i.destats.nf2 == sol_o.destats.nf2
412+
@test sol_i.destats.naccept == sol_o.destats.naccept
413+
@test 19 <= sol_i.destats.naccept <= 21
414+
@test abs(sol_i.destats.nf - 7 * sol_i.destats.naccept) < 4
415+
# adaptive time step
416+
sol_i = solve(ode_i, alg)
417+
sol_o = solve(ode_o, alg)
418+
@test_broken sol_i.t sol_o.t
419+
@test_broken sol_i.u sol_o.u
420+
end
421+
422+
@testset "DPRKN4" begin
423+
alg = DPRKN4()
424+
dt = 0.5
425+
# fixed time step
426+
sol_i = solve(ode_i, alg, adaptive = false, dt = dt)
427+
sol_o = solve(ode_o, alg, adaptive = false, dt = dt)
428+
@test sol_i.t sol_o.t
429+
@test sol_i.u sol_o.u
430+
@test sol_i.destats.nf == sol_o.destats.nf
431+
@test sol_i.destats.nf2 == sol_o.destats.nf2
432+
@test sol_i.destats.naccept == sol_o.destats.naccept
433+
@test 19 <= sol_i.destats.naccept <= 21
434+
@test abs(sol_i.destats.nf - 4 * sol_i.destats.naccept) < 4
435+
# adaptive time step
436+
sol_i = solve(ode_i, alg)
437+
sol_o = solve(ode_o, alg)
438+
@test sol_i.t sol_o.t
439+
@test sol_i.u sol_o.u
440+
end
441+
442+
@testset "DPRKN5" begin
443+
alg = DPRKN5()
444+
dt = 0.5
445+
# fixed time step
446+
sol_i = solve(ode_i, alg, adaptive = false, dt = dt)
447+
sol_o = solve(ode_o, alg, adaptive = false, dt = dt)
448+
@test sol_i.t sol_o.t
449+
@test sol_i.u sol_o.u
450+
@test sol_i.destats.nf == sol_o.destats.nf
451+
@test sol_i.destats.nf2 == sol_o.destats.nf2
452+
@test sol_i.destats.naccept == sol_o.destats.naccept
453+
@test 19 <= sol_i.destats.naccept <= 21
454+
@test abs(sol_i.destats.nf - 6 * sol_i.destats.naccept) < 4
455+
# adaptive time step
456+
sol_i = solve(ode_i, alg)
457+
sol_o = solve(ode_o, alg)
458+
@test sol_i.t sol_o.t
459+
@test sol_i.u sol_o.u
460+
end
461+
462+
@testset "DPRKN6" begin
463+
alg = DPRKN6()
464+
dt = 0.5
465+
# fixed time step
466+
sol_i = solve(ode_i, alg, adaptive = false, dt = dt)
467+
sol_o = solve(ode_o, alg, adaptive = false, dt = dt)
468+
@test sol_i.t sol_o.t
469+
@test_broken sol_i.u sol_o.u
470+
@test sol_i.destats.nf == sol_o.destats.nf
471+
@test sol_i.destats.nf2 == sol_o.destats.nf2
472+
@test sol_i.destats.naccept == sol_o.destats.naccept
473+
@test 19 <= sol_i.destats.naccept <= 21
474+
@test abs(sol_i.destats.nf - 6 * sol_i.destats.naccept) < 4
475+
# adaptive time step
476+
sol_i = solve(ode_i, alg)
477+
sol_o = solve(ode_o, alg)
478+
@test_broken sol_i.t sol_o.t
479+
@test_broken sol_i.u sol_o.u
480+
end
481+
482+
@testset "DPRKN6FM" begin
483+
alg = DPRKN6FM()
484+
dt = 0.5
485+
# fixed time step
486+
sol_i = solve(ode_i, alg, adaptive = false, dt = dt)
487+
sol_o = solve(ode_o, alg, adaptive = false, dt = dt)
488+
@test sol_i.t sol_o.t
489+
@test sol_i.u sol_o.u
490+
@test sol_i.destats.nf == sol_o.destats.nf
491+
@test sol_i.destats.nf2 == sol_o.destats.nf2
492+
@test sol_i.destats.naccept == sol_o.destats.naccept
493+
@test 19 <= sol_i.destats.naccept <= 21
494+
@test abs(sol_i.destats.nf - 6 * sol_i.destats.naccept) < 4
495+
# adaptive time step
496+
sol_i = solve(ode_i, alg)
497+
sol_o = solve(ode_o, alg)
498+
@test_broken sol_i.t sol_o.t
499+
@test_broken sol_i.u sol_o.u
500+
end
501+
502+
@testset "DPRKN8" begin
503+
alg = DPRKN8()
504+
dt = 0.5
505+
# fixed time step
506+
sol_i = solve(ode_i, alg, adaptive = false, dt = dt)
507+
sol_o = solve(ode_o, alg, adaptive = false, dt = dt)
508+
@test sol_i.t sol_o.t
509+
@test sol_i.u sol_o.u
510+
@test sol_i.destats.nf == sol_o.destats.nf
511+
@test sol_i.destats.nf2 == sol_o.destats.nf2
512+
@test sol_i.destats.naccept == sol_o.destats.naccept
513+
@test 19 <= sol_i.destats.naccept <= 21
514+
@test abs(sol_i.destats.nf - 9 * sol_i.destats.naccept) < 4
515+
# adaptive time step
516+
sol_i = solve(ode_i, alg)
517+
sol_o = solve(ode_o, alg)
518+
@test_broken sol_i.t sol_o.t
519+
@test_broken sol_i.u sol_o.u
520+
end
521+
522+
@testset "DPRKN12" begin
523+
alg = DPRKN12()
524+
dt = 0.5
525+
# fixed time step
526+
sol_i = solve(ode_i, alg, adaptive = false, dt = dt)
527+
sol_o = solve(ode_o, alg, adaptive = false, dt = dt)
528+
@test sol_i.t sol_o.t
529+
@test sol_i.u sol_o.u
530+
@test sol_i.destats.nf == sol_o.destats.nf
531+
@test sol_i.destats.nf2 == sol_o.destats.nf2
532+
@test sol_i.destats.naccept == sol_o.destats.naccept
533+
@test 19 <= sol_i.destats.naccept <= 21
534+
@test abs(sol_i.destats.nf - 17 * sol_i.destats.naccept) < 4
535+
# adaptive time step
536+
sol_i = solve(ode_i, alg)
537+
sol_o = solve(ode_o, alg)
538+
@test_broken sol_i.t sol_o.t
539+
@test_broken sol_i.u sol_o.u
540+
end
541+
end

0 commit comments

Comments
 (0)