stodes: Steady-state solver for ordinary differential equations (ODE)...

Description Usage Arguments Details Value Author(s) References See Also Examples

View source: R/stodes.R

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

1
2
3
4
5
6
7
8
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, 
       forcings = NULL, initforc = NULL, fcontrol = NULL,
       spmethod = "yale", control = NULL, times = time, ...)

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, times

time for which steady-state is wanted; the default is times=0. (note- since version 1.7, 'times' has been added as an alias to 'time').

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.

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 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 assumed to be reached.

sparsetype

the sparsity structure of the Jacobian, one of "sparseint", "sparseusr", "sparsejan", "sparsereturn", ..., The sparsity can be estimated internally by stodes (first and last option) or given by the user (other two). See details. Note: setting sparsetype equal to "sparsereturn" will not solve for steady-state but solely return the ian and jan vectors.

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

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 the elements jan as used in the stodes code. See details. In all other cases, ignored. If inz is NULL, the sparsity will be determined by stodes.

lrw

the length of the work array of the 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 that lrw should be increased. 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).

In case of an error induced by a too small value of lrw, its value can be assessed by the attributes()$dims value.

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.

spmethod

the sparse method to be used, one of "yale", "ilut", "ilutp". The default uses the yale sparse matrix solver; the other use preconditioned GMRES (generalised minimum residual method) solvers from FORTRAN package sparsekit. ilut stands for incomplete LU factorisation with trheshold (or tolerances, droptol); the "p" iin ilutp stands for pivoting.

control

only used if spmethod not equal to "yale", a list with the control options of the preconditioned solvers. The default is list( droptol = 1e-3, permtol = 1e-3, fillin = 10, lenplufac = 2). droptol is the tolerance in ilut, ilutp to decide when to drop a value. permtol is used in ilutp, to decide whether or not to permute variables. See Saad 1994, the manual of sparskit and Saad 2003, chapter 10 for details.

forcings

only used if ‘dllname’ is specified: a vector with the forcing function values, or 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.

This feature is here for compatibility with models defined in compiled code from package deSolve; see deSolve's 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 deSolve's package vignette "compiledCode".

fcontrol

A list of control parameters for the forcing functions. See deSolve's package vignette "compiledCode".

...

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 several choices for the sparsity specification, selected by argument sparsetype.

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 package vignette for details on how to write models in compiled code.

When the spmethod equals ilut or ilutp, a number of parameters can be specified in argument control. They are:

fillin, the fill-in parameter. Each row of L and each row of U will have a maximum of lfil elements (excluding the diagonal element). lfil must be >= 0.

droptol, sets the threshold for dropping small terms in the factorization.

When ilutp is chosen the following arguments can also be specified:

permtol = tolerance ratio used to determne whether or not to permute two columns. At step i columns i and j are permuted when abs(a(i,j))*permtol .gt. abs(a(i,i)) [0 –> never permute; good values 0.1 to 0.01]

lenplufac = sets the working array - increase its value if a warning.

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).

In case the argument sparsetype is set to "sparsereturn", then two extra attributes will be returned, i.e. ian and jan. These can then be used to speed up subsequent calculations - see last example.

Author(s)

Karline Soetaert <karline.soetaert@nioz.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.

When spmethod = "yale" then 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.

else the functions ilut and ilutp from sparsekit package are used:

Yousef Saad, 1994. SPARSKIT: a basic tool kit for sparse matrix computations. VERSION 2

Yousef Saad, 2003. Iterative methods for Sparse Linear Systems. Society for Industrial and Applied Mathematics.

See Also

steady, for a general interface to most of the steady-state solvers

steady.band, to find the steady-state of ODE models with a banded Jacobian

steady.1D, steady.2D, steady.3D, steady-state solvers for 1-D, 2-D and 3-D partial differential equations.

stode, iterative steady-state solver for ODEs with full or banded Jacobian.

runsteady, steady-state solver by dynamically running to steady-state

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
## =======================================================================
## 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
print(system.time(
  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")

# the size of lrw is in the attributes()$dims vector.     
attributes(ss)     

## =======================================================================
## deriving sparsity structure and speeding up calculations
## =======================================================================

sparse <- stodes(runif(N), func = model, parms = NULL,
                 sparsetype = "sparsereturn", decay = 5, ing = 20)

ian <- attributes(sparse)$ian
jan <- attributes(sparse)$jan
nnz <- length(jan)
inz <- c(ian, jan)

print(system.time(
s2   <- stodes(runif(N), func = model, parms = NULL, positive = TRUE,
               sparsetype = "sparsejan", inz = inz, decay = 5, ing = 20)
))

# Can also be used with steady.1D, by setting jactype = "sparsejan".
# The advantage is this allows easy plotting...

print(system.time(
s2b   <- steady.1D(runif(N), func = model, parms = NULL, method = "stodes", 
                   nspec = 1, jactype = "sparsejan", inz = inz, 
               decay = 5, ing = 20, verbose = FALSE)
))

plot(s2b, grid = dist, xyswap = TRUE, type = "l", lwd = 2,
     xlab = "Nonlocal exchange", ylab = "sediment depth",
     main = "steady 1-D, sparse jacobian imposed")

rootSolve documentation built on Sept. 23, 2021, 3 a.m.