MRKISS 2025-09-10
A tiny library with zero dependencies that aims to make it easy to use & experiment with explicit Runge-Kutta methods.
Loading...
Searching...
No Matches
mrkiss_solvers_nt.f90 File Reference

RK Solvers. More...

Go to the source code of this file.

Data Types

interface  mrkiss_solvers_nt::deq_iface
 Type for ODE dydt functions. More...
interface  mrkiss_solvers_nt::stepp_iface
 Type step processing subroutine. More...
interface  mrkiss_solvers_nt::sdf_iface
 Type SDF function on a solution point. More...

Modules

module  mrkiss_solvers_nt
 Solvers.

Functions/Subroutines

One Step Solvers
subroutine, public mrkiss_solvers_nt::step_all (status, y_deltas, dy, deq, y, param, a, b, c, t_delta)
 Compute one step of a embedded Runge-Kutta method expressed as a Butcher Tableau.
subroutine, public mrkiss_solvers_nt::step_one (status, y_delta, dy, deq, y, param, a, b, c, t_delta)
 Compute one step of a Runge-Kutta method expressed as a Butcher Tableau.
subroutine, public mrkiss_solvers_nt::step_richardson (status, y_delta, dy, deq, y, param, a, b, c, p, t_delta)
 Compute one Richardson Extrapolation Step.
One Step Test Solvers
subroutine, public mrkiss_solvers_nt::step_rk4 (status, y_delta, dy, deq, y, param, t_delta)
 Compute one step of RK4 – the classical 4th order Runge-Kutta method (mrkiss_erk_kutta_4).
subroutine, public mrkiss_solvers_nt::step_rkf45 (status, y_deltas, dy, deq, y, param, t_delta)
 Compute one step of RKF45 (mrkiss_eerk_fehlberg_4_5).
subroutine, public mrkiss_solvers_nt::step_dp54 (status, y_deltas, dy, deq, y, param, t_delta)
 Compute one step of DP45 (mrkiss_eerk_dormand_prince_5_4)
Multistep Solvers
subroutine, public mrkiss_solvers_nt::fixed_t_steps (status, istats, solution, deq, y, param, a, b, c, p_o, max_pts_o, t_delta_o, t_end_o, t_max_o)
 Take multiple fixed time steps with a simple Runge-Kutta method.
subroutine, public mrkiss_solvers_nt::fixed_y_steps (status, istats, solution, deq, y, param, a, b, c, y_delta_len_targ, t_delta_max, t_delta_min_o, y_delta_len_tol_o, max_bisect_o, no_bisect_error_o, y_delta_len_idxs_o, max_pts_o, y_sol_len_max_o, t_max_o)
 Take multiple adaptive steps with constant \(\mathbf{\Delta{y}}\) length using a simple Runge-Kutta method.
subroutine, public mrkiss_solvers_nt::sloppy_fixed_y_steps (status, istats, solution, deq, y, param, a, b, c, y_delta_len_targ, t_delta_ini, t_max_o, t_delta_min_o, t_delta_max_o, y_delta_len_idxs_o, adj_short_o, max_pts_o, y_sol_len_max_o)
 Take multiple adaptive steps adjusting for the length of \(\mathbf{\Delta{y}}\) using a simple Runge-Kutta method.
subroutine, public mrkiss_solvers_nt::adaptive_steps (status, istats, solution, deq, y, param, a, b, c, p, t_max_o, t_end_o, t_delta_ini_o, t_delta_min_o, t_delta_max_o, t_delta_fac_min_o, t_delta_fac_max_o, t_delta_fac_fdg_o, error_tol_abs_o, error_tol_rel_o, max_pts_o, max_bisect_o, no_bisect_error_o, sdf_o, sdf_tol_o, stepp_o)
 Take multiple adaptive steps with an embedded Runge-Kutta method using relatively traditional step size controls.
Multistep Meta Solvers
subroutine, public mrkiss_solvers_nt::fixed_t_steps_between (status, istats, solution, deq, y, param, a, b, c, steps_per_pnt, p_o)
 Take a solution with \(t\) values pre-populated, and take steps_per_pnt Runge-Kutta steps between each \(t\) value.
subroutine, public mrkiss_solvers_nt::interpolate_solution (status, istats, solution, src_solution, deq, param, num_src_pts_o, linear_interp_o)
 Create an interpolated solution from a source solution.

Detailed Description

RK Solvers.

Author
Mitch Richling http://www.mitchr.me/
Keywords
runge kutta embedded butcher tableau
Standards
F2023
See also
https://github.com/richmit/MRKISS

Definition in file mrkiss_solvers_nt.f90.