Modern Fortran implementation of the DDEABM Adams-Bashforth algorithm.
The original SLATEC code is a public domain work of the US Government. The modifications are Copyright (c) 2014-2018, Jacob Williams.
Note
This module depends on roots-fortran.
When compiled with FPM, this will automatically be downloaded (see fpm.toml
)
Note
The default real kind (wp
) can be
changed using optional preprocessor flags.
This library was built with real kind:
real(kind=real64)
[8 bytes]
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer, | public, | parameter | :: | ddeabm_rk | = | real64 |
real kind used by this module [8 bytes] |
integer, | private, | parameter | :: | wp | = | ddeabm_rk |
local copy of |
real(kind=wp), | private, | parameter | :: | d1mach2 | = | huge(1.0_wp) |
the largest magnitude |
real(kind=wp), | private, | parameter | :: | d1mach4 | = | epsilon(1.0_wp) |
the largest relative spacing |
function to determine if a root is bracketed
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_with_event_class), | intent(inout) | :: | me | |||
real(kind=wp), | intent(in) | :: | t1 |
first time point |
||
real(kind=wp), | intent(in) | :: | t2 |
second time point |
||
real(kind=wp), | intent(in), | dimension(:) | :: | x1 |
state at t1 |
|
real(kind=wp), | intent(in), | dimension(:) | :: | x2 |
state at t1 |
|
real(kind=wp), | intent(in) | :: | g1 |
event function value at t1 |
||
real(kind=wp), | intent(in) | :: | g2 |
event function value at t2 |
if the root is bracketed
function to determine if a root is bracketed (vector event version)
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_with_event_class_vec), | intent(inout) | :: | me | |||
integer, | intent(in), | optional | :: | ig |
the event function number |
|
real(kind=wp), | intent(in) | :: | t1 |
first time point |
||
real(kind=wp), | intent(in) | :: | t2 |
second time point |
||
real(kind=wp), | intent(in), | dimension(:) | :: | x1 |
state at t1 |
|
real(kind=wp), | intent(in), | dimension(:) | :: | x2 |
state at t1 |
|
real(kind=wp), | intent(in) | :: | g1 |
event function ig value at t1 |
||
real(kind=wp), | intent(in) | :: | g2 |
event function ig value at t2 |
if the root is bracketed
Interface to the ddeabm_class derivative function.
Known as df
in the documentation. provided by the user
to define the system of first order differential equations
which is to be solved. for the given values of t and the
vector x()=(x(1),x(2),...,x(neq)), the subroutine must
evaluate the neq components of the system of differential
equations dx/dt=df(t,x) and store the derivatives in the
array xdot(), that is, xdot(i) = dx(i)/dt for
equations i=1,...,neq.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_class), | intent(inout) | :: | me | |||
real(kind=wp), | intent(in) | :: | t |
time |
||
real(kind=wp), | intent(in), | dimension(:) | :: | x |
state |
|
real(kind=wp), | intent(out), | dimension(:) | :: | xdot |
derivative of state () |
Interface to the step reporting function.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_class), | intent(inout) | :: | me | |||
real(kind=wp), | intent(in) | :: | t |
time |
||
real(kind=wp), | intent(in), | dimension(:) | :: | x |
state |
Interface to the ddeabm_with_event_class scalar event function.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_with_event_class), | intent(inout) | :: | me | |||
real(kind=wp), | intent(in) | :: | t |
time |
||
real(kind=wp), | intent(in), | dimension(:) | :: | x |
state |
|
real(kind=wp), | intent(out) | :: | g |
event function: . The event is located when . |
Interface to the ddeabm_with_event_class_vec scalar event function.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_with_event_class_vec), | intent(inout) | :: | me | |||
real(kind=wp), | intent(in) | :: | t |
time |
||
real(kind=wp), | intent(in), | dimension(:) | :: | x |
state |
|
real(kind=wp), | intent(out), | dimension(:) | :: | g |
vector event function: . The event is located when any . |
|
integer, | intent(in), | optional | :: | ig |
the event function to compute:
If |
The main DDEABM integration class.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer, | private | :: | neq | = | 0 |
the number of (first order) differential equations to be integrated (>=0) |
|
integer, | private | :: | maxnum | = | 500 |
the expense of solving the problem is monitored by counting the
number of steps attempted. when this exceeds |
|
procedure(deriv_func), | private, | pointer | :: | df | => | null() |
to define the system of first order differential equations
which is to be solved. for the given values of |
procedure(report_func), | private, | pointer | :: | report | => | null() |
a function for reporting the integration steps (optional).
this can be associated in ddeabm_initialize or
ddeabm_with_event_initialize, and is enabled
using the |
logical, | private | :: | scalar_tols | = | .true. | ||
real(kind=wp), | private, | allocatable, dimension(:) | :: | rtol |
the user input rel tol |
||
real(kind=wp), | private, | allocatable, dimension(:) | :: | atol |
the user input abs tol |
||
real(kind=wp), | private, | allocatable, dimension(:) | :: | rtol_tmp |
rel tol used internally |
||
real(kind=wp), | private, | allocatable, dimension(:) | :: | atol_tmp |
abs tol used internally |
||
integer, | private, | dimension(4) | :: | info | = | 0 |
info array |
integer, | private | :: | initial_step_mode | = | 1 |
how to choose the initial step |
|
real(kind=wp), | private | :: | initial_step_size | = | 0.0_wp |
when |
|
integer, | private | :: | icount | = | 0 |
replaced |
|
real(kind=wp), | private | :: | tprev | = | 0.0_wp |
replaced |
|
real(kind=wp), | private, | dimension(:), allocatable | :: | two |
array |
||
real(kind=wp), | private, | dimension(:), allocatable | :: | gstr |
array |
||
real(kind=wp), | private, | allocatable, dimension(:) | :: | ypout |
|
||
real(kind=wp), | private, | allocatable, dimension(:) | :: | yp |
|
||
real(kind=wp), | private, | allocatable, dimension(:) | :: | yy |
|
||
real(kind=wp), | private, | allocatable, dimension(:) | :: | wt |
|
||
real(kind=wp), | private, | allocatable, dimension(:) | :: | p |
|
||
real(kind=wp), | private, | allocatable, dimension(:, :) | :: | phi |
|
||
real(kind=wp), | private, | dimension(12) | :: | alpha | = | 0.0_wp | |
real(kind=wp), | private, | dimension(12) | :: | beta | = | 0.0_wp | |
real(kind=wp), | private, | dimension(12) | :: | psi | = | 0.0_wp | |
real(kind=wp), | private, | dimension(12) | :: | v | = | 0.0_wp | |
real(kind=wp), | private, | dimension(12) | :: | w | = | 0.0_wp | |
real(kind=wp), | private, | dimension(13) | :: | sig | = | 0.0_wp | |
real(kind=wp), | private, | dimension(13) | :: | g | = | 0.0_wp | |
real(kind=wp), | private, | dimension(11) | :: | gi | = | 0.0_wp | |
real(kind=wp), | private | :: | h | = | 0.0_wp |
the step size |
|
real(kind=wp), | private | :: | eps | = | 0.0_wp |
if the tolerances have been increased by the
code ( |
|
real(kind=wp), | private | :: | x | = | 0.0_wp |
contains the current value of the
independent variable, i.e. the farthest point
integration has reached. this will be different
from |
|
real(kind=wp), | private | :: | xold | = | 0.0_wp | ||
real(kind=wp), | private | :: | hold | = | 0.0_wp | ||
real(kind=wp), | private | :: | told | = | 0.0_wp | ||
real(kind=wp), | private | :: | delsgn | = | 0.0_wp | ||
real(kind=wp), | private | :: | tstop | = | 0.0_wp |
(see |
|
real(kind=wp), | private | :: | twou | = | 0.0_wp |
machine dependent parameter |
|
real(kind=wp), | private | :: | fouru | = | 0.0_wp |
machine dependent parameter |
|
logical, | private | :: | start | = | .false. |
indicator for dsteps code |
|
logical, | private | :: | phase1 | = | .false. |
indicator for dsteps code |
|
logical, | private | :: | nornd | = | .false. |
indicator for dsteps code |
|
logical, | private | :: | stiff | = | .false. |
indicator for stiffness detection |
|
logical, | private | :: | intout | = | .false. |
indicator for intermediate-output |
|
integer, | private | :: | ns | = | 0 | ||
integer, | private | :: | kord | = | 0 | ||
integer, | private | :: | kold | = | 0 | ||
integer, | private | :: | init | = | 0 |
initialization indicator |
|
integer, | private | :: | ksteps | = | 0 |
counter for attempted steps |
|
integer, | private | :: | kle4 | = | 0 |
step counter for stiffness detection |
|
integer, | private | :: | iquit | = | 0 |
termination flag |
|
integer, | private | :: | kprev | = | 0 | ||
integer, | private | :: | ivc | = | 0 | ||
integer, | private, | dimension(10) | :: | iv | = | 0 | |
integer, | private | :: | kgi | = | 0 | ||
logical, | private | :: | error | = | .false. |
will be set in |
procedure, public, non_overridable :: initialize => ddeabm_initialize | |
procedure, public, non_overridable :: integrate => ddeabm_wrapper | |
procedure, public, non_overridable :: destroy => destroy_ddeabm | |
procedure, public, non_overridable :: first_call => ddeabm_new_problem | |
procedure, public, non_overridable :: interpolate => ddeabm_interp | ../../ state interpolation function |
procedure, public, non_overridable :: stop_integration => ddeabm_stop_integration | ../../ user can call this in |
procedure, private, non_overridable :: ddeabm | |
procedure, private, non_overridable :: ddes | |
procedure, private, non_overridable :: dhstrt | |
procedure, private, non_overridable :: dsteps |
A version of ddeabm_class with event location (root finding).
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
real(kind=wp), | private | :: | tol | = | 0.0_wp |
tolerance for root finding |
|
real(kind=wp), | private | :: | t_saved | = | 0.0_wp |
time of last successful step (used to continue after a root finding) |
|
real(kind=wp), | private, | dimension(:), allocatable | :: | x_saved |
state of last successful step. size |
||
procedure(event_func), | private, | pointer | :: | gfunc | => | null() |
event function: g(t,x)=0 at event |
procedure(bracket_func), | private, | pointer | :: | bracket | => | null() |
a function for determining if the root is bracketed. the default just checks for a sign change in the event function. the user can specify this to use other conditions. |
procedure, public, non_overridable :: initialize => ddeabm_initialize | |
procedure, public, non_overridable :: integrate => ddeabm_wrapper | |
procedure, public, non_overridable :: destroy => destroy_ddeabm | |
procedure, public, non_overridable :: first_call => ddeabm_new_problem | |
procedure, public, non_overridable :: interpolate => ddeabm_interp | ../../ state interpolation function |
procedure, public, non_overridable :: stop_integration => ddeabm_stop_integration | ../../ user can call this in |
procedure, public, non_overridable :: initialize_event => ddeabm_with_event_initialize | |
procedure, public, non_overridable :: integrate_to_event => ddeabm_with_event_wrapper | main routine for integration to an event (note that integrate can also still be used from ddeabm_class) |
A version of ddeabm_class with event location (root finding) of a vector function (it will stop on the first root it encounters)
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer, | private | :: | n_g_eqns | = | 0 |
number of |
|
real(kind=wp), | private, | dimension(:), allocatable | :: | tol |
tolerance for root finding.
this is a vector (size |
||
real(kind=wp), | private | :: | t_saved | = | 0.0_wp |
time of last successful step (used to continue after a root finding) |
|
real(kind=wp), | private, | dimension(:), allocatable | :: | x_saved |
state of last successful step. size |
||
procedure(event_func_vec), | private, | pointer | :: | gfunc | => | null() |
event function: g(t,x)=0 at event |
procedure(bracket_func_vec), | private, | pointer | :: | bracket | => | null() |
a function for determining if the root is bracketed. the default just checks for a sign change in the event function. the user can specify this to use other conditions. |
procedure, public, non_overridable :: initialize => ddeabm_initialize | |
procedure, public, non_overridable :: integrate => ddeabm_wrapper | |
procedure, public, non_overridable :: destroy => destroy_ddeabm | |
procedure, public, non_overridable :: first_call => ddeabm_new_problem | |
procedure, public, non_overridable :: interpolate => ddeabm_interp | ../../ state interpolation function |
procedure, public, non_overridable :: stop_integration => ddeabm_stop_integration | ../../ user can call this in |
procedure, public, non_overridable :: initialize_event => ddeabm_with_event_initialize_vec | |
procedure, public, non_overridable :: integrate_to_event => ddeabm_with_event_wrapper_vec | main routine for integration to an event (note that integrate can also still be used from ddeabm_class) |
Compute the maximum norm of the first n
elements of vector v
.
Replacement for the original SLATEC routine.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=wp), | intent(in), | dimension(:) | :: | v | ||
integer, | intent(in) | :: | n |
Call this to abort the integration if there is an error.
(Can be called in df
by the user.) The step
will be considered invalid and not returned.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_class), | intent(inout) | :: | me |
Call this to indicate that a new problem is being solved.
It sets info(1) = 0
(see ddeabm documentation).
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_class), | intent(inout) | :: | me |
Initialize the ddeabm_class, and set the variables that cannot be changed during a problem.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_class), | intent(inout) | :: | me | |||
integer, | intent(in) | :: | neq |
number of equations to be integrated |
||
integer, | intent(in) | :: | maxnum |
maximum number of integration steps allowed |
||
procedure(deriv_func) | :: | df |
derivative function dx/dt |
|||
real(kind=wp), | intent(in), | dimension(:) | :: | rtol |
relative tolerance for integration (see ddeabm) |
|
real(kind=wp), | intent(in), | dimension(:) | :: | atol |
absolution tolerance for integration (see ddeabm) |
|
procedure(report_func), | optional | :: | report |
reporting function |
||
integer, | intent(in), | optional | :: | initial_step_mode |
how to choose the initial step |
|
real(kind=wp), | intent(in), | optional | :: | initial_step_size |
for |
Initialize ddeabm_with_event_class class, and set the variables that cannot be changed during a problem.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_with_event_class), | intent(inout) | :: | me | |||
integer, | intent(in) | :: | neq |
number of equations to be integrated |
||
integer, | intent(in) | :: | maxnum |
maximum number of integration steps allowed |
||
procedure(deriv_func) | :: | df |
derivative function |
|||
real(kind=wp), | intent(in), | dimension(:) | :: | rtol |
relative tolerance for integration (see ddeabm) |
|
real(kind=wp), | intent(in), | dimension(:) | :: | atol |
absolution tolerance for integration (see ddeabm) |
|
procedure(event_func) | :: | g |
Event function . This should be a smooth function than can have values and . When (within the tolerance), then a root has been located and the integration will stop. |
|||
real(kind=wp), | intent(in) | :: | root_tol |
tolerance for the root finding |
||
procedure(report_func), | optional | :: | report |
reporting function |
||
procedure(bracket_func), | optional | :: | bracket |
root bracketing function. if not present, the default is used. |
||
integer, | intent(in), | optional | :: | initial_step_mode |
how to choose the initial step |
|
real(kind=wp), | intent(in), | optional | :: | initial_step_size |
for |
Initialize ddeabm_with_event_class_vec class, and set the variables that cannot be changed during a problem.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_with_event_class_vec), | intent(inout) | :: | me | |||
integer, | intent(in) | :: | neq |
number of equations to be integrated |
||
integer, | intent(in) | :: | maxnum |
maximum number of integration steps allowed |
||
procedure(deriv_func) | :: | df |
derivative function |
|||
real(kind=wp), | intent(in), | dimension(:) | :: | rtol |
relative tolerance for integration (see ddeabm) |
|
real(kind=wp), | intent(in), | dimension(:) | :: | atol |
absolution tolerance for integration (see ddeabm) |
|
integer, | intent(in) | :: | ng |
number of event functions in |
||
procedure(event_func_vec) | :: | g |
Event function . This should be a smooth function than can have values and . When (within the tolerance), then a root has been located and the integration will stop. |
|||
real(kind=wp), | intent(in), | dimension(:) | :: | root_tol |
tolerance for the root finding.
This should be sized |
|
procedure(report_func), | optional | :: | report |
reporting function |
||
procedure(bracket_func_vec), | optional | :: | bracket |
root bracketing function. if not present, the default is used. |
||
integer, | intent(in), | optional | :: | initial_step_mode |
how to choose the initial step |
|
real(kind=wp), | intent(in), | optional | :: | initial_step_size |
for |
Destructor for ddeabm_class.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_class), | intent(out) | :: | me |
Wrapper routine for ddeabm.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_class), | intent(inout) | :: | me | |||
real(kind=wp), | intent(inout) | :: | t | |||
real(kind=wp), | intent(inout), | dimension(me%neq) | :: | y | ||
real(kind=wp), | intent(in) | :: | tout |
time at which a solution is desired. |
||
real(kind=wp), | intent(in), | optional | :: | tstop |
for some problems it may not be permissible to integrate
past a point |
|
integer, | intent(out) | :: | idid |
indicator reporting what the code did. you must monitor this integer variable to decide what action to take next. |
||
integer, | intent(in), | optional | :: | integration_mode |
Step mode:
1 - normal integration from |
|
real(kind=wp), | intent(in), | optional | :: | tstep |
Fixed time step to use for |
Wrapper routine for ddeabm, with event finding.
It will integrate until g(t,x)=0
or t=tmax
(whichever comes first).
Note that a root at the initial time is ignored (user should check for
this manually)
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_with_event_class), | intent(inout) | :: | me | |||
real(kind=wp), | intent(inout) | :: | t | |||
real(kind=wp), | intent(inout), | dimension(me%neq) | :: | y | ||
real(kind=wp), | intent(in) | :: | tmax |
max time at which a solution is desired.
(if root not located, it will integrate to |
||
real(kind=wp), | intent(in), | optional | :: | tstop |
for some problems it may not be permissible to integrate
past a point |
|
integer, | intent(out) | :: | idid |
indicator reporting what the code did.
you must monitor this integer variable to
decide what action to take next.
|
||
real(kind=wp), | intent(out) | :: | gval |
value of the event function |
||
integer, | intent(in), | optional | :: | integration_mode |
Step mode:
1 - normal integration from |
|
real(kind=wp), | intent(in), | optional | :: | tstep |
Fixed time step to use for reporting and evaluation of event function. If not present, then default integrator steps are used. |
|
logical, | intent(in), | optional | :: | continue |
to continue after a previous event location. This option can be used after a previous call to this routine has found an event, to continue integration to the next event without having to restart the integration. It will not report the initial point (which would have been reported as the last point of the previous call). |
Interpolation function. Can be used for dense output after a step. It calls the low-level routine dintp.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_class), | intent(inout) | :: | me | |||
real(kind=wp), | intent(in) | :: | tc |
point at which solution is desired |
||
real(kind=wp), | intent(out), | dimension(me%neq) | :: | yc |
interpolated state at |
Wrapper routine for ddeabm, with event finding for multiple functions.
It will integrate until any g(t,x)=0
or t=tmax
(whichever comes first).
Note that a root at the initial time is ignored (user should check for
this manually)
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_with_event_class_vec), | intent(inout) | :: | me | |||
real(kind=wp), | intent(inout) | :: | t | |||
real(kind=wp), | intent(inout), | dimension(me%neq) | :: | y | ||
real(kind=wp), | intent(in) | :: | tmax |
max time at which a solution is desired.
(if root not located, it will integrate to |
||
real(kind=wp), | intent(in), | optional | :: | tstop |
for some problems it may not be permissible to integrate
past a point |
|
integer, | intent(out) | :: | idid |
indicator reporting what the code did.
you must monitor this integer variable to
decide what action to take next.
|
||
real(kind=wp), | intent(out), | dimension(:) | :: | gval |
value of the event functions |
|
integer, | intent(in), | optional | :: | integration_mode |
Step mode:
1 - normal integration from |
|
real(kind=wp), | intent(in), | optional | :: | tstep |
Fixed time step to use for reporting and evaluation of event function. If not present, then default integrator steps are used. |
|
logical, | intent(in), | optional | :: | continue |
to continue after a previous event location. This option can be used after a previous call to this routine has found an event, to continue integration to the next event without having to restart the integration. It will not report the initial point (which would have been reported as the last point of the previous call). |
solve an initial value problem in ordinary differential equations using an adams-bashforth method.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_class), | intent(inout) | :: | me | |||
integer, | intent(in) | :: | neq |
the number of (first order) differential equations to be integrated. (neq >= 1) |
||
real(kind=wp), | intent(inout) | :: | t |
value of the independent variable. on input, set it to the initial point of the integration. the code changes its value. on output, the solution was successfully advanced to the output value of t. |
||
real(kind=wp), | intent(inout), | dimension(neq) | :: | y | ||
real(kind=wp), | intent(in) | :: | tout | |||
integer, | intent(inout), | dimension(4) | :: | info | ||
real(kind=wp), | intent(inout), | dimension(:) | :: | rtol | ||
real(kind=wp), | intent(inout), | dimension(:) | :: | atol | ||
integer, | intent(out) | :: | idid |
ddeabm merely allocates storage for ddes to relieve the user of the inconvenience of a long call list. consequently ddes is used as described in the comments for ddeabm.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_class), | intent(inout) | :: | me | |||
integer, | intent(in) | :: | neq | |||
real(kind=wp), | intent(inout) | :: | t | |||
real(kind=wp), | intent(inout), | dimension(neq) | :: | y | ||
real(kind=wp), | intent(in) | :: | tout | |||
integer, | intent(inout), | dimension(4) | :: | info | ||
real(kind=wp), | intent(inout), | dimension(:) | :: | rtol | ||
real(kind=wp), | intent(inout), | dimension(:) | :: | atol | ||
integer, | intent(inout) | :: | idid | |||
real(kind=wp), | intent(inout), | dimension(neq) | :: | ypout | ||
real(kind=wp), | intent(inout), | dimension(neq) | :: | yp | ||
real(kind=wp), | intent(inout), | dimension(neq) | :: | yy | ||
real(kind=wp), | intent(inout), | dimension(neq) | :: | wt | ||
real(kind=wp), | intent(inout), | dimension(neq) | :: | p | ||
real(kind=wp), | intent(inout), | dimension(neq, 16) | :: | phi | ||
real(kind=wp), | intent(inout), | dimension(12) | :: | alpha | ||
real(kind=wp), | intent(inout), | dimension(12) | :: | beta | ||
real(kind=wp), | intent(inout), | dimension(12) | :: | psi | ||
real(kind=wp), | intent(inout), | dimension(12) | :: | v | ||
real(kind=wp), | intent(inout), | dimension(12) | :: | w | ||
real(kind=wp), | intent(inout), | dimension(13) | :: | sig | ||
real(kind=wp), | intent(inout), | dimension(13) | :: | g | ||
real(kind=wp), | intent(inout), | dimension(11) | :: | gi | ||
real(kind=wp), | intent(inout) | :: | h | |||
real(kind=wp), | intent(inout) | :: | eps | |||
real(kind=wp), | intent(inout) | :: | x | |||
real(kind=wp), | intent(inout) | :: | xold | |||
real(kind=wp), | intent(inout) | :: | hold | |||
real(kind=wp), | intent(inout) | :: | told | |||
real(kind=wp), | intent(inout) | :: | delsgn | |||
real(kind=wp), | intent(inout) | :: | tstop | |||
real(kind=wp), | intent(inout) | :: | twou | |||
real(kind=wp), | intent(inout) | :: | fouru | |||
logical, | intent(inout) | :: | start | |||
logical, | intent(inout) | :: | phase1 | |||
logical, | intent(inout) | :: | nornd | |||
logical, | intent(inout) | :: | stiff | |||
logical, | intent(inout) | :: | intout | |||
integer, | intent(inout) | :: | ns | |||
integer, | intent(inout) | :: | kord | |||
integer, | intent(inout) | :: | kold | |||
integer, | intent(inout) | :: | init | |||
integer, | intent(inout) | :: | ksteps | |||
integer, | intent(inout) | :: | kle4 | |||
integer, | intent(inout) | :: | iquit | |||
integer, | intent(inout) | :: | kprev | |||
integer, | intent(inout) | :: | ivc | |||
integer, | intent(inout), | dimension(10) | :: | iv | ||
integer, | intent(inout) | :: | kgi |
Computes a starting step size to be used in solving initial value problems in ordinary differential equations.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_class), | intent(inout) | :: | me | |||
integer, | intent(in) | :: | neq |
the number of (first order) differential equations to be integrated. |
||
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(neq) | :: | y |
the vector of initial values of the |
|
real(kind=wp), | intent(in), | dimension(neq) | :: | yprime |
the vector of derivatives of the |
|
real(kind=wp), | intent(in), | dimension(neq) | :: | etol |
the vector of error tolerances corresponding to
the |
|
integer, | intent(in) | :: | morder |
the order of the formula which will be used by the initial value method for taking the first integration step. |
||
real(kind=wp), | intent(in) | :: | small |
a small positive machine dependent constant
which is used for protecting against computations with
numbers which are too small relative to the precision of
floating point arithmetic. |
||
real(kind=wp), | intent(in) | :: | big |
a large positive machine dependent constant which is used for preventing machine overflows. a reasonable choice is to set big to (approximately) the square root of the largest double precision number which can be held in the machine. |
||
real(kind=wp), | intent(inout), | dimension(neq) | :: | spy |
work array of length |
|
real(kind=wp), | intent(inout), | dimension(neq) | :: | pv |
work array of length |
|
real(kind=wp), | intent(inout), | dimension(neq) | :: | yp |
work array of length |
|
real(kind=wp), | intent(inout), | dimension(neq) | :: | sf |
work array of length |
|
real(kind=wp), | intent(out) | :: | h |
appropriate starting step size to be attempted by the differential equation method. |
approximate the solution at xout
by evaluating the
polynomial computed in dsteps at xout
. must be used in
conjunction with dsteps.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=wp), | intent(in) | :: | x | |||
real(kind=wp), | intent(in), | dimension(neqn) | :: | y | ||
real(kind=wp), | intent(in) | :: | xout |
point at which solution is desired |
||
real(kind=wp), | intent(out), | dimension(neqn) | :: | yout |
solution at xout |
|
real(kind=wp), | intent(out), | dimension(neqn) | :: | ypout |
derivative of solution at xout |
|
integer, | intent(in) | :: | neqn | |||
integer, | intent(in) | :: | kold | |||
real(kind=wp), | intent(in), | dimension(neqn, 16) | :: | phi | ||
integer, | intent(in) | :: | ivc | |||
integer, | intent(in), | dimension(10) | :: | iv | ||
integer, | intent(in) | :: | kgi | |||
real(kind=wp), | intent(in), | dimension(11) | :: | gi | ||
real(kind=wp), | intent(in), | dimension(12) | :: | alpha | ||
real(kind=wp), | intent(in), | dimension(13) | :: | og | ||
real(kind=wp), | intent(in), | dimension(12) | :: | ow | ||
real(kind=wp), | intent(in) | :: | ox | |||
real(kind=wp), | intent(in), | dimension(neqn) | :: | oy |
integrate a system of first order ordinary differential equations one step.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_class), | intent(inout) | :: | me | |||
integer, | intent(in) | :: | neqn |
number of equations to be integrated |
||
real(kind=wp), | intent(inout), | dimension(neqn) | :: | y |
solution vector at x |
|
real(kind=wp), | intent(inout) | :: | x |
independent variable |
||
real(kind=wp), | intent(inout) | :: | h |
appropriate step size for next step. normally determined by code |
||
real(kind=wp), | intent(inout) | :: | eps |
local error tolerance |
||
real(kind=wp), | intent(inout), | dimension(neqn) | :: | wt |
vector of weights for error criterion |
|
logical, | intent(inout) | :: | start |
logical variable set .true. for first step, .false. otherwise |
||
real(kind=wp), | intent(inout) | :: | hold |
step size used for last successful step |
||
integer, | intent(inout) | :: | k |
appropriate order for next step (determined by code) |
||
integer, | intent(inout) | :: | kold |
order used for last successful step |
||
logical, | intent(inout) | :: | crash |
logical variable set .true. when no step can be taken, .false. otherwise. |
||
real(kind=wp), | intent(inout), | dimension(neqn, 16) | :: | phi |
required for the interpolation subroutine dintp |
|
real(kind=wp), | intent(inout), | dimension(neqn) | :: | p |
required for the interpolation subroutine dintp |
|
real(kind=wp), | intent(inout), | dimension(neqn) | :: | yp |
derivative of solution vector at x after successful step |
|
real(kind=wp), | intent(inout), | dimension(12) | :: | psi | ||
real(kind=wp), | intent(inout), | dimension(12) | :: | alpha |
required for the interpolation subroutine dintp |
|
real(kind=wp), | intent(inout), | dimension(12) | :: | beta | ||
real(kind=wp), | intent(inout), | dimension(13) | :: | sig | ||
real(kind=wp), | intent(inout), | dimension(12) | :: | v | ||
real(kind=wp), | intent(inout), | dimension(12) | :: | w |
required for the interpolation subroutine dintp |
|
real(kind=wp), | intent(inout), | dimension(13) | :: | g |
required for the interpolation subroutine dintp |
|
logical, | intent(inout) | :: | phase1 | |||
integer, | intent(inout) | :: | ns |
number of dsteps taken with size h |
||
logical, | intent(inout) | :: | nornd | |||
integer, | intent(inout) | :: | ksteps |
counter on attempted steps |
||
real(kind=wp), | intent(inout) | :: | twou |
2.*u where u is machine unit roundoff quantity |
||
real(kind=wp), | intent(inout) | :: | fouru |
4.*u where u is machine unit roundoff quantity |
||
real(kind=wp), | intent(inout) | :: | xold |
required for the interpolation subroutine dintp |
||
integer, | intent(inout) | :: | kprev | |||
integer, | intent(inout) | :: | ivc |
required for the interpolation subroutine dintp |
||
integer, | intent(inout), | dimension(10) | :: | iv |
required for the interpolation subroutine dintp |
|
integer, | intent(inout) | :: | kgi |
required for the interpolation subroutine dintp |
||
real(kind=wp), | intent(inout), | dimension(11) | :: | gi |
required for the interpolation subroutine dintp |
Report an error message.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
character(len=*), | intent(in) | :: | subrou |
subroutine where error occurred |
||
character(len=*), | intent(in) | :: | messg |
error message |
||
integer, | intent(in) | :: | nerr |
error code |
||
integer, | intent(in) | :: | level |
[Not used here] * -1: warning message (once), * 0: warning message, * 1: recoverable error, * 2: fatal error |