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
time sequence for which output is wanted; the first
value of times must be the initial time; if only one step is
to be taken; set times = NULL.
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:
func <- function(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
must 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
lsodes() is called. See package vignette "compiledCode"
for more details.
parms
vector or list of parameters used in func or
jacfunc.
rtol
relative error tolerance, either a scalar or an array as
long as y. See details.
atol
absolute error tolerance, either a scalar or an array as
long as y. See details.
jacvec
if not NULL, an R function that computes a
column of the Jacobian of the system of differential equations
dydot(i)/dy(j),
or a string giving the name of a function or
subroutine in ‘dllname’ that computes the column of the
Jacobian (see vignette "compiledCode" for more about this option).
The R
calling sequence for jacvec is identical to that of
func, but with extra parameter j, denoting the column
number. Thus, jacvec should be called as: jacvec =
func(t, y, j, parms) and jacvec should return a vector
containing column j of the Jacobian, i.e. its i-th value is
dydot(i)/dy(j).
If this function is absent, lsodes will
generate the Jacobian by differences.
sparsetype
the sparsity structure of the Jacobian, one of
"sparseint" or "sparseusr", "sparsejan", ...,
The sparsity can be estimated internally by lsodes (first option)
or given by the user (last two). See details.
nnz
the number of nonzero elements in the sparse Jacobian (if
this is unknown, use an estimate).
inz
if sparsetype equal to "sparseusr", a two-columned matrix
with the (row, column) indices to the nonzero elements in the sparse
Jacobian. If sparsetype = "sparsejan", a vector with the elements
ian followed by he elements jan as used in the lsodes code. See details.
In all other cases, ignored.
rootfunc
if not NULL, an R function that computes the
function whose root has to be estimated or a string giving the name
of a function or subroutine in ‘dllname’ that computes the root
function. The R calling sequence for rootfunc is identical
to that of func. rootfunc should return a vector with
the function values whose root is sought.
verbose
if TRUE: full output to the screen, e.g. will
print the diagnostiscs of the integration - see details.
nroot
only used if ‘dllname’ is specified: the number of
constraint functions whose roots are desired during the integration;
if rootfunc is an R-function, the solver estimates the number
of roots.
tcrit
if not NULL, then lsodes cannot integrate
past tcrit. The FORTRAN routine lsodes overshoots its
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.
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
times, to avoid that the simulation possibly ignores
short-term events. If 0, no maximal size is specified.
hini
initial step size to be attempted; if 0, the initial step
size is determined by the solver.
ynames
logical, if FALSE names of state variables are
not passed to function func; this may speed up the simulation
especially for multi-D models.
maxord
the maximum order to be allowed. 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.
maxsteps
maximal number of steps per output interval taken by the
solver.
lrw
the length of the real work array rwork; due to the
sparsicity, this cannot be readily predicted. If NULL, a
guess will be made, and if not sufficient, lsodes will return
with a message indicating the size of rwork actually required.
Therefore, some experimentation may be necessary to estimate the
value of lrw.
For instance, if you get the error:
DLSODES- RWORK length is insufficient to proceed.
Length needed is .ge. LENRW (=I1), exceeds LRW (=I2)
In above message, I1 = 27627 I2 = 25932
set lrw equal to 27627 or a higher value
liw
the length of the integer work array iwork; due to the
sparsicity, this cannot be readily predicted. If NULL, a guess will
be made, and if not sufficient, lsodes will return with a
message indicating the size of iwork actually required. Therefore,
some experimentation may be necessary to estimate the value of
liw.
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 calculated 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.
These names will be used to label the output matrix.
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(times), max(times)] is done by taking the value at
the closest data extreme.
See forcings or package vignette "compiledCode".
initforc
if not NULL, the name of the forcing function
initialisation function, as provided in
‘dllname’. It MUST be present if forcings has been given a
value.
See forcings or package vignette "compiledCode".
fcontrol
A list of control parameters for the forcing functions.
See forcings or vignette compiledCode.
events
A list 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 func and
jacfunc allowing this to be a generic function.
Details
The work is done by the FORTRAN subroutine lsodes, whose
documentation should be consulted for details (it is included as
comments in the source file ‘src/opkdmain.f’). The implementation
is based on the November, 2003 version of lsodes, from Netlib.
lsodes is applied for stiff problems, where the Jacobian has a
sparse structure.
There are several choices depending on whether jacvec
is specified and depending on the setting of sparsetype.
If function jacvec is present, then it should return the j-th
column of the Jacobian matrix.
There are also several choices for the sparsity specification, selected by
argument sparsetype.
sparsetype = "sparseint". The sparsity is estimated
by the solver, based on numerical differences.
In this case, it is advisable to provide an estimate of the number
of non-zero elements in the Jacobian (nnz).
This value can be approximate; upon return the number of nonzero
elements actually required will be known (1st element of attribute
dims).
In this case, inz need not be specified.
sparsetype = "sparseusr". The sparsity is determined by
the user. In this case, inz should be a matrix, containing indices
(row, column) to the nonzero elements in the Jacobian matrix.
The number of nonzeros nnz will be set equal to the number of rows
in inz.
sparsetype = "sparsejan". The sparsity is also determined by
the user.
In this case, inz should be a vector, containting the ian and
jan elements of the sparse storage format, as used in the sparse solver.
Elements of ian should be the first n+1 elements of this vector, and
contain the starting locations in jan of columns 1.. n.
jan contains the row indices of the nonzero locations of
the Jacobian, reading in columnwise order.
The number of nonzeros nnz will be set equal to the length of inz - (n+1).
sparsetype = "1D", "2D", "3D".
The sparsity is estimated by the solver, based on numerical differences.
Assumes finite differences in a 1D, 2D or 3D regular grid - used by
functions ode.1D, ode.2D, ode.3D.
Similar are "2Dmap", and "3Dmap", which also include a
mapping variable (passed in nnz).
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 ‘doc/examples/dynload’ subdirectory
of the deSolve package directory.
lsodes 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, lsodes may
return false roots, or return the same root at two or more
nearly equal values of time.
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 ‘lsodes’
returns with an unrecoverable error. If y has a names
attribute, it will be used to label the columns of the output value.
Author(s)
Karline Soetaert <karline.soetaert@nioz.nl>
References
Alan C. Hindmarsh, ODEPACK, A Systematized Collection of ODE Solvers,
in Scientific Computing, R. S. Stepleman et al. (Eds.), North-Holland,
Amsterdam, 1983, pp. 55-64.
S. C. Eisenstat, M. C. Gursky, M. H. Schultz, and A. H. Sherman, Yale
Sparse Matrix Package: I. The Symmetric Codes,
Int. J. Num. Meth. Eng., 18 (1982), pp. 1145-1151.
S. C. Eisenstat, M. C. Gursky, M. H. Schultz, and A. H. Sherman, Yale
Sparse Matrix Package: II. The Nonsymmetric Codes, Research Report
No. 114, Dept. of Computer Sciences, Yale University, 1977.
See Also
rk,
rk4 and euler for
Runge-Kutta integrators.
lsoda, lsode,
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.
Examples
## Various ways to solve the same model.
## =======================================================================
## The example from lsodes source code
## A chemical model
## =======================================================================
n <- 12
y <- rep(1, n)
dy <- rep(0, n)
times <- c(0, 0.1*(10^(0:4)))
rtol <- 1.0e-4
atol <- 1.0e-6
parms <- c(rk1 = 0.1, rk2 = 10.0, rk3 = 50.0, rk4 = 2.5, rk5 = 0.1,
rk6 = 10.0, rk7 = 50.0, rk8 = 2.5, rk9 = 50.0, rk10 = 5.0,
rk11 = 50.0, rk12 = 50.0,rk13 = 50.0, rk14 = 30.0,
rk15 = 100.0,rk16 = 2.5, rk17 = 100.0,rk18 = 2.5,
rk19 = 50.0, rk20 = 50.0)
#
chemistry <- function (time, Y, pars) {
with (as.list(pars), {
dy[1] <- -rk1 *Y[1]
dy[2] <- rk1 *Y[1] + rk11*rk14*Y[4] + rk19*rk14*Y[5] -
rk3 *Y[2]*Y[3] - rk15*Y[2]*Y[12] - rk2*Y[2]
dy[3] <- rk2 *Y[2] - rk5 *Y[3] - rk3*Y[2]*Y[3] -
rk7*Y[10]*Y[3] + rk11*rk14*Y[4] + rk12*rk14*Y[6]
dy[4] <- rk3 *Y[2]*Y[3] - rk11*rk14*Y[4] - rk4*Y[4]
dy[5] <- rk15*Y[2]*Y[12] - rk19*rk14*Y[5] - rk16*Y[5]
dy[6] <- rk7 *Y[10]*Y[3] - rk12*rk14*Y[6] - rk8*Y[6]
dy[7] <- rk17*Y[10]*Y[12] - rk20*rk14*Y[7] - rk18*Y[7]
dy[8] <- rk9 *Y[10] - rk13*rk14*Y[8] - rk10*Y[8]
dy[9] <- rk4 *Y[4] + rk16*Y[5] + rk8*Y[6] +
rk18*Y[7]
dy[10] <- rk5 *Y[3] + rk12*rk14*Y[6] + rk20*rk14*Y[7] +
rk13*rk14*Y[8] - rk7 *Y[10]*Y[3] - rk17*Y[10]*Y[12] -
rk6 *Y[10] - rk9*Y[10]
dy[11] <- rk10*Y[8]
dy[12] <- rk6 *Y[10] + rk19*rk14*Y[5] + rk20*rk14*Y[7] -
rk15*Y[2]*Y[12] - rk17*Y[10]*Y[12]
return(list(dy))
})
}
## =======================================================================
## application 1. lsodes estimates the structure of the Jacobian
## and calculates the Jacobian by differences
## =======================================================================
out <- lsodes(func = chemistry, y = y, parms = parms, times = times,
atol = atol, rtol = rtol, verbose = TRUE)
## =======================================================================
## application 2. the structure of the Jacobian is input
## lsodes calculates the Jacobian by differences
## this is not so efficient...
## =======================================================================
## elements of Jacobian that are not zero
nonzero <- matrix(nc = 2, byrow = TRUE, data = c(
1, 1, 2, 1, # influence of sp1 on rate of change of others
2, 2, 3, 2, 4, 2, 5, 2, 12, 2,
2, 3, 3, 3, 4, 3, 6, 3, 10, 3,
2, 4, 3, 4, 4, 4, 9, 4, # d (dyi)/dy4
2, 5, 5, 5, 9, 5, 12, 5,
3, 6, 6, 6, 9, 6, 10, 6,
7, 7, 9, 7, 10, 7, 12, 7,
8, 8, 10, 8, 11, 8,
3,10, 6,10, 7,10, 8,10, 10,10, 12,10,
2,12, 5,12, 7,12, 10,12, 12,12)
)
## when run, the default length of rwork is too small
## lsodes will tell the length actually needed
# out2 <- lsodes(func = chemistry, y = y, parms = parms, times = times,
# inz = nonzero, atol = atol,rtol = rtol) #gives warning
out2 <- lsodes(func = chemistry, y = y, parms = parms, times = times,
sparsetype = "sparseusr", inz = nonzero,
atol = atol, rtol = rtol, verbose = TRUE, lrw = 353)
## =======================================================================
## application 3. lsodes estimates the structure of the Jacobian
## the Jacobian (vector) function is input
## =======================================================================
chemjac <- function (time, Y, j, pars) {
with (as.list(pars), {
PDJ <- rep(0,n)
if (j == 1){
PDJ[1] <- -rk1
PDJ[2] <- rk1
} else if (j == 2) {
PDJ[2] <- -rk3*Y[3] - rk15*Y[12] - rk2
PDJ[3] <- rk2 - rk3*Y[3]
PDJ[4] <- rk3*Y[3]
PDJ[5] <- rk15*Y[12]
PDJ[12] <- -rk15*Y[12]
} else if (j == 3) {
PDJ[2] <- -rk3*Y[2]
PDJ[3] <- -rk5 - rk3*Y[2] - rk7*Y[10]
PDJ[4] <- rk3*Y[2]
PDJ[6] <- rk7*Y[10]
PDJ[10] <- rk5 - rk7*Y[10]
} else if (j == 4) {
PDJ[2] <- rk11*rk14
PDJ[3] <- rk11*rk14
PDJ[4] <- -rk11*rk14 - rk4
PDJ[9] <- rk4
} else if (j == 5) {
PDJ[2] <- rk19*rk14
PDJ[5] <- -rk19*rk14 - rk16
PDJ[9] <- rk16
PDJ[12] <- rk19*rk14
} else if (j == 6) {
PDJ[3] <- rk12*rk14
PDJ[6] <- -rk12*rk14 - rk8
PDJ[9] <- rk8
PDJ[10] <- rk12*rk14
} else if (j == 7) {
PDJ[7] <- -rk20*rk14 - rk18
PDJ[9] <- rk18
PDJ[10] <- rk20*rk14
PDJ[12] <- rk20*rk14
} else if (j == 8) {
PDJ[8] <- -rk13*rk14 - rk10
PDJ[10] <- rk13*rk14
PDJ[11] <- rk10
} else if (j == 10) {
PDJ[3] <- -rk7*Y[3]
PDJ[6] <- rk7*Y[3]
PDJ[7] <- rk17*Y[12]
PDJ[8] <- rk9
PDJ[10] <- -rk7*Y[3] - rk17*Y[12] - rk6 - rk9
PDJ[12] <- rk6 - rk17*Y[12]
} else if (j == 12) {
PDJ[2] <- -rk15*Y[2]
PDJ[5] <- rk15*Y[2]
PDJ[7] <- rk17*Y[10]
PDJ[10] <- -rk17*Y[10]
PDJ[12] <- -rk15*Y[2] - rk17*Y[10]
}
return(PDJ)
})
}
out3 <- lsodes(func = chemistry, y = y, parms = parms, times = times,
jacvec = chemjac, atol = atol, rtol = rtol)
## =======================================================================
## application 4. The structure of the Jacobian (nonzero elements) AND
## the Jacobian (vector) function is input
## =======================================================================
out4 <- lsodes(func = chemistry, y = y, parms = parms, times = times,
lrw = 351, sparsetype = "sparseusr", inz = nonzero,
jacvec = chemjac, atol = atol, rtol = rtol,
verbose = TRUE)
# The sparsejan variant
# note: errors in inz may cause R to break, so this is not without danger...
# out5 <- lsodes(func = chemistry, y = y, parms = parms, times = times,
# jacvec = chemjac, atol = atol, rtol = rtol, sparsetype = "sparsejan",
# inz = c(1,3,8,13,17,21,25,29,32,32,38,38,43, # ian
# 1,2, 2,3,4,5,12, 2,3,4,6,10, 2,3,4,9, 2,5,9,12, 3,6,9,10, # jan
# 7,9,10,12, 8,10,11, 3,6,7,8,10,12, 2,5,7,10,12), lrw = 343)