
lsode
provides an interface to the FORTRAN ODE
solver of the same name, written by Alan C. Hindmarsh and Andrew
H. Sherman.
It combines parts of the code lsodar
and can thus find the root
of at least one of a set of constraint functions g(i) of the independent
and dependent variables. This can be used to stop the simulation or to
trigger events, i.e. a sudden change in one of the state variables.
The system of ODE's is written as an Rfunction 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.
lsode
is very similar to vode
, but uses a
fixed-step-interpolate method rather than the variable-coefficient
method in vode
. 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.lsode(y, times, func, parms, rtol = 1e-6, atol = 1e-6,
jacfunc = NULL, jactype = "fullint", mf = NULL, rootfunc = NULL,
verbose = FALSE, nroot = 0, 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,...)
y
has a name attribute, the names will be used to label the output
matrix.times
must be the initial time; if only one step is
to be taken; set times
= NULL
.
func
or
jacfunc
.y
. See details.y
. See details.NULL
, an Rfunction that computes the
Jacobian of the system of differential equations
$\partial\dot{y}_i/\partial y_j$, or
a string giving the name of a function or subroutine in
"fullint"
, "fullusr"
, "bandusr"
or
"bandint"
- either full or banded and estimated internally or
by user; overruled if mf
is not NULL
jactype
- provides more options than jactype
- see
details.NULL
, an Rfunction that computes the
function whose root has to be estimated or a string giving the name
of a function or subroutine in
diagnostiscs
of the integration - see details.rootfunc
is an R-function, the solver estimates the number
of roots.NULL
, then lsode
cannot integrate
past tcrit
. The FORTRAN routine lsode
overshoots its
targets (times points in the vector times
), and interpolates
values for the desired hmin
if you don't know why!hmax
is set to the largest difference in
times
, to avoid that the simulation possibly ignores
short-term events. If 0, no maximal size is specifFALSE
names of state variables are not
passed to function func
; this may speed up the simulation especially
for multi-D models.NULL
uses the default,
i.e. order 12 if implicit Adams method (meth = 1), order 5 if BDF
method (meth = 2). Reduce maxord to save storage space.func
and
jacfunc
. See package vignette "compiledCode"
.NULL
, the name of the initialisation function
(which initialises values of parameters), as provided in
"compiledCode"
.initfunc
is in the dll: the
parameters passed to the initialiser, to initialise the common
blocks (FORTRAN) or global variables (C, C++).func
and jacfunc
.func
and jacfunc
.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 checkenout
> 0: the names of output variables calculated in the
compiled function func
, present in the shared library.
These names will be used to label the output matrix.times
), max(times
)] is doneNULL
, the name of the forcing function
initialisation function, as provided in
forcings
has been given a
value.
See forcings compiledCode
.func
and
jacfunc
allowing this to be a generic function.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 `lsode'
returns with an unrecoverable error. If y
has a names
attribute, it will be used to label the columns of the output value.lsode
, whose
documentation should be consulted for details (it is included as
comments in the source file lsode
, 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 [object Object],[object Object],[object Object],[object Object]
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.
mf
is a positive two-digit integer, mf
=
(10*METH + MITER), where
[object Object],[object Object]
If MITER = 1 or 4, the user must supply a subroutine jacfunc
.
Inspection of the example below shows how to specify both a banded and
full Jacobian.
The input parameters rtol
, and atol
determine the
error control performed by the solver. 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 deSolve
package directory.
lsode
can find the root of at least one of a set of constraint functions
rootfunc
of the independent and dependent variables. It then returns the
solution at the root if that occurs sooner than the specified stop
condition, and otherwise returns the solution according the specified
stop condition.
Caution: Because of numerical errors in the function
rootfun
due to roundoff and integration error, lsode
may
return false roots, or return the same root at two or more
nearly equal values of time
.
rk
,rk4
andeuler
for
Runge-Kutta integrators.lsoda
,lsodes
,lsodar
,vode
,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.
## =======================================================================
## Example 1:
## Various ways to solve the same model.
## =======================================================================
## the model, 5 state variables
f1 <- function (t, y, parms) {
ydot <- vector(len = 5)
ydot[1] <- 0.1*y[1] -0.2*y[2]
ydot[2] <- -0.3*y[1] +0.1*y[2] -0.2*y[3]
ydot[3] <- -0.3*y[2] +0.1*y[3] -0.2*y[4]
ydot[4] <- -0.3*y[3] +0.1*y[4] -0.2*y[5]
ydot[5] <- -0.3*y[4] +0.1*y[5]
return(list(ydot))
}
## the Jacobian, written as a full matrix
fulljac <- function (t, y, parms) {
jac <- matrix(nrow = 5, ncol = 5, byrow = TRUE,
data = c(0.1, -0.2, 0 , 0 , 0 ,
-0.3, 0.1, -0.2, 0 , 0 ,
0 , -0.3, 0.1, -0.2, 0 ,
0 , 0 , -0.3, 0.1, -0.2,
0 , 0 , 0 , -0.3, 0.1))
return(jac)
}
## the Jacobian, written in banded form
bandjac <- function (t, y, parms) {
jac <- matrix(nrow = 3, ncol = 5, byrow = TRUE,
data = c( 0 , -0.2, -0.2, -0.2, -0.2,
0.1, 0.1, 0.1, 0.1, 0.1,
-0.3, -0.3, -0.3, -0.3, 0))
return(jac)
}
## initial conditions and output times
yini <- 1:5
times <- 1:20
## default: stiff method, internally generated, full Jacobian
out <- lsode(yini, times, f1, parms = 0, jactype = "fullint")
## stiff method, user-generated full Jacobian
out2 <- lsode(yini, times, f1, parms = 0, jactype = "fullusr",
jacfunc = fulljac)
## stiff method, internally-generated banded Jacobian
## one nonzero band above (up) and below(down) the diagonal
out3 <- lsode(yini, times, f1, parms = 0, jactype = "bandint",
bandup = 1, banddown = 1)
## stiff method, user-generated banded Jacobian
out4 <- lsode(yini, times, f1, parms = 0, jactype = "bandusr",
jacfunc = bandjac, bandup = 1, banddown = 1)
## non-stiff method
out5 <- lsode(yini, times, f1, parms = 0, mf = 10)
## =======================================================================
## Example 2:
## diffusion on a 2-D grid
## partially specified Jacobian
## =======================================================================
diffusion2D <- function(t, Y, par) {
y <- matrix(nrow = n, ncol = n, data = Y)
dY <- r*y # production
## diffusion in X-direction; boundaries = 0-concentration
Flux <- -Dx * rbind(y[1,],(y[2:n,]-y[1:(n-1),]),-y[n,])/dx
dY <- dY - (Flux[2:(n+1),]-Flux[1:n,])/dx
## diffusion in Y-direction
Flux <- -Dy * cbind(y[,1],(y[,2:n]-y[,1:(n-1)]),-y[,n])/dy
dY <- dY - (Flux[,2:(n+1)]-Flux[,1:n])/dy
return(list(as.vector(dY)))
}
## parameters
dy <- dx <- 1 # grid size
Dy <- Dx <- 1 # diffusion coeff, X- and Y-direction
r <- 0.025 # production rate
times <- c(0, 1)
n <- 50
y <- matrix(nrow = n, ncol = n, 0)
pa <- par(ask = FALSE)
## initial condition
for (i in 1:n) {
for (j in 1:n) {
dst <- (i - n/2)^2 + (j - n/2)^2
y[i, j] <- max(0, 1 - 1/(n*n) * (dst - n)^2)
}
}
filled.contour(y, color.palette = terrain.colors)
## =======================================================================
## jacfunc need not be estimated exactly
## a crude approximation, with a smaller bandwidth will do.
## Here the half-bandwidth 1 is used, whereas the true
## half-bandwidths are equal to n.
## This corresponds to ignoring the y-direction coupling in the ODEs.
## =======================================================================
print(system.time(
for (i in 1:20) {
out <- lsode(func = diffusion2D, y = as.vector(y), times = times,
parms = NULL, jactype = "bandint", bandup = 1, banddown = 1)
filled.contour(matrix(nrow = n, ncol = n, out[2,-1]), zlim = c(0,1),
color.palette = terrain.colors, main = i)
y <- out[2, -1]
}
))
par(ask = pa)
Run the code above in your browser using DataLab