stodes {rootSolve}R Documentation

Steady-state solver for ordinary differential equations (ODE) with a sparse jacobian.

Description

Estimates the steady-state condition for a system of ordinary differential equations (ODE) in the form:

dy/dt = f(t,y)

and where the jacobian matrix df/dy has an arbitrary sparse structure.

Uses a newton-raphson method, implemented in Fortran.

The system of ODE's is written as an R function or defined in compiled code that has been dynamically loaded.

Usage

stodes(y, time=0, func, parms=NULL, 
rtol=1e-6, atol=1e-8, ctol=1e-8, sparsetype="sparseint",verbose=FALSE,
nnz=NULL, inz=NULL, lrw=NULL, ngp=NULL, positive=FALSE, maxiter=100, 
ynames=TRUE, dllname=NULL, initfunc=dllname, initpar=parms, 
rpar=NULL, ipar=NULL, nout=0, outnames = NULL, ...)

Arguments

y the initial guess of (state) values for the ode system, a vector. If y has a name attribute, the names will be used to label the output matrix.
time time for which steady-state is wanted; the default is time=0
func either a user-supplied function that computes the values of the derivatives in the ode system (the model definition) at time time, or a character string giving the name of a compiled function in a dynamically loaded shared library.
If func is a user-supplied function, it must be called as: yprime = func(t, y, parms). t is the time point at which the steady-state is wanted, 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 of parameters (which may have a names attribute).
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 (possibly with a names attribute) are global values that are required as output. If func is a string, then dllname must give the name of the shared library (without extension) which must be loaded before stodes() is called. see Details for more information.
parms other parameters passed to func
rtol relative error tolerance, either a scalar or a vector, one value for each y
atol absolute error tolerance, either a scalar or a vector, one value for each y
ctol if between two iterations, the maximal change in y is less than this amount, steady-state is reached
sparsetype the sparsity pattern, to date only "sparseint", sparse jacobian, estimated internally by stodes
verbose if TRUE: full output to the screen, e.g. will output the steady-state settings
nnz the number of nonzero elements in the sparse Jacobian (if this is unknown, use an estimate); If NULL, a guess will be made, and if not sufficient, stodes will return with a message indicating the size actually required.
If a solution is found, the minimal value of nnz actually required is returned by the solver (1st element of attribute dims)
inz (row,column) indices to the nonzero elements in the sparse Jacobian. If this is NULL, the sparsity will be determined by stodes
lrw the length of the work array of solver; due to the sparsicity, this cannot be readily predicted. If NULL, a guess will be made, and if not sufficient, stodes will return with a message indicating the size actually required. Therefore, some experimentation may be necessary to estimate the value of lrw
If a solution is found, the minimal value of lrw actually required is returned by the solver (3rd element of attribute dims)
ngp number of groups of independent state variables. Due to the sparsicity, this cannot be readily predicted. If NULL, a guess will be made, and if not sufficient, stodes will return with a message indicating the size actually required. Therefore, some experimentation may be necessary to estimate the value of ngp
If a solution is found, the minimal value of ngp actually required is returned by the solver (2nd element of attribute dims
positive either a logical or a vector with indices of the state variables that have to be non-negative; if TRUE, the state variables are forced to be non-negative numbers
maxiter maximal number of iterations during one call to the solver
ynames if FALSE: names of state variables are not passed to function func ; this may speed up the simulation especially for multi-D models
dllname a string giving the name of the shared library (without extension) that contains all the compiled function or subroutine definitions referred to in func.
initfunc if not NULL, the name of the initialisation function (which initialises values of parameters), as provided in ‘dllname’. See details.
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
ipar only when ‘dllname’ is specified: a vector with integer values passed to the dll-functions whose names are specified by func
nout only used if ‘dllname’ is specified: the number of output variables calculated in the compiled function func, present in the shared library
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

Details

The work is done by a Fortran 77 routine that implements the Newton-Raphson method.

stodes is to be used for problems, where the Jacobian has a sparse structure.

There are two choices for the sparsity specification, depending on whether inz is present.

Either way, the Jacobian itself is always generated by the solver (i.e. there is no provision to provide an analytic Jacobian).
This is done by perturbing simulataneously a combination of state variables that do not affect each other. This significantly reduces computing time. The number of groups with independent state variables can be given by ngp

The input parameters rtol, atol and ctol determine the error control performed by the solver. See help for stode for details.

Models may be defined in compiled C or Fortran code, as well as in R. See help for stode for details.

Value

A list containing

y A vector with the state variable values from the last iteration during estimation of steady-state condition of the system of equations. If y has a names attribute, it will be used to label the output values.
... the number of "global" values returned

The output will have the attribute steady, which returns TRUE, if steady-state has been reached and the attribute precis with an estimate of the precision attained during each iteration, the mean absolute rate of change (sum(abs(dy))/n).

Author(s)

Karline Soetaert <k.soetaert@nioo.knaw.nl>

References

For a description of the Newton-Raphson method, e.g.
Press, WH, Teukolsky, SA, Vetterling, WT, Flannery, BP, 1996. Numerical Recipes in FORTRAN. The Art of Scientific computing. 2nd edition. Cambridge University Press.

The algorithm uses linear algebra routines from the Yale sparse matrix package:
Eisenstat, S.C., Gursky, M.C., Schultz, M.H., Sherman, A.H., 1982. Yale Sparse Matrix Package. i. The symmetric codes. Int. J. Num. meth. Eng. 18, 1145-1151.

See Also

stode, steady-state solver with full or banded jacobian.

runsteady, for steady-state estimation by dynamically running till the derivatives become 0.

steady.band, for steady-state estimation, when the jacobian matrix is banded, and where the state variables need NOT to be rearranged

steady.1D, for steady-state estimation, when the jacobian matrix is banded, and where the state variables need to be rearranged

Examples

# 1000 simultaneous equations
model <- function (time,OC,parms,decay,ing)
{
 # model describing C in a sediment,
 # Upper boundary = imposed flux, lower boundary = zero-gradient
 Flux  <- v * c(OC[1] ,OC) +              # advection
          -Kz*diff(c(OC[1],OC,OC[N]))/dx  # diffusion;
 Flux[1]<- flux     # imposed flux

 # Rate of change= Flux gradient and first-order consumption
 dOC   <- -diff(Flux)/dx - decay*OC

 # Fraction of OC in first 5 layers is translocated to mean depth
 # (layer N/2)
 dOC[1:5]  <- dOC[1:5] - ing*OC[1:5]
 dOC[N/2]  <- dOC[N/2] + ing*sum(OC[1:5])
 list(dOC)
}
v    <- 0.1    # cm/yr
flux <- 10
dx   <- 0.01
N    <- 1000
dist <- seq(dx/2,by=dx,len=N)
Kz   <- 1  #bioturbation (diffusion), cm2/yr
ss   <- stodes(runif(N),func=model,parms=NULL,
               positive=TRUE, decay=5,ing=20,verbose=TRUE)

plot(ss$y[1:N],dist,ylim=rev(range(dist)),type="l",lwd=2,
     xlab="Nonlocal exchange",ylab="sediment depth",
     main="stodes, sparse jacobian")

[Package rootSolve version 1.3 Index]