rk {deSolve} | R Documentation |
Solving initial value problems for non-stiff systems of first-order ordinary differential equations (ODEs).
The R function rk
is a top-level function that provides
interfaces to a collection of common explicit one-step solvers of the
Runge-Kutta family with fixed or variable time steps.
The system of ODE's is written as an R function (which may, of
course, use .C
, .Fortran
,
.Call
, etc., to call foreign code) or be defined in
compiled code that has been dynamically loaded. A vector of
parameters is passed to the ODEs, so the solver may be used as part of
a modeling package for ODEs, or for parameter estimation using any
appropriate modeling tool for non-linear models in R such as
optim
, nls
, nlm
or
nlme
rk(y, times, func, parms, rtol = 1e-6, atol = 1e-6, verbose = FALSE, tcrit = NULL, hmin = 0, hmax = NULL, hini = hmax, ynames = TRUE, method = rkMethod("rk45dp7", ... ), maxsteps = 5000, dllname = NULL, initfunc = dllname, initpar = parms, rpar = NULL, ipar = NULL, nout = 0, outnames = NULL, ...)
y |
the initial (state) values for the ODE system. If y
has a name attribute, the names will be used to label the output
matrix.
|
times |
times at which explicit estimates for y are
desired. The first value in times must be the initial time.
|
func |
either an R-function that computes the values of the
derivatives in the ODE system (the model definition) at time
t, or a character string giving the name of a compiled function in a
dynamically loaded shared library.
If func is an R-function, it must be defined as:
yprime = func(t, y, parms, ...) . t is the current
time point in the integration, y is the current estimate of
the variables in the ODE system. If the initial values y has
a names attribute, the names will be available inside func .
parms is a vector or list of parameters; ... (optional) are
any other arguments passed to the function.
The return value of func should be a list, whose first
element is a vector containing the derivatives of y with
respect to time , and whose next elements are global values
that are required at each point in times . The derivatives
should be specified in the same order as the state variables y .
If func is
a string, then dllname must give the name of the shared
library (without extension) which must be loaded before
vode() is called. See package vignette "compiledCode"
for more details.
|
parms |
vector or list of parameters used in func |
rtol |
relative error tolerance, either a scalar or an array as
long as y . Only applicable to methods with variable time
step, see details.
|
atol |
absolute error tolerance, either a scalar or an array as
long as y . Only applicable to methods with variable time
step, see details.
|
tcrit |
if not NULL , then rk cannot integrate past
tcrit . The solver routines may overshoot their targets (times
points in the vector times ), and interpolates values for the
desired time points. If there is a time beyond which integration
should not proceed (perhaps because of a singularity), that should
be provided in tcrit .
|
verbose |
a logical value that, when TRUE, triggers more verbose output from the ODE solver. |
hmin |
an optional minimum value of the integration stepsize. In
special situations this parameter may speed up computations with the
cost of precision. Don't use hmin if you don't know why!
|
hmax |
an optional maximum value of the integration stepsize. If
not specified, hmax is set to the maximum of hini and
the largest difference in times , to avoid that the simulation
possibly ignores short-term events. If 0, no maximal size is
specified. Note that hmin and hmax are ignored by
fixed step methods like "rk4" or "euler" .
|
hini |
initial step size to be attempted; if 0, the initial step
size is determined automatically by solvers with flexible time step.
Setting hini = 0 for fixed step methods forces setting of
internal time steps identically to external time steps provided by
times .
|
ynames |
if FALSE : names of state variables are not passed
to function func ; this may speed up the simulation especially
for large models.
|
method |
the integrator to use. This can either be a string
constant naming one of the pre-defined methods or a call to function
rkMethod specifying a user-defined method. The most
common methods are the fixed-step methods "euler" ,
"rk2" , "rk4" or the variable step methods
"rk23bs" , "rk34f" , "rk45f" or
"rk45dp7" .
|
maxsteps |
average maximal number of steps per output interval
taken by the solver. This argument is defined such to ensure
compatibility with the Livermore-solvers, but the maximum number of
steps in total is calculated as length(times) * maxsteps .
|
dllname |
a string giving the name of the shared library
(without extension) that contains all the compiled function or
subroutine definitions refered to in func and
jacfunc . See package vignette "compiledCode" .
|
initfunc |
if not NULL , the name of the initialisation function
(which initialises values of parameters), as provided in
‘dllname’. See package vignette "compiledCode" .
|
initpar |
only when ‘dllname’ is specified and an
initialisation function initfunc is in the dll: the
parameters passed to the initialiser, to initialise the common
blocks (fortran) or global variables (C, C++).
|
rpar |
only when ‘dllname’ is specified: a vector with
double precision values passed to the dll-functions whose names are
specified by func and jacfunc .
|
ipar |
only when ‘dllname’ is specified: a vector with
integer values passed to the dll-functions whose names are specified
by func and jacfunc .
|
nout |
only used if dllname is specified and the model is
defined in compiled code: the number of output variables calculated
in the compiled function func , present in the shared
library. Note: it is not automatically checked whether this is
indeed the number of output variables calculed in the dll - you have
to perform this check in the code. See package vignette
"compiledCode" .
|
outnames |
only used if ‘dllname’ is specified and
nout > 0: the names of output variables calculated in the
compiled function func , present in the shared library.
|
... |
additional arguments passed to func allowing this
to be a generic function.
|
The Runge-Kutta solvers are primarily provided for didactic reasons.
For most practical cases, solvers of the Livermore family (ODEPACK;
lsoda
, lsode
, lsodes
,
lsodar
, vode
, daspk
) are
superior and more thoroughly tested; in most cases they are also
faster).
In addition to this, some of the Livermore solvers are also suitable for stiff ODEs, differential algebraic equations (DAEs), or partial differential equations (PDEs).
Function rk
is a generalized implementation that can be used to
evaluate different solvers of the Runge-Kutta family. A pre-defined
set of common method parameters is in function rkMethod
which also allows to supply user-defined Butcher tables.
The input parameters rtol
, and atol
determine the error
control performed by the solver. The solver will control the vector
of estimated local errors in y, according to an inequality of
the form max-norm of ( e/ewt ) <= 1, where
ewt is a vector of positive error weights. The values of
rtol
and atol
should all be non-negative. The form of
ewt is:
rtol * abs(y) + atol
where multiplication of two vectors is element-by-element.
Models can be defined in R as a user-supplied
R-function, that must be called as: yprime = func(t, y,
parms)
. t
is the current time point in the integration,
y
is the current estimate of the variables in the ODE system.
The return value of func
should be a list, whose first element
is a vector containing the derivatives of y with respect to time, and
whose second element contains output variables that are required at
each point in time. Examples are given below.
A matrix with up to as many rows as elements in times
and
as many columns as elements in y
plus the number of "global"
values returned in the next elements of the return from func
,
plus and additional column for the time value. There will be a row
for each element in times
unless the solver returns with an
unrecoverable error. If y
has a names attribute, it will be
used to label the columns of the output value.
The output will have the attributes istate
, and rstate
,
two vectors with several useful elements, whose interpretation is
compatible with lsoda
:
el 1: |
0 for normal return, -2 means excess accuracy requested (tolerances too small). |
el 12: |
the number of steps taken for the problem so far. |
el 13: |
the number of function evaluations for the problem so far. |
el 15: |
the order of the method. |
Arguments rpar
and ipar
are provided for compatibility
with lsoda
, but not yet thoroughly tested.
Thomas Petzoldt thomas.petzoldt@tu-dresden.de
Butcher, J. C. (1987) The numerical analysis of ordinary differential equations, Runge-Kutta and general linear methods, Wiley, Chichester and New York.
Engeln-Muellges, G. and Reutter, F. (1996) Numerik Algorithmen: Entscheidungshilfe zur Auswahl und Nutzung. VDI Verlag, Duesseldorf.
Press, W. H., Teukolsky, S. A., Vetterling, W. T. and Flannery, B. P. (2007) Numerical Recipes in C. Cambridge University Press.
rkMethod
for a list of available Runge-Kutta
parameter sets,
rk4
and euler
for special
versions without interpolation (and less overhead),
lsoda
, lsode
,
lsodes
, lsodar
, vode
,
daspk
for solvers of the Livermore family,
ode
for a general interface to most of the ODE solvers,
ode.band
for solving models with a banded
Jacobian,
ode.1D
for integrating 1-D models,
ode.2D
for integrating 2-D models,
ode.3D
for integrating 3-D models,
diagnostics
to print diagnostic messages.
## ========================================================= ## Example: Resource-producer-consumer Lotka-Volterra model ## ========================================================= ## Note: ## parameters are a list, names accessible via "with" function ## (see also ode and lsoda examples) lvmodel <- function(t, x, parms) { S <- x[1] # substrate P <- x[2] # producer K <- x[3] # consumer with(parms, { import <- approx(signal$times, signal$import, t)$y dS <- import - b * S * P + g * K dP <- c * S * P - d * K * P dK <- e * P * K - f * K res <- c(dS, dP, dK) list(res) }) } ## vector of timesteps times <- seq(0, 100, length = 101) ## external signal with rectangle impulse signal <- as.data.frame(list(times = times, import = rep(0,length(times)))) signal$import[signal$times >= 10 & signal$times <= 11] <- 0.2 ## Parameters for steady state conditions parms <- list(b = 0.0, c = 0.1, d = 0.1, e = 0.1, f = 0.1, g = 0.0) ## Start values for steady state y <- xstart <- c(S = 1, P = 1, K = 1) system.time({ ## Euler method out1 <- as.data.frame(rk(xstart, times, lvmodel, parms, hini = 0.1, method = "euler")) ## classical Runge-Kutta 4th order out2 <- as.data.frame(rk(xstart, times, lvmodel, parms, hini = 1, method = "rk4")) ## Dormand-Prince method of order 5(4) out3 <- as.data.frame(rk(xstart, times, lvmodel, parms, hmax = 1, method = "rk45dp7")) }) mf <- par(mfrow = c(2,2)) plot (out1$time, out1$S, type = "l", ylab = "Substrate") lines(out2$time, out2$S, col = "red", lty = "dotted", lwd = 2) lines(out3$time, out3$S, col = "green", lty = "dotted") plot (out1$time, out1$P, type = "l", ylab = "Producer") lines(out2$time, out2$P, col = "red", lty = "dotted") lines(out3$time, out3$P, col = "green", lty = "dotted") plot (out1$time, out1$K, type = "l", ylab = "Consumer") lines(out2$time, out2$K, col = "red", lty = "dotted", lwd = 2) lines(out3$time, out3$K, col = "green", lty = "dotted") plot (out1$P, out1$K, type = "l", xlab = "Producer", ylab = "Consumer") lines(out2$P, out2$K, col = "red", lty = "dotted", lwd = 2) lines(out3$P, out3$K, col = "green", lty = "dotted") legend("center", legend = c("euler", "rk4", "rk45dp7"), lty = c(1, 3, 3), lwd = c(1, 2, 1), col = c("black", "red", "green")) par(mfrow = mf)