Main integration class for variable step size Runge-Kutta methods
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
type(stepsize_class), | private | :: | stepsize_method |
the method for varying the step size |
|||
real(kind=wp), | private, | dimension(:), allocatable | :: | rtol |
relative tolerance ( |
||
real(kind=wp), | private, | dimension(:), allocatable | :: | atol |
absolute tolerance ( |
||
integer, | private | :: | hinit_method | = | 1 |
if automatically computing the inital step size, which
method to use. 1 = |
|
integer, | private | :: | num_rejected_steps | = | 0 |
number of rejected steps |
|
real(kind=wp), | private | :: | last_accepted_step_size | = | zero |
the last accepted step size |
user-callable method to stop the integration
User-callable method to stop the integration.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(rk_class), | intent(inout) | :: | me |
get status code and message
Get the status of an integration.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(rk_class), | intent(in) | :: | me | |||
integer, | intent(out), | optional | :: | istatus |
status code ( |
|
character(len=:), | intent(out), | optional, | allocatable | :: | message |
status message |
Returns true if there was an error. Can use rk_class_status to get more info.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(rk_class), | intent(in) | :: | me |
Returns the properties of the method.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(rk_class), | intent(in) | :: | me |
properties of the method
the step routine for the rk method
rk step function for the variable-step methods.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(rk_variable_step_class), | intent(inout) | :: | me | |||
real(kind=wp), | intent(in) | :: | t |
initial time |
||
real(kind=wp), | intent(in), | dimension(me%n) | :: | x |
initial state vector |
|
real(kind=wp), | intent(in) | :: | h |
time step |
||
real(kind=wp), | intent(out), | dimension(me%n) | :: | xf |
final state vector |
|
real(kind=wp), | intent(out), | dimension(me%n) | :: | xerr |
truncation error estimate |
initialize the class (set n,f, and report)
Initialize the rk_variable_step_class.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(rk_variable_step_class), | intent(inout) | :: | me | |||
integer, | intent(in) | :: | n |
number of equations |
||
procedure(deriv_func) | :: | f |
derivative function |
|||
real(kind=wp), | intent(in), | optional, | dimension(:) | :: | rtol |
relative tolerance (if size=1,
then same tol used for all
equations). If not present, a default
of |
real(kind=wp), | intent(in), | optional, | dimension(:) | :: | atol |
absolute tolerance (if size=1,
then same tol used for all
equations). If not present, a default
of |
type(stepsize_class), | intent(in), | optional | :: | stepsize_method |
method for varying the step size |
|
integer, | intent(in), | optional | :: | hinit_method |
which method (1 or 2) to use for
automatic initial step size
computation.
1 = use |
|
procedure(report_func), | optional | :: | report |
for reporting the steps |
||
procedure(event_func), | optional | :: | g |
for stopping at an event |
||
logical, | intent(in), | optional | :: | stop_on_errors |
stop the program for any errors (default is False) |
|
integer, | intent(in), | optional | :: | max_number_of_steps |
max number of steps allowed |
|
integer, | intent(in), | optional | :: | report_rate |
how often to call report function:
|
|
class(root_method), | intent(in), | optional | :: | solver |
the root-finding method to use for even finding.
if not present, then |
Main integration routine for the rk_variable_step_class.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(rk_variable_step_class), | intent(inout) | :: | me | |||
real(kind=wp), | intent(in) | :: | t0 |
initial time |
||
real(kind=wp), | intent(in), | dimension(:) | :: | x0 |
initial state |
|
real(kind=wp), | intent(in) | :: | h |
initial abs(time step) |
||
real(kind=wp), | intent(in) | :: | tf |
final time |
||
real(kind=wp), | intent(out), | dimension(:) | :: | xf |
final state |
Event-finding integration routine for the rk_variable_step_class. Integrates until g(t,x)=0, or until t=tf (whichever happens first).
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(rk_variable_step_class), | intent(inout) | :: | me | |||
real(kind=wp), | intent(in) | :: | t0 |
initial time |
||
real(kind=wp), | intent(in), | dimension(me%n) | :: | x0 |
initial state |
|
real(kind=wp), | intent(in) | :: | h |
abs(time step) |
||
real(kind=wp), | intent(in) | :: | tmax |
max final time if event not located |
||
real(kind=wp), | intent(in) | :: | tol |
function tolerance for root finding |
||
real(kind=wp), | intent(out) | :: | tf |
actual final time reached |
||
real(kind=wp), | intent(out), | dimension(me%n) | :: | xf |
final state (at tf) |
|
real(kind=wp), | intent(out) | :: | gf |
g value at tf |
Return some info about the integration.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(rk_variable_step_class), | intent(in) | :: | me | |||
integer, | intent(out), | optional | :: | num_steps |
number of steps taken |
|
integer, | intent(out), | optional | :: | num_rejected_steps |
number of rejected steps |
|
real(kind=wp), | intent(out), | optional | :: | last_accepted_step_size |
the last accepted step size
|
for automatically computing the initial step size [this is from DDEABM]
Computes a starting step size to be used in solving initial value problems in ordinary differential equations.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(rk_variable_step_class), | intent(inout) | :: | me | |||
real(kind=wp), | intent(in) | :: | a |
the initial point of integration. |
||
real(kind=wp), | intent(in) | :: | b |
a value of the independent variable used to define
the direction of integration. a reasonable choice is to
set |
||
real(kind=wp), | intent(in), | dimension(me%n) | :: | y |
the vector of initial values of the |
|
real(kind=wp), | intent(in), | dimension(me%n) | :: | yprime |
the vector of derivatives of the |
|
real(kind=wp), | intent(in), | dimension(me%n) | :: | etol |
the vector of error tolerances corresponding to
the |
|
real(kind=wp), | intent(out) | :: | h |
appropriate starting step size to be attempted by the differential equation method. |
for automatically computing the initial step size [this is from DOP853]
computation of an initial step size guess
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(rk_variable_step_class), | intent(inout) | :: | me | |||
real(kind=wp), | intent(in) | :: | x | |||
real(kind=wp), | intent(in), | dimension(:) | :: | y |
dimension(n) |
|
real(kind=wp), | intent(in) | :: | posneg |
posneg = sign(1.0_wp,xend-x) |
||
real(kind=wp), | intent(in), | dimension(:) | :: | f0 |
dimension(n) |
|
real(kind=wp), | intent(in) | :: | hmax | |||
real(kind=wp), | intent(in), | dimension(:) | :: | atol | ||
real(kind=wp), | intent(in), | dimension(:) | :: | rtol |
Begin a rk_variable_step_class integration.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(rk_variable_step_class), | intent(inout) | :: | me |
Compute the initial step size.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(rk_variable_step_class), | intent(inout) | :: | me | |||
real(kind=wp), | intent(in) | :: | t0 |
initial time |
||
real(kind=wp), | intent(in) | :: | tf |
final time |
||
real(kind=wp), | dimension(me%n) | :: | x0 |
initial state |
||
real(kind=wp), | intent(in) | :: | h0 |
user-input initial step size (if zero, then one is computed) |
step size to use
returns p
, the order of the method
Returns the order of the RK method
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(rk_variable_step_class), | intent(in) | :: | me |
order of the method
type,extends(rk_class),abstract,public :: rk_variable_step_class !! Main integration class for variable step size Runge-Kutta methods private type(stepsize_class) :: stepsize_method !! the method for varying the step size real(wp),dimension(:),allocatable :: rtol !! relative tolerance (`size(n)`) real(wp),dimension(:),allocatable :: atol !! absolute tolerance (`size(n)`) integer :: hinit_method = 1 !! if automatically computing the inital step size, which !! method to use. 1 = `hstart`, 2 = `hinit`. integer :: num_rejected_steps = 0 !! number of rejected steps real(wp) :: last_accepted_step_size = zero !! the last accepted step size `dt` from the integration !! (positive or negative) contains private procedure(step_func_variable),deferred :: step !! the step routine for the rk method procedure,public :: initialize => initialize_variable_step !! initialize the class (set n,f, and report) procedure,public :: integrate => integrate_variable_step procedure,public :: integrate_to_event => integrate_to_event_variable_step procedure,public :: info => info_variable_step procedure :: hstart !! for automatically computing the initial step size [this is from DDEABM] procedure :: hinit !! for automatically computing the initial step size [this is from DOP853] procedure :: begin_integration => begin_integration_rk_variable_step_class procedure :: compute_initial_step procedure :: order !! returns `p`, the order of the method end type rk_variable_step_class