vode {deSolve} | R Documentation |
Solver for Ordinary Differential Equations (ODE)
Description
Solves the initial value problem for stiff or nonstiff systems of ordinary differential equations (ODE) in the form:
dy/dt = f(t,y)
The R function vode
provides an interface to the FORTRAN ODE
solver of the same name, written by Peter N. Brown, Alan C. Hindmarsh
and George D. Byrne.
The system of ODE's is written as an R function or be defined in compiled code that has been dynamically loaded.
In contrast to lsoda
, the user has to specify whether or
not the problem is stiff and choose the appropriate solution method.
vode
is very similar to lsode
, but uses a
variable-coefficient method rather than the fixed-step-interpolate
methods in lsode
. In addition, in vode it is possible
to choose whether or not a copy of the Jacobian is saved for reuse in
the corrector iteration algorithm; In lsode
, a copy is not
kept.
Usage
vode(y, times, func, parms, rtol = 1e-6, atol = 1e-6,
jacfunc = NULL, jactype = "fullint", mf = NULL, verbose = FALSE,
tcrit = NULL, hmin = 0, hmax = NULL, hini = 0, ynames = TRUE,
maxord = NULL, bandup = NULL, banddown = NULL, maxsteps = 5000,
dllname = NULL, initfunc = dllname, initpar = parms, rpar = NULL,
ipar = NULL, nout = 0, outnames = NULL, forcings=NULL,
initforc = NULL, fcontrol=NULL, events=NULL, lags = NULL,...)
Arguments
y |
the initial (state) values for the ODE system. If |
times |
time sequence for which output is wanted; the first
value of |
func |
either an R-function that computes the values of the
derivatives in the ODE system (the model definition) at time
If The return value of If |
parms |
vector or list of parameters used in |
rtol |
relative error tolerance, either a scalar or an array as
long as |
atol |
absolute error tolerance, either a scalar or an array as
long as |
jacfunc |
if not In some circumstances, supplying
If the Jacobian is a full matrix, If the Jacobian is banded, |
jactype |
the structure of the Jacobian, one of
|
mf |
the "method flag" passed to function vode - overrules
|
verbose |
if TRUE: full output to the screen, e.g. will
print the |
tcrit |
if not |
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 largest difference in
|
hini |
initial step size to be attempted; if 0, the initial step size is determined by the solver. |
ynames |
logical; if |
maxord |
the maximum order to be allowed. |
bandup |
number of non-zero bands above the diagonal, in case the Jacobian is banded. |
banddown |
number of non-zero bands below the diagonal, in case the Jacobian is banded. |
maxsteps |
maximal number of steps per output interval taken by the solver. |
dllname |
a string giving the name of the shared library
(without extension) that contains all the compiled function or
subroutine definitions refered to in |
initfunc |
if not |
initpar |
only when ‘dllname’ is specified and an
initialisation function |
rpar |
only when ‘dllname’ is specified: a vector with
double precision values passed to the dll-functions whose names are
specified by |
ipar |
only when ‘dllname’ is specified: a vector with
integer values passed to the dll-functions whose names are specified
by |
nout |
only used if |
outnames |
only used if ‘dllname’ is specified and
|
forcings |
only used if ‘dllname’ is specified: a list with
the forcing function data sets, each present as a two-columned matrix,
with (time,value); interpolation outside the interval
[min( See forcings or package vignette |
initforc |
if not |
fcontrol |
A list of control parameters for the forcing functions.
forcings or package vignette |
events |
A matrix or data frame that specifies events, i.e. when the value of a state variable is suddenly changed. See events for more information. |
lags |
A list that specifies timelags, i.e. the number of steps that has to be kept. To be used for delay differential equations. See timelags, dede for more information. |
... |
additional arguments passed to |
Details
Before using the integrator vode
, the user has to decide
whether or not the problem is stiff.
If the problem is nonstiff, use method flag mf
= 10, which
selects a nonstiff (Adams) method, no Jacobian used.
If the problem is stiff, there are four standard choices which can be
specified with jactype
or mf
.
The options for jactype are
- jac = "fullint":
a full Jacobian, calculated internally by vode, corresponds to
mf
= 22,- jac = "fullusr":
a full Jacobian, specified by user function
jacfunc
, corresponds tomf
= 21,- jac = "bandusr":
a banded Jacobian, specified by user function
jacfunc
; the size of the bands specified bybandup
andbanddown
, corresponds tomf
= 24,- jac = "bandint":
a banded Jacobian, calculated by vode; the size of the bands specified by
bandup
andbanddown
, corresponds tomf
= 25.
More options are available when specifying mf directly.
The legal values of mf
are 10, 11, 12, 13, 14, 15, 20, 21, 22,
23, 24, 25, -11, -12, -14, -15, -21, -22, -24, -25.
mf
is a signed two-digit integer, mf = JSV*(10*METH +
MITER)
, where
- JSV = SIGN(mf)
indicates the Jacobian-saving strategy: JSV = 1 means a copy of the Jacobian is saved for reuse in the corrector iteration algorithm. JSV = -1 means a copy of the Jacobian is not saved.
- METH
indicates the basic linear multistep method: METH = 1 means the implicit Adams method. METH = 2 means the method based on backward differentiation formulas (BDF-s).
- MITER
indicates the corrector iteration method: MITER = 0 means functional iteration (no Jacobian matrix is involved).
MITER = 1 means chord iteration with a user-supplied full (NEQ by NEQ) Jacobian.
MITER = 2 means chord iteration with an internally generated (difference quotient) full Jacobian (using NEQ extra calls to
func
per df/dy value).MITER = 3 means chord iteration with an internally generated diagonal Jacobian approximation (using 1 extra call to
func
per df/dy evaluation).MITER = 4 means chord iteration with a user-supplied banded Jacobian.
MITER = 5 means chord iteration with an internally generated banded Jacobian (using ML+MU+1 extra calls to
func
per df/dy evaluation).
If MITER = 1 or 4, the user must supply a subroutine jacfunc
.
The example for integrator lsode
demonstrates how to
specify both a banded and full Jacobian.
The input parameters rtol
, and atol
determine the
error control performed by the solver. If the request for
precision exceeds the capabilities of the machine, vode will return an
error code. See lsoda
for details.
The diagnostics of the integration can be printed to screen
by calling diagnostics
. If verbose
= TRUE
,
the diagnostics will written to the screen at the end of the integration.
See vignette("deSolve") for an explanation of each element in the vectors containing the diagnostic properties and how to directly access them.
Models may be defined in compiled C or FORTRAN code, as well as
in an R-function. See package vignette "compiledCode"
for details.
More information about models defined in compiled code is in the package vignette ("compiledCode"); information about linking forcing functions to compiled code is in forcings.
Examples in both C and FORTRAN are in the ‘dynload’ subdirectory
of the deSolve
package directory.
Value
A matrix of class deSolve
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 FORTRAN routine ‘vode’
returns with an unrecoverable error. If y
has a names
attribute, it will be used to label the columns of the output value.
Note
From version 1.10.4, the default of atol
was changed from 1e-8 to 1e-6,
to be consistent with the other solvers.
Author(s)
Karline Soetaert <karline.soetaert@nioz.nl>
References
P. N. Brown, G. D. Byrne, and A. C. Hindmarsh, 1989. VODE: A Variable
Coefficient ODE Solver, SIAM J. Sci. Stat. Comput., 10, pp. 1038-1051.
Also, LLNL Report UCRL-98412, June 1988.
doi:10.1137/0910062
G. D. Byrne and A. C. Hindmarsh, 1975. A Polyalgorithm for the Numerical Solution of Ordinary Differential Equations. ACM Trans. Math. Software, 1, pp. 71-96. doi:10.1145/355626.355636
A. C. Hindmarsh and G. D. Byrne, 1977. EPISODE: An Effective Package for the Integration of Systems of Ordinary Differential Equations. LLNL Report UCID-30112, Rev. 1.
G. D. Byrne and A. C. Hindmarsh, 1976. EPISODEB: An Experimental Package for the Integration of Systems of Ordinary Differential Equations with Banded Jacobians. LLNL Report UCID-30132, April 1976.
A. C. Hindmarsh, 1983. ODEPACK, a Systematized Collection of ODE Solvers. in Scientific Computing, R. S. Stepleman et al., eds., North-Holland, Amsterdam, pp. 55-64.
K. R. Jackson and R. Sacks-Davis, 1980. An Alternative Implementation of Variable Step-Size Multistep Formulas for Stiff ODEs. ACM Trans. Math. Software, 6, pp. 295-318. doi:10.1145/355900.355903
Netlib: https://netlib.org
See Also
-
rk
, -
lsoda
,lsode
,lsodes
,lsodar
,daspk
for other 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.
Examples
## =======================================================================
## ex. 1
## The famous Lorenz equations: chaos in the earth's atmosphere
## Lorenz 1963. J. Atmos. Sci. 20, 130-141.
## =======================================================================
chaos <- function(t, state, parameters) {
with(as.list(c(state)), {
dx <- -8/3 * x + y * z
dy <- -10 * (y - z)
dz <- -x * y + 28 * y - z
list(c(dx, dy, dz))
})
}
state <- c(x = 1, y = 1, z = 1)
times <- seq(0, 100, 0.01)
out <- vode(state, times, chaos, 0)
plot(out, type = "l") # all versus time
plot(out[,"x"], out[,"y"], type = "l", main = "Lorenz butterfly",
xlab = "x", ylab = "y")
## =======================================================================
## ex. 2
## SCOC model, in FORTRAN - to see the FORTRAN code:
## browseURL(paste(system.file(package="deSolve"),
## "/doc/examples/dynload/scoc.f",sep=""))
## example from Soetaert and Herman, 2009, chapter 3. (simplified)
## =======================================================================
## Forcing function data
Flux <- matrix(ncol = 2, byrow = TRUE, data = c(
1, 0.654, 11, 0.167, 21, 0.060, 41, 0.070, 73, 0.277, 83, 0.186,
93, 0.140,103, 0.255, 113, 0.231,123, 0.309,133, 1.127,143, 1.923,
153,1.091,163, 1.001, 173, 1.691,183, 1.404,194, 1.226,204, 0.767,
214,0.893,224, 0.737, 234, 0.772,244, 0.726,254, 0.624,264, 0.439,
274,0.168,284, 0.280, 294, 0.202,304, 0.193,315, 0.286,325, 0.599,
335,1.889,345, 0.996, 355, 0.681,365, 1.135))
parms <- c(k = 0.01)
meanDepo <- mean(approx(Flux[,1], Flux[,2], xout = seq(1, 365, by = 1))$y)
Yini <- c(y = as.double(meanDepo/parms))
times <- 1:365
out <- vode(Yini, times, func = "scocder",
parms = parms, dllname = "deSolve",
initforc = "scocforc", forcings = Flux,
initfunc = "scocpar", nout = 2,
outnames = c("Mineralisation", "Depo"))
matplot(out[,1], out[,c("Depo", "Mineralisation")],
type = "l", col = c("red", "blue"), xlab = "time", ylab = "Depo")
## Constant interpolation of forcing function - left side of interval
fcontrol <- list(method = "constant")
out2 <- vode(Yini, times, func = "scocder",
parms = parms, dllname = "deSolve",
initforc = "scocforc", forcings = Flux, fcontrol = fcontrol,
initfunc = "scocpar", nout = 2,
outnames = c("Mineralisation", "Depo"))
matplot(out2[,1], out2[,c("Depo", "Mineralisation")],
type = "l", col = c("red", "blue"), xlab = "time", ylab = "Depo")
## Constant interpolation of forcing function - middle of interval
fcontrol <- list(method = "constant", f = 0.5)
out3 <- vode(Yini, times, func = "scocder",
parms = parms, dllname = "deSolve",
initforc = "scocforc", forcings = Flux, fcontrol = fcontrol,
initfunc = "scocpar", nout = 2,
outnames = c("Mineralisation", "Depo"))
matplot(out3[,1], out3[,c("Depo", "Mineralisation")],
type = "l", col = c("red", "blue"), xlab = "time", ylab = "Depo")
plot(out, out2, out3)