A version of ddeabm_class with event location (root finding).
Call initialize_event to set up the integration.
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. |
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 |
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 |
Destructor for ddeabm_class.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(ddeabm_class), | intent(out) | :: | 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 |
state interpolation function
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 |
user can call this in df
routine to stop the integration
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 |
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 |
main routine for integration to an event (note that integrate can also still be used from ddeabm_class)
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). |
type, extends(ddeabm_class), public :: ddeabm_with_event_class !! A version of [[ddeabm_class]] with event location (root finding). !! !! Call *initialize_event* to set up the integration. private real(wp) :: tol = 0.0_wp !! tolerance for root finding real(wp) :: t_saved = 0.0_wp !! time of last successful step !! (used to continue after a root finding) real(wp), dimension(:), allocatable :: x_saved !! state of last successful step. size `(neq)`. !! (used to continue after a root finding) procedure(event_func), pointer :: gfunc => null() !! event function: g(t,x)=0 at event procedure(bracket_func), 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. contains private procedure, non_overridable, public :: initialize_event => ddeabm_with_event_initialize procedure, non_overridable, public :: 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]]) end type ddeabm_with_event_class