Skip to content
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

Added Algorithm FineRKN4 #1975

Closed
wants to merge 7 commits into from
Closed

Conversation

HenryLangner
Copy link
Contributor

@HenryLangner HenryLangner commented Jun 19, 2023

Added the 4th order Nyström-Method presented by J. M. Fine in the article "Low Order Practical Runge-Kutta-Nyström Methods" (Also linked here: #677). Works with adaptive step-size. The Method passes the implemented convergence tests (harm. Oscillator):

using OrdinaryDiffEq, Test, RecursiveArrayTools, DiffEqDevTools, Statistics

u0 = fill(0.0, 2);
v0 = ones(2);
function f1_harmonic(dv, v, u, p, t)
    dv .= -u
end;
function f2_harmonic(du, v, u, p, t)
    du .= v
end;
function harmonic_analytic(y0, p, x)
    v0, u0 = y0.x
    ArrayPartition(-u0 * sin(x) + v0 * cos(x), u0 * cos(x) + v0 * sin(x))
end;
ff_harmonic = DynamicalODEFunction(f1_harmonic, f2_harmonic; analytic = harmonic_analytic);
prob = DynamicalODEProblem(ff_harmonic, v0, u0, (0.0, 5.0));

dts = 1.0 ./ 2.0 .^ (5:-1:0);

println("##### In Place #####")
sim = test_convergence(dts, prob, FineRKN4(), dense_errors = true);
println(@test sim.𝒪est[:l2]≈5 rtol=1e-1)
println(@test sim.𝒪est[:L2]≈4 rtol=1e-1)
println("with:")
println(sim.𝒪est)
println("##### Out of Place ##### ")

u0 = 0.0;
v0 = 1.0;
function f1_harmonic_nip(v, u, p, t)
    -u
end;
function f2_harmonic_nip(v, u, p, t)
    v
end;

ff_harmonic_nip = DynamicalODEFunction(f1_harmonic_nip, f2_harmonic_nip;
                                       analytic = harmonic_analytic);
prob = DynamicalODEProblem(ff_harmonic_nip, v0, u0, (0.0, 5.0));

csim = test_convergence(dts, prob, FineRKN4(), dense_errors = true);
println(@test sim.𝒪est[:l2]≈5 rtol=1e-1)
println(@test sim.𝒪est[:L2]≈4 rtol=1e-1)
println("with:")
println(sim.𝒪est)

Yields:
grafik

There is one more PR on the way which adds adaptive step-size to FineRKN5.
Once these two PRs are merged I will add convergence test for FineRKN4, FineRKN5 and Nystrom4 regarding the damped_oscillator introduced in #1973.

CC @ranocha

Removed all code related to adaptive-step-size feature for ```FineRKN5```.
This feature is beeing pushed in another PR.
@ranocha ranocha mentioned this pull request Jun 19, 2023
20 tasks
Copy link
Member

@ranocha ranocha left a comment

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Thanks! I just have a minor comment

src/perform_step/rkn_perform_step.jl Outdated Show resolved Hide resolved
@ranocha
Copy link
Member

ranocha commented Jun 19, 2023

Do you have some tests with adaptive step size control here?

@HenryLangner
Copy link
Contributor Author

Do you have some tests with adaptive step size control here?

I just added them!

Copy link
Member

@ranocha ranocha left a comment

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Thanks!

@HenryLangner HenryLangner deleted the FineRKN4 branch July 3, 2023 11:50
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Projects
None yet
Development

Successfully merging this pull request may close these issues.

2 participants