Title: | Solve Delay Differential Equations |
---|---|
Description: | Solves ordinary and delay differential equations, where the objective function is written in either R or C. Suitable only for non-stiff equations, the solver uses a 'Dormand-Prince' method that allows interpolation of the solution at any point. This approach is as described by Hairer, Norsett and Wanner (1993) <ISBN:3540604529>. Support is also included for iterating difference equations. |
Authors: | Rich FitzJohn [aut, cre], Wes Hinsley [aut], Imperial College of Science, Technology and Medicine [cph] |
Maintainer: | Rich FitzJohn <[email protected]> |
License: | MIT + file LICENSE |
Version: | 1.0.7 |
Built: | 2024-11-22 05:37:51 UTC |
Source: | https://github.com/mrc-ide/dde |
Solve a difference (or recurrence) equation by iterating it a number of times.
difeq(y, steps, target, parms, ..., n_out = 0L, n_history = 0L, grow_history = FALSE, return_history = n_history > 0, dllname = "", parms_are_real = TRUE, ynames = names(y), outnames = NULL, return_by_column = TRUE, return_initial = TRUE, return_step = TRUE, return_output_with_y = TRUE, restartable = FALSE, return_minimal = FALSE) difeq_continue(obj, steps, y = NULL, ..., copy = FALSE, parms = NULL, return_history = NULL, return_by_column = NULL, return_initial = NULL, return_step = NULL, return_output_with_y = NULL, restartable = NULL) yprev(step, i = NULL)
difeq(y, steps, target, parms, ..., n_out = 0L, n_history = 0L, grow_history = FALSE, return_history = n_history > 0, dllname = "", parms_are_real = TRUE, ynames = names(y), outnames = NULL, return_by_column = TRUE, return_initial = TRUE, return_step = TRUE, return_output_with_y = TRUE, restartable = FALSE, return_minimal = FALSE) difeq_continue(obj, steps, y = NULL, ..., copy = FALSE, parms = NULL, return_history = NULL, return_by_column = NULL, return_initial = NULL, return_step = NULL, return_output_with_y = NULL, restartable = NULL) yprev(step, i = NULL)
y |
The initial state of the system. Must be a numeric
vector (and will be passed through |
steps |
A vector of steps to return the system at. The
first step is taken as step zero, and the solution will
be recorded at every other step in the vector. So to step a
system from time zero to times 1, 2, 3, ..., n use 0:n. Must be
integer values and will be passed through |
target |
The target function to advance. This can either be
an R function taking arguments |
parms |
Parameters to pass through to the difference function |
... |
Dummy arguments - nothing is allowed here, but this means that all further arguments must be specified by name (not order) so I can easily reorder them later on. |
n_out |
The number of output variables (in addition to the
difference equation variables). If given, then an R function
must return an attribute |
n_history |
The number of iterations of history to save during the simulation. By default, no history is saved. |
grow_history |
Logical indicating if history should be grown
during the simulation. If |
return_history |
Logical indicating if history is to be
returned. By default, history is returned if |
dllname |
Name of the shared library (without extension) to
find the function |
parms_are_real |
Logical, indicating if |
ynames |
Logical, indicating if the output should be named
following the names of the input vector |
outnames |
An optional character vector, used when
|
return_by_column |
Logical, indicating if the output should be
returned organised by column (rather than row). This incurs a
slight cost for transposing the matrices. If you can work with
matrices that are transposed relative to |
return_initial |
Logical, indicating if the output should
include the initial conditions. Specifying |
return_step |
Logical, indicating if a row (or column if
|
return_output_with_y |
Logical, indicating if the output
should be bound together with the returned matrix |
restartable |
Logical, indicating if the problem should be
restartable. If |
return_minimal |
Shorthand option - if set to |
obj |
An object to continue from; this must be the results of
running a simulation with the option |
copy |
Logical, indicating if the pointer should be copied
before continuing. If |
step |
The step to access (not that this is not an offset,
but the actual step; within your target function you'd write
things like |
i |
index within the state vector |
# Here is a really simple equation that just increases by 'p' each # time (p is the parameter vector and could be any R structure). rhs <- function(i, y, p) y + p y0 <- 1 t <- 0:10 p <- 5 dde::difeq(y0, t, rhs, p)
# Here is a really simple equation that just increases by 'p' each # time (p is the parameter vector and could be any R structure). rhs <- function(i, y, p) y + p y0 <- 1 t <- 0:10 p <- 5 dde::difeq(y0, t, rhs, p)
Solve a replicate set of difference (or recurrence) equation by
iterating it a number of times. This is a wrapper around
difeq
that does not (yet) do anything clever to
avoid many allocations.
difeq_replicate(n, y, steps, target, parms, ..., n_out = 0L, n_history = 0L, grow_history = FALSE, return_history = n_history > 0, dllname = "", parms_are_real = TRUE, ynames = NULL, outnames = NULL, return_by_column = TRUE, return_initial = TRUE, return_step = TRUE, return_output_with_y = TRUE, restartable = FALSE, return_minimal = FALSE, as_array = TRUE)
difeq_replicate(n, y, steps, target, parms, ..., n_out = 0L, n_history = 0L, grow_history = FALSE, return_history = n_history > 0, dllname = "", parms_are_real = TRUE, ynames = NULL, outnames = NULL, return_by_column = TRUE, return_initial = TRUE, return_step = TRUE, return_output_with_y = TRUE, restartable = FALSE, return_minimal = FALSE, as_array = TRUE)
n |
Number of replicates. It is an error to request zero replicates. |
y |
The initial state of the system. Must be either a
numeric vector or a |
steps |
A vector of steps to return the system at. The
first step is taken as step zero, and the solution will
be recorded at every other step in the vector. So to step a
system from time zero to times 1, 2, 3, ..., n use 0:n. Must be
integer values and will be passed through |
target |
The target function to advance. This can either be
an R function taking arguments |
parms |
Parameters to pass through to the difference function |
... |
Dummy arguments - nothing is allowed here, but this means that all further arguments must be specified by name (not order) so I can easily reorder them later on. |
n_out |
The number of output variables (in addition to the
difference equation variables). If given, then an R function
must return an attribute |
n_history |
The number of iterations of history to save during the simulation. By default, no history is saved. |
grow_history |
Logical indicating if history should be grown
during the simulation. If |
return_history |
Logical indicating if history is to be
returned. By default, history is returned if |
dllname |
Name of the shared library (without extension) to
find the function |
parms_are_real |
Logical, indicating if |
ynames |
Logical, indicating if the output should be named
following the names of the input vector |
outnames |
An optional character vector, used when
|
return_by_column |
Logical, indicating if the output should be
returned organised by column (rather than row). This incurs a
slight cost for transposing the matrices. If you can work with
matrices that are transposed relative to |
return_initial |
Logical, indicating if the output should
include the initial conditions. Specifying |
return_step |
Logical, indicating if a row (or column if
|
return_output_with_y |
Logical, indicating if the output
should be bound together with the returned matrix |
restartable |
Logical, indicating if the problem should be
restartable. If |
return_minimal |
Shorthand option - if set to |
as_array |
(Defunct) Logical, indicating if the output should
be converted into an array. If |
It is not currently possible to replicate over a set of parameters at once yet; the same parameter set will be used for all replications.
The details of how replication is done here are all considered
implementation details and are up for change in the future - in
particular if the models are run in turn or simultaneously (and
the effect that has on the random number stream). Logic around
naming output may change in future too; note that varying names in
the y
here will have some unexpected behaviours.
# Here is a really simple equation that does a random walk with # steps that are normally distributed: rhs <- function(i, y, p) y + runif(1) y0 <- 1 t <- 0:10 p <- 5 dde::difeq_replicate(10, y0, t, rhs, p)
# Here is a really simple equation that does a random walk with # steps that are normally distributed: rhs <- function(i, y, p) y + runif(1) y0 <- 1 t <- 0:10 p <- 5 dde::difeq_replicate(10, y0, t, rhs, p)
Integrate an ODE or DDE with dopri.
dopri(y, times, func, parms, ..., n_out = 0L, output = NULL, rtol = 1e-06, atol = 1e-06, step_size_min = 0, step_size_max = Inf, step_size_initial = 0, step_max_n = 100000L, step_size_min_allow = FALSE, tcrit = NULL, event_time = NULL, event_function = NULL, method = "dopri5", stiff_check = 0, verbose = FALSE, callback = NULL, n_history = 0, grow_history = FALSE, return_history = n_history > 0, dllname = "", parms_are_real = TRUE, ynames = names(y), outnames = NULL, return_by_column = TRUE, return_initial = TRUE, return_time = TRUE, return_output_with_y = TRUE, return_statistics = FALSE, restartable = FALSE, return_minimal = FALSE) dopri5(y, times, func, parms, ...) dopri853(y, times, func, parms, ...) dopri_continue(obj, times, y = NULL, ..., copy = FALSE, parms = NULL, tcrit = NULL, return_history = NULL, return_by_column = NULL, return_initial = NULL, return_statistics = NULL, return_time = NULL, return_output_with_y = NULL, restartable = NULL) ylag(t, i = NULL)
dopri(y, times, func, parms, ..., n_out = 0L, output = NULL, rtol = 1e-06, atol = 1e-06, step_size_min = 0, step_size_max = Inf, step_size_initial = 0, step_max_n = 100000L, step_size_min_allow = FALSE, tcrit = NULL, event_time = NULL, event_function = NULL, method = "dopri5", stiff_check = 0, verbose = FALSE, callback = NULL, n_history = 0, grow_history = FALSE, return_history = n_history > 0, dllname = "", parms_are_real = TRUE, ynames = names(y), outnames = NULL, return_by_column = TRUE, return_initial = TRUE, return_time = TRUE, return_output_with_y = TRUE, return_statistics = FALSE, restartable = FALSE, return_minimal = FALSE) dopri5(y, times, func, parms, ...) dopri853(y, times, func, parms, ...) dopri_continue(obj, times, y = NULL, ..., copy = FALSE, parms = NULL, tcrit = NULL, return_history = NULL, return_by_column = NULL, return_initial = NULL, return_statistics = NULL, return_time = NULL, return_output_with_y = NULL, restartable = NULL) ylag(t, i = NULL)
y |
Initial conditions for the integration |
times |
Times where output is needed. Unlike |
func |
Function to integrate. Can be an R function of
arguments |
parms |
Parameters to pass through to the derivatives. |
... |
Dummy arguments - nothing is allowed here, but this means that all further arguments must be specified by name (not order) so I can easily reorder them later on. |
n_out |
Number of "output" variables (not differential
equation variables) to compute via the routine |
output |
The output routine; either an R function taking
arguments |
rtol |
The per-step relative tolerance. The total accuracy will be less than this. |
atol |
The per-step absolute tolerance. |
step_size_min |
The minimum step size. The actual minimum
used will be the largest of the absolute value of this
|
step_size_max |
The largest step size. By default there is
no maximum step size (Inf) so the solver can take as large a
step as it wants to. If you have short-lived fluctuations in
your rhs that the solver may skip over by accident, then specify
a smaller maximum step size here (or use |
step_size_initial |
The initial step size. By default the integrator will guess the step size automatically, but one can be given here instead. |
step_max_n |
The maximum number of steps allowed. If the
solver takes more steps than this it will throw an error. Note
the number of evaluations of |
step_size_min_allow |
Logical, indicating if when a step size
is driven down to |
tcrit |
An optional vector of critical times that the solver
must stop at (rather than interpolating over). This can include
an end time that we can't go past, or points within the
integration that must be stopped at exactly (for example cases
where the derivatives change abruptly). Note that this differs
from the interpretation of this parameter in deSolve; there
|
event_time |
Vector of times to fire events listed in
|
event_function |
Function to fire at events. For R models
( |
method |
The integration method to use, as a string. The
supported methods are |
stiff_check |
How often to check that the problem has become
stiff. If zero, then the problem is never checked, and if
positive then the problem is checked every |
verbose |
Be verbose, and print information about each step.
This may be useful for learning about models that misbehave.
Valid values are |
callback |
Callback function that can be used to make verbose
output more useful. This can be used to return more information
about the evaluation as it proceeds, generally as information
printed to the screen. The function must accept arguments
|
n_history |
Number of history points to retain. This needs to be greater than zero for delay differential equations to work. Alternatively, this may be greater than zero to return model outputs that can be inspected later. |
grow_history |
Logical indicating if history should be grown
during the simulation. If |
return_history |
Logical indicating if history should be
returned alongside the output or discarded. By default, history
is retained if |
dllname |
Name of the shared library (without extension) to
find the function |
parms_are_real |
Logical, indicating if |
ynames |
Logical, indicating if the output should be named
following the names of the input vector |
outnames |
An optional character vector, used when
|
return_by_column |
Logical, indicating if the output should be
returned organised by column (rather than row). This incurs a
slight cost for transposing the matrices. If you can work with
matrices that are transposed relative to |
return_initial |
Logical, indicating if the output should
include the initial conditions. Specifying |
return_time |
Logical, indicating if a row (or column if
|
return_output_with_y |
Logical, indicating if the output
should be bound together with the returned matrix |
return_statistics |
Logical, indicating if statistics about
the run should be included. If |
restartable |
Logical, indicating if the problem should be
restartable. If |
return_minimal |
Shorthand option - if set to |
obj |
An object to continue from; this must be the results of
running an integration with the option |
copy |
Logical, indicating if the pointer should be copied
before continuing. If |
t |
The time to access (not that this is not an offset,
but the actual time; within your target function you'd write
things like |
i |
index within the state vector |
Like deSolve::lsoda
, this function has many
arguments. This is far from ideal, and I would welcome any
approach for simplifying it a bit.
The options return_by_column
, return_initial
,
return_time
, return_output_with_y
exist because
these options all carry out modifications of the data at the end
of solving the ODE and this can incur a small but measurable cost.
When solving an ODE repeatedly (e.g., in the context of an MCMC or
optimisation) it may be useful to do as little as possible. For
simple problems this can save around 5-10% of the total
computational time (especially the transpose). The shorthand
option return_minimal
will set all to FALSE
when
used.
At present the return value is transposed relative to deSolve. This might change in future.
Debugging a failed integration can be difficult, but dopri
provides a couple of tools to get more information about where a
failure might have occurred. Most simply, one can pass
verbose = TRUE
which will print information about the
time and the step size at each point just before the step is
stated. Passing in verbose = dde:::VERBOSE_EVAL
will
print information just before every evaluation of the target
function (there are several evaluations per step).
However, this does not provide information about the state just
before failure. To get that, one must provide a callback
function - this is an R function that will be called just before
a step or evaluation (based on the value of the verbose
argument) in place of the default print. Define a callback
function with arguments t
, h
and y
where
t
is the time (beginning of a step or location of an
evaluation), h
is the step size (or NA
for an
evaluation) and y
is the state at the point of the step
or evaluation. Your callback function can do anything - you can
print to the screen (using cat
or message
), you
can store results using a closure and <<-
or you could
conditionally use a browser()
call to debug
interactively. However, it is not possible for the callback to
affect the solution short of throwing an error and interrupting
it. See the Examples for an example of use.
dopri_interpolate
which can be used to
efficiently sample from output of dopri
, and the package
vignette which shows in more detail how to solve delay
differential equations and to use compiled objective functions.
# The lorenz attractor: lorenz <- function(t, y, p) { sigma <- p[[1L]] R <- p[[2L]] b <- p[[3L]] c(sigma * (y[[2L]] - y[[1L]]), R * y[[1L]] - y[[2L]] - y[[1L]] * y[[3L]], -b * y[[3L]] + y[[1L]] * y[[2L]]) } p <- c(10, 28, 8 / 3) y0 <- c(10, 1, 1) tt <- seq(0, 100, length.out = 40000) y <- dde::dopri(y0, tt, lorenz, p, return_time = FALSE) plot(y[, c(1, 3)], type = "l", lwd = 0.5, col = "#00000066") # If we want to print progress as the integration progresses we can # use the verbose argument: y <- dde::dopri(y0, c(0, 0.1), lorenz, p, verbose = TRUE) # Or print the y values too using a callback: callback <- function(t, h, y) { message(sprintf("t: %f, h: %e, y: [%s]", t, h, paste(format(y, 5), collapse = ", "))) } y <- dde::dopri(y0, c(0, 0.1), lorenz, p, verbose = TRUE, callback = callback)
# The lorenz attractor: lorenz <- function(t, y, p) { sigma <- p[[1L]] R <- p[[2L]] b <- p[[3L]] c(sigma * (y[[2L]] - y[[1L]]), R * y[[1L]] - y[[2L]] - y[[1L]] * y[[3L]], -b * y[[3L]] + y[[1L]] * y[[2L]]) } p <- c(10, 28, 8 / 3) y0 <- c(10, 1, 1) tt <- seq(0, 100, length.out = 40000) y <- dde::dopri(y0, tt, lorenz, p, return_time = FALSE) plot(y[, c(1, 3)], type = "l", lwd = 0.5, col = "#00000066") # If we want to print progress as the integration progresses we can # use the verbose argument: y <- dde::dopri(y0, c(0, 0.1), lorenz, p, verbose = TRUE) # Or print the y values too using a callback: callback <- function(t, h, y) { message(sprintf("t: %f, h: %e, y: [%s]", t, h, paste(format(y, 5), collapse = ", "))) } y <- dde::dopri(y0, c(0, 0.1), lorenz, p, verbose = TRUE, callback = callback)
Interpolate the Dormand-Prince output after an integration. This only interpolates the core integration variables and not any additional output variables.
dopri_interpolate(h, t)
dopri_interpolate(h, t)
h |
The interpolation history. This can be the output
running |
t |
The times at which interpolated output is required. These times must fall within the included history (i.e., the times that the original simulation was run) or an error will be thrown. |
This decouples the integration of the equations and the generation of output; it is not necessary for use of the package, but may come in useful where you need to do (for example) root finding on the time course of a problem, or generate minimal output in some cases and interrogate the solution more deeply in others. See the examples and the package vignette for a full worked example.
Rich FitzJohn
# Here is the Lorenz attractor implemented as an R function lorenz <- function(t, y, p) { sigma <- p[[1L]] R <- p[[2L]] b <- p[[3L]] c(sigma * (y[[2L]] - y[[1L]]), R * y[[1L]] - y[[2L]] - y[[1L]] * y[[3L]], -b * y[[3L]] + y[[1L]] * y[[2L]]) } # Standard parameters and a reasonable starting point: p <- c(10, 28, 8 / 3) y0 <- c(10, 1, 1) # Run the integration for times [0, 50] and return minimal output, # but *do* record and return history. y <- dopri(y0, c(0, 50), lorenz, p, n_history = 5000, return_history = TRUE, return_time = FALSE, return_initial = FALSE, return_by_column = FALSE) # Very little output is returned (just 3 numbers being the final # state of the system), but the "history" attribute is fairly # large matrix of history information. It is not printed though # as its contents should not be relied on. What does matter is # the range of supported times printed (i.e., [0, 50]) and the # number of entries (~2000). y # Generate an interpolated set of variables using this; first for # 1000 steps over the full range: tt <- seq(0, 50, length.out = 1000) yy <- dopri_interpolate(y, tt) plot(yy[, c(1, 3)], type = "l") # Then for 50000 tt <- seq(0, 50, length.out = 50000) yy <- dopri_interpolate(y, tt) plot(yy[, c(1, 3)], type = "l")
# Here is the Lorenz attractor implemented as an R function lorenz <- function(t, y, p) { sigma <- p[[1L]] R <- p[[2L]] b <- p[[3L]] c(sigma * (y[[2L]] - y[[1L]]), R * y[[1L]] - y[[2L]] - y[[1L]] * y[[3L]], -b * y[[3L]] + y[[1L]] * y[[2L]]) } # Standard parameters and a reasonable starting point: p <- c(10, 28, 8 / 3) y0 <- c(10, 1, 1) # Run the integration for times [0, 50] and return minimal output, # but *do* record and return history. y <- dopri(y0, c(0, 50), lorenz, p, n_history = 5000, return_history = TRUE, return_time = FALSE, return_initial = FALSE, return_by_column = FALSE) # Very little output is returned (just 3 numbers being the final # state of the system), but the "history" attribute is fairly # large matrix of history information. It is not printed though # as its contents should not be relied on. What does matter is # the range of supported times printed (i.e., [0, 50]) and the # number of entries (~2000). y # Generate an interpolated set of variables using this; first for # 1000 steps over the full range: tt <- seq(0, 50, length.out = 1000) yy <- dopri_interpolate(y, tt) plot(yy[, c(1, 3)], type = "l") # Then for 50000 tt <- seq(0, 50, length.out = 50000) yy <- dopri_interpolate(y, tt) plot(yy[, c(1, 3)], type = "l")