dvode_module Module

Module for the refactored DVODE.

See also

  • The original DVODE.f (1989-2002) can be found here

History

  • Translated into modern Fortran by Jacob Williams.

Todo

Replace GOTO statements


Uses

  • module~~dvode_module~~UsesGraph module~dvode_module dvode_module iso_fortran_env iso_fortran_env module~dvode_module->iso_fortran_env module~dvode_blas_module dvode_blas_module module~dvode_module->module~dvode_blas_module module~dvode_kinds_module dvode_kinds_module module~dvode_module->module~dvode_kinds_module module~dvode_linpack_module dvode_linpack_module module~dvode_module->module~dvode_linpack_module module~dvode_blas_module->module~dvode_kinds_module module~dvode_kinds_module->iso_fortran_env module~dvode_linpack_module->module~dvode_blas_module module~dvode_linpack_module->module~dvode_kinds_module

Variables

Type Visibility Attributes Name Initial
integer, private, parameter :: wp = dvode_wp

real kind to use

real(kind=wp), private, parameter :: zero = 0.0_wp
real(kind=wp), private, parameter :: epmach = epsilon(1.0_wp)

machine epsilon

integer, private, parameter :: iumach = output_unit

standard output unit number

integer, private, parameter, dimension(2) :: mord = [12, 5]

Abstract Interfaces

abstract interface

  • private function f_dvnorm(me, n, v, w)

    Arguments

    Type IntentOptional Attributes Name
    class(dvode_t), intent(inout) :: me
    integer, intent(in) :: n
    real(kind=wp), intent(in) :: v(n)
    real(kind=wp), intent(in) :: w(n)

    Return Value real(kind=wp)

abstract interface

  • private subroutine f_func(me, neq, t, y, ydot)

    Arguments

    Type IntentOptional Attributes Name
    class(dvode_t), intent(inout) :: me
    integer :: neq
    real(kind=wp) :: t
    real(kind=wp) :: y(neq)
    real(kind=wp) :: ydot(neq)

abstract interface

  • private subroutine f_jac(me, neq, t, y, ml, mu, pd, nrowpd)

    Arguments

    Type IntentOptional Attributes Name
    class(dvode_t), intent(inout) :: me
    integer :: neq
    real(kind=wp) :: t
    real(kind=wp) :: y(neq)
    integer :: ml
    integer :: mu
    real(kind=wp) :: pd(nrowpd,neq)
    integer :: nrowpd

abstract interface

  • private subroutine f_dewset(me, n, itol, rtol, atol, ycur, ewt)

    Arguments

    Type IntentOptional Attributes Name
    class(dvode_t), intent(inout) :: me
    integer, intent(in) :: n
    integer, intent(in) :: itol
    real(kind=wp), intent(in) :: rtol(*)
    real(kind=wp), intent(in) :: atol(*)
    real(kind=wp), intent(in) :: ycur(n)
    real(kind=wp), intent(out) :: ewt(n)

Derived Types

type, public ::  dvode_data_t

Internal variables for DVODE. variables which are communicated between subroutines in the dvode package, or which are to be saved between calls to dvode.

Components

Type Visibility Attributes Name Initial
real(kind=wp), private :: acnrm = zero

weighted r.m.s. norm of accumulated correction vectors.

real(kind=wp), private :: ccmxj = zero

threshhold on drc for updating the jacobian. (see drc.)

real(kind=wp), private :: conp = zero

the saved value of tq(5).

real(kind=wp), private :: crate = zero

estimated corrector convergence rate constant.

real(kind=wp), private :: drc = zero

relative change in h*rl1 since last dvjac call.

real(kind=wp), private :: el(13) = zero

real array of integration coefficients. see dvset.

real(kind=wp), private :: eta = zero

saved tentative ratio of new to old h.

real(kind=wp), private :: etamax = zero

saved maximum value of eta to be allowed.

real(kind=wp), private :: h = zero

the step size.

real(kind=wp), private :: hmin = zero

the minimum absolute value of the step size h to be used.

real(kind=wp), private :: hmxi = zero

inverse of the maximum absolute value of h to be used. hmxi = 0.0 is allowed and corresponds to an infinite hmax.

real(kind=wp), private :: hnew = zero

the step size to be attempted on the next step.

real(kind=wp), private :: hscal = zero

step size h used in scaling of nordsieck array yh in dvjust. (if iord = +1, dvjust assumes that hscal = tau(1).) see references 1 and 2 for details.

real(kind=wp), private :: prl1 = zero

the saved value of rl1.

real(kind=wp), private :: rc = zero

ratio of current h*rl1 to value on last dvjac call.

real(kind=wp), private :: rl1 = zero

the reciprocal of the coefficient el(1).

real(kind=wp), private :: tau(13) = zero

real vector of past nq step sizes, length 13.

real(kind=wp), private :: tq(5) = zero

a real vector of length 5 in which dvset stores constants used for the convergence test, the error test, and the selection of h at a new order.

real(kind=wp), private :: tn = zero

the independent variable, updated on each step taken.

real(kind=wp), private :: uround = zero

the machine unit roundoff. the smallest positive real number such that 1.0 + uround /= 1.0

integer, private :: icf = 0

integer flag for convergence failure in dvnlsd:

Read more…
integer, private :: init = 0

saved integer flag indicating whether initialization of the problem has been done (init = 1) or not.

integer, private :: ipup = 0

saved flag to signal updating of newton matrix.

Read more…
integer, private :: jcur = 0

output flag from dvjac showing jacobian status:

Read more…
integer, private :: jstart = 0

integer flag used as input to dvstep:

Read more…
integer, private :: jsv = 0

integer flag for jacobian saving, = sign(mf).

integer, private :: kflag = 0

a completion code from dvstep with the following meanings:

Read more…
integer, private :: kuth = 0

input flag to dvstep showing whether h was reduced by the driver. kuth = 1 if h was reduced, = 0 otherwise.

integer, private :: l = 0

integer variable, nq + 1, current order plus one.

integer, private :: lmax = 0

maxord + 1 (used for dimensioning).

integer, private :: lyh = 0

saved integer pointers to segments of rwork and iwork.

integer, private :: lewt = 0

saved integer pointers to segments of rwork and iwork.

integer, private :: lacor = 0

saved integer pointers to segments of rwork and iwork.

integer, private :: lsavf = 0

saved integer pointers to segments of rwork and iwork.

integer, private :: lwm = 0

saved integer pointers to segments of rwork and iwork.

integer, private :: liwm = 0

saved integer pointers to segments of rwork and iwork.

integer, private :: locjs = 0

a pointer to the saved jacobian, whose storage starts at wm(locjs), if jsv = 1.

integer, private :: maxord = 0

the maximum order of integration method to be allowed.

integer, private :: meth = 0

the method flags. see mf.

integer, private :: miter = 0

the method flags. see mf.

integer, private :: msbj = 0

the maximum number of steps between j evaluations, = 50.

integer, private :: mxhnil = 0

saved value of optional input mxhnil.

integer, private :: mxstep = 0

saved value of optional input mxstep.

integer, private :: n = 0

the number of first-order odes, = neq.

integer, private :: newh = 0

saved integer to flag change of h.

integer, private :: newq = 0

the method order to be used on the next step.

integer, private :: nhnil = 0

saved counter for occurrences of t + h = t.

integer, private :: nq = 0

integer variable, the current integration method order.

integer, private :: nqnyh = 0

saved value of nq*nyh.

integer, private :: nqwait = 0

a counter controlling the frequency of order changes. an order change is about to be considered if nqwait = 1.

integer, private :: nslj = 0

the number of steps taken as of the last jacobian update.

integer, private :: nslp = 0

saved value of nst as of last newton matrix update.

integer, private :: nyh = 0

saved value of the initial value of neq.

real(kind=wp), private :: hu = zero

the step size in t last used.

integer, private :: ncfn = 0

number of nonlinear convergence failures so far.

integer, private :: netf = 0

the number of error test failures of the integrator so far.

integer, private :: nfe = 0

the number of f evaluations for the problem so far.

integer, private :: nje = 0

the number of jacobian evaluations so far.

integer, private :: nlu = 0

the number of matrix lu decompositions so far.

integer, private :: nni = 0

number of nonlinear iterations so far.

integer, private :: nqu = 0

the method order last used.

integer, private :: nst = 0

the number of steps taken for the problem so far.

type, public ::  dvode_t

Main DVODE class.

Read more…

Components

Type Visibility Attributes Name Initial
type(dvode_data_t), private :: dat

internal data (formerly in common blocks)

real(kind=wp), private :: etaq = zero
real(kind=wp), private :: etaqm1 = zero
integer, private :: lunit = -1

logical unit number for messages. the default is obtained by a call to iumach (may be machine-dependent).

integer, private :: mesflg = 1

print control flag:

Read more…
procedure(f_func), private, pointer :: f => null()
procedure(f_jac), private, pointer :: jac => null()
procedure(f_dewset), private, pointer :: dewset => null()
procedure(f_dvnorm), private, pointer :: dvnorm => null()

Type-Bound Procedures

procedure, public :: initialize ../../

this routine must be called first before using the solver.

procedure, public :: solve => dvode ../../

main solver routine.

procedure, public :: xsetun

set the logical unit number, lun, for output of messages from dvode, if the default is not desired. the default value of lun is output_unit. this call may be made at any time and will take effect immediately.

procedure, public :: xsetf

set a flag to control the printing of messages by dvode:

Read more…
procedure, public :: dvsrco

saves and restores the contents of the internal variables used by dvode. dvsrco is useful if one is interrupting a run and restarting later, or alternating between two or more problems solved with dvode.

procedure, public :: dvindy

provide derivatives of y, of various orders, at a specified point t, if desired. it may be called only after a successful return from dvode.

procedure, private :: dvhin
procedure, private :: dvstep
procedure, private :: dvset
procedure, private :: dvjust
procedure, private :: dvnlsd
procedure, private :: dvjac
procedure, private :: dvsol
procedure, private :: ixsav
procedure, private :: xerrwd

Functions

private function dvnorm_default(me, n, v, w)

weighted root-mean-square vector norm.

Read more…

Arguments

Type IntentOptional Attributes Name
class(dvode_t), intent(inout) :: me
integer, intent(in) :: n
real(kind=wp), intent(in) :: v(n)
real(kind=wp), intent(in) :: w(n)

Return Value real(kind=wp)

private function ixsav(me, ipar, ivalue, iset)

save and recall error message control parameters.

Read more…

Arguments

Type IntentOptional Attributes Name
class(dvode_t), intent(inout) :: me
integer, intent(in) :: ipar

parameter indicator (1 for lunit, 2 for mesflg).

integer :: ivalue

the value to be set for the parameter, if iset = .true.

logical, intent(in) :: iset

logical flag to indicate whether to read or write. if iset = .true., the parameter will be given the value ivalue. if iset = .false., the parameter will be unchanged, and ivalue is a dummy argument.

Return Value integer


Subroutines

private subroutine initialize(me, f, jac, dewset, dvnorm)

Set the function pointers. This must be called before dvode is called.

Read more…

Arguments

Type IntentOptional Attributes Name
class(dvode_t), intent(inout) :: me
procedure(f_func) :: f

the name of the user-supplied subroutine defining the ode system. the system must be put in the first-order form dy/dt = f(t,y), where f is a vector-valued function of the scalar t and the vector y. subroutine f is to compute the function f. it is to have the form subroutine f (neq, t, y, ydot) real(wp) t, y(neq), ydot(neq) where neq, t, and y are input, and the array ydot = f(t,y) is output. y and ydot are arrays of length neq. subroutine f should not alter y(1),...,y(neq). f must be declared external in the calling program.

Read more…
procedure(f_jac), optional :: jac

the name of the user-supplied routine (miter = 1 or 4) to compute the jacobian matrix, df/dy, as a function of the scalar t and the vector y. it is to have the form f_jac, where neq, t, y, ml, mu, and nrowpd are input and the array pd is to be loaded with partial derivatives (elements of the jacobian matrix) in the output. pd must be given a first dimension of nrowpd. t and y have the same meaning as in subroutine f.

Read more…
procedure(f_dewset), optional :: dewset

the following subroutine is called just before each internal integration step, and sets the array of error weights, ewt, as described under itol/rtol/atol above: subroutine dewset (neq, itol, rtol, atol, ycur, ewt) where neq, itol, rtol, and atol are as in the dvode call sequence, ycur contains the current dependent variable vector, and ewt is the array of weights set by dewset.

Read more…
procedure(f_dvnorm), optional :: dvnorm

the following is a real function routine which computes the weighted root-mean-square norm of a vector v: d = dvnorm (n, v, w) where:

Read more…

public subroutine dvode(me, neq, y, t, tout, itol, rtol, atol, itask, istate, iopt, rwork, lrw, iwork, liw, mf)

dvode: variable-coefficient ordinary differential equation solver, with fixed-leading-coefficient implementation.

Read more…

Arguments

Type IntentOptional Attributes Name
class(dvode_t), intent(inout) :: me
integer, intent(in) :: neq

the size of the ode system (number of first order ordinary differential equations). neq may not be increased during the problem, but can be decreased (with istate = 3 in the input).

real(kind=wp), intent(inout) :: y(*)

a real array for the vector of dependent variables, of length neq or more. used for both input and output on the first call (istate = 1), and only for output on other calls. on the first call, y must contain the vector of initial values. in the output, y contains the computed solution evaluated at t. if desired, the y array may be used for other purposes between calls to the solver.

Read more…
real(kind=wp), intent(inout) :: t

the independent variable. in the input, t is used only on the first call, as the initial point of the integration. in the output, after each call, t is the value at which a computed solution y is evaluated (usually the same as tout). on an error return, t is the farthest point reached.

real(kind=wp), intent(in) :: tout

the next value of t at which a computed solution is desired.

Read more…
integer, intent(in) :: itol

an indicator for the type of error control. 1 or 2 according as atol is a scalar or array. see description under atol.

real(kind=wp), intent(in) :: rtol(*)

a relative error tolerance parameter, either a scalar or an array of length neq. see description under atol.

real(kind=wp), intent(in) :: atol(*)

an absolute error tolerance parameter, either a scalar or an array of length neq.

Read more…
integer, intent(in) :: itask

an index specifying the task to be performed. input only. itask has the following values and meanings:

Read more…
integer, intent(inout) :: istate

an index used for input and output to specify the the state of the calculation.

Read more…
integer, intent(in) :: iopt

an integer flag to specify whether or not any optional input is being used on this call. the optional input is listed separately below.

Read more…
real(kind=wp) :: rwork(lrw)

a real working array. the length of rwork must be at least 20 + nyh*(maxord + 1) + 3*neq + lwm where:

Read more…
integer, intent(in) :: lrw

the length of the array rwork, as declared by the user. (this will be checked by the solver.)

integer :: iwork(liw)

an integer work array. the length of iwork must be at least:

Read more…
integer, intent(in) :: liw

the length of the array iwork, as declared by the user. (this will be checked by the solver.)

integer, intent(in) :: mf

method flag. standard values are:

Read more…

private subroutine dvhin(me, n, t0, y0, ydot, tout, uround, ewt, itol, atol, y, temp, h0, niter, ier)

this routine computes the step size, h0, to be attempted on the first step, when the user has not supplied a value for this.

Read more…

Arguments

Type IntentOptional Attributes Name
class(dvode_t), intent(inout) :: me
integer, intent(in) :: n

size of ode system

real(kind=wp), intent(in) :: t0

initial value of independent variable

real(kind=wp), intent(in) :: y0(*)

vector of initial conditions

real(kind=wp), intent(in) :: ydot(*)

vector of initial first derivatives

real(kind=wp), intent(in) :: tout

first output value of independent variable

real(kind=wp), intent(in) :: uround

machine unit roundoff

real(kind=wp), intent(in) :: ewt(*)

error weights and tolerance parameters as described in the driver routine

integer, intent(in) :: itol

error weights and tolerance parameters as described in the driver routine

real(kind=wp), intent(in) :: atol(*)

error weights and tolerance parameters as described in the driver routine

real(kind=wp), intent(inout) :: y(n)

work array of length n

real(kind=wp), intent(inout) :: temp(n)

work array of length n

real(kind=wp), intent(out) :: h0

step size to be attempted

integer, intent(out) :: niter

number of iterations (and of f evaluations) to compute h0

integer, intent(out) :: ier

the error flag, returned with the value:

Read more…

private subroutine dvindy(me, t, k, yh, ldyh, dky, iflag)

dvindy computes interpolated values of the k-th derivative of the dependent variable vector y, and stores it in dky. this routine is called within the package with k = 0 and t = tout, but may also be called by the user for any k up to the current order. (see detailed instructions in the usage documentation.)

Read more…

Arguments

Type IntentOptional Attributes Name
class(dvode_t), intent(inout) :: me
real(kind=wp), intent(in) :: t

value of independent variable where answers are desired (normally the same as the t last returned by dvode). for valid results, t must lie between tcur - hu and tcur. (see optional output for tcur and hu.)

integer, intent(in) :: k

integer order of the derivative desired. k must satisfy 0 <= k <= nqcur, where nqcur is the current order (see optional output). the capability corresponding to k = 0, i.e. computing y(t), is already provided by dvode directly. since nqcur >= 1, the first derivative dy/dt is always available with dvindy.

real(kind=wp), intent(in) :: yh(ldyh,*)

the history array yh

integer, intent(in) :: ldyh

column length of yh, equal to the initial value of neq.

real(kind=wp), intent(out) :: dky(*)

a real array of length neq containing the computed value of the k-th derivative of y(t).

integer, intent(out) :: iflag

integer flag, returned as 0 if k and t were legal, -1 if k was illegal, and -2 if t was illegal. on an error return, a message is also written.

private subroutine dvstep(me, y, yh, ldyh, yh1, ewt, savf, vsav, acor, wm, iwm)

dvstep performs one step of the integration of an initial value problem for a system of ordinary differential equations.

Read more…

Arguments

Type IntentOptional Attributes Name
class(dvode_t), intent(inout) :: me
real(kind=wp), intent(inout) :: y(*)

an array of length n used for the dependent variable vector.

real(kind=wp), intent(inout) :: yh(ldyh,*)

an ldyh by lmax array containing the dependent variables and their approximate scaled derivatives, where lmax = maxord + 1. yh(i,j+1) contains the approximate j-th derivative of y(i), scaled by h**j/factorial(j) (j = 0,1,...,nq). on entry for the first step, the first two columns of yh must be set from the initial values.

integer, intent(in) :: ldyh

a constant integer >= n, the first dimension of yh. n is the number of odes in the system.

real(kind=wp), intent(inout) :: yh1(*)

a one-dimensional array occupying the same space as yh.

real(kind=wp), intent(in) :: ewt(*)

an array of length n containing multiplicative weights for local error measurements. local errors in y(i) are compared to 1.0/ewt(i) in various error tests.

real(kind=wp) :: savf(*)

an array of working storage, of length n. also used for input of yh(*,maxord+2) when jstart = -1 and maxord < the current order nq.

real(kind=wp) :: vsav(*)

a work array of length n passed to subroutine dvnlsd.

real(kind=wp) :: acor(*)

a work array of length n, used for the accumulated corrections. on a successful return, acor(i) contains the estimated one-step local error in y(i).

real(kind=wp) :: wm(*)

real work array associated with matrix operations in dvnlsd.

integer :: iwm(*)

integer work array associated with matrix operations in dvnlsd.

private subroutine dvset(me)

dvset is called by dvstep and sets coefficients for use there.

Read more…

Arguments

Type IntentOptional Attributes Name
class(dvode_t), intent(inout) :: me

private subroutine dvjust(me, yh, ldyh, iord)

this subroutine adjusts the yh array on reduction of order, and also when the order is increased for the stiff option (meth = 2).

Arguments

Type IntentOptional Attributes Name
class(dvode_t), intent(inout) :: me
real(kind=wp), intent(inout) :: yh(ldyh,*)
integer, intent(in) :: ldyh

leading dimension of yh

integer, intent(in) :: iord

an integer flag used when meth = 2 to indicate an order increase (iord = +1) or an order decrease (iord = -1).

private subroutine dvnlsd(me, y, yh, ldyh, vsav, savf, ewt, acor, iwm, wm, nflag)

subroutine dvnlsd is a nonlinear system solver, which uses functional iteration or a chord (modified newton) method. for the chord method direct linear algebraic system solvers are used. subroutine dvnlsd then handles the corrector phase of this integration package.

Arguments

Type IntentOptional Attributes Name
class(dvode_t), intent(inout) :: me
real(kind=wp), intent(inout) :: y(*)

the dependent variable, a vector of length n

real(kind=wp), intent(inout) :: yh(ldyh,*)

the nordsieck (taylor) array, ldyh by lmax, input and output. on input, it contains predicted values.

integer, intent(in) :: ldyh

a constant >= n, the first dimension of yh

real(kind=wp) :: vsav(*)

unused work array.

real(kind=wp) :: savf(*)

a work array of length n.

real(kind=wp), intent(in) :: ewt(*)

an error weight vector of length n

real(kind=wp), intent(inout) :: acor(*)

a work array of length n, used for the accumulated corrections to the predicted y vector.

integer, intent(inout) :: iwm(*)

integer work array associated with matrix operations in chord iteration (miter /= 0).

real(kind=wp), intent(inout) :: wm(*)

real work array associated with matrix operations in chord iteration (miter /= 0).

integer, intent(inout) :: nflag

input/output flag, with values and meanings as follows:

Read more…

private subroutine dvjac(me, y, yh, ldyh, ewt, ftem, savf, wm, iwm, ierpj)

dvjac is called by dvnlsd to compute and process the matrix p = i - h*rl1*j , where j is an approximation to the jacobian.

Read more…

Arguments

Type IntentOptional Attributes Name
class(dvode_t), intent(inout) :: me
real(kind=wp), intent(inout) :: y(*)

vector containing predicted values on entry.

real(kind=wp), intent(in) :: yh(ldyh,*)

the nordsieck array, an ldyh by lmax array.

integer, intent(in) :: ldyh

a constant >= n, the first dimension of yh.

real(kind=wp), intent(in) :: ewt(*)

an error weight vector of length n.

real(kind=wp), intent(in) :: ftem(*)
real(kind=wp), intent(in) :: savf(*)

array containing f evaluated at predicted y.

real(kind=wp), intent(inout) :: wm(*)

real work space for matrices. in the output, it contains the inverse diagonal matrix if miter = 3 and the lu decomposition of p if miter is 1, 2 , 4, or 5. storage of matrix elements starts at wm(3). storage of the saved jacobian starts at wm(locjs). wm also contains the following matrix-related data: wm(1) = sqrt(uround), used in numerical jacobian step. wm(2) = h*rl1, saved for later use if miter = 3.

integer, intent(inout) :: iwm(*)

integer work space containing pivot information, starting at iwm(31), if miter is 1, 2, 4, or 5. iwm also contains band parameters ml = iwm(1) and mu = iwm(2) if miter is 4 or 5.

integer, intent(out) :: ierpj

output error flag, = 0 if no trouble, 1 if the p matrix is found to be singular.

private subroutine dacopy(nrow, ncol, a, nrowa, b, nrowb)

this routine copies one rectangular array, a, to another, b, where a and b may have different row dimensions, nrowa and nrowb. the data copied consists of nrow rows and ncol columns.

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nrow
integer, intent(in) :: ncol
real(kind=wp), intent(in) :: a(nrowa,ncol)
integer, intent(in) :: nrowa
real(kind=wp), intent(out) :: b(nrowb,ncol)
integer, intent(in) :: nrowb

private subroutine dvsol(me, wm, iwm, x, iersl)

This routine manages the solution of the linear system arising from a chord iteration. it is called if miter /= 0:

Read more…

Arguments

Type IntentOptional Attributes Name
class(dvode_t), intent(inout) :: me
real(kind=wp), intent(inout) :: wm(*)

real work space containing the inverse diagonal matrix if miter = 3 and the lu decomposition of the matrix otherwise. storage of matrix elements starts at wm(3). wm also contains the following matrix-related data: wm(1) = sqrt(uround) (not used here), wm(2) = hrl1, the previous value of h*rl1, used if miter = 3.

integer :: iwm(*)

integer work space containing pivot information, starting at iwm(31), if miter is 1, 2, 4, or 5. iwm also contains band parameters ml = iwm(1) and mu = iwm(2) if miter is 4 or 5.

real(kind=wp), intent(inout) :: x(*)

the right-hand side vector on input, and the solution vector on output, of length n.

integer, intent(out) :: iersl

output flag. iersl = 0 if no trouble occurred. iersl = 1 if a singular matrix arose with miter = 3.

private subroutine dvsrco(me, sav, job)

this routine saves or restores (depending on job) the contents of the dvode internal variables.

Arguments

Type IntentOptional Attributes Name
class(dvode_t), intent(inout) :: me
type(dvode_data_t), intent(inout) :: sav
integer, intent(in) :: job

flag indicating to save or restore the data:

Read more…

private subroutine dewset_default(me, n, itol, rtol, atol, ycur, ewt)

Set error weight vector.

Read more…

Arguments

Type IntentOptional Attributes Name
class(dvode_t), intent(inout) :: me
integer, intent(in) :: n
integer, intent(in) :: itol
real(kind=wp), intent(in) :: rtol(*)
real(kind=wp), intent(in) :: atol(*)
real(kind=wp), intent(in) :: ycur(n)
real(kind=wp), intent(out) :: ewt(n)

private subroutine xerrwd(me, msg, nmes, nerr, level, ni, i1, i2, nr, r1, r2)

write error message with values.

Read more…

Arguments

Type IntentOptional Attributes Name
class(dvode_t), intent(inout) :: me
character(len=*), intent(in) :: msg

the message (character array).

integer, intent(in) :: nmes

the length of msg (number of characters).

integer, intent(in) :: nerr

the error number (not used).

integer, intent(in) :: level

the error level:

Read more…
integer, intent(in) :: ni

number of integers (0, 1, or 2) to be printed with message.

integer, intent(in) :: i1

integer to be printed, depending on ni.

integer, intent(in) :: i2

integer to be printed, depending on ni.

integer, intent(in) :: nr

number of reals (0, 1, or 2) to be printed with message.

real(kind=wp), intent(in) :: r1

real to be printed, depending on nr.

real(kind=wp), intent(in) :: r2

real to be printed, depending on nr.

private subroutine xsetf(me, mflag)

reset the error print control flag.

Read more…

Arguments

Type IntentOptional Attributes Name
class(dvode_t), intent(inout) :: me
integer :: mflag

private subroutine xsetun(me, lun)

reset the logical unit number for error messages.

Read more…

Arguments

Type IntentOptional Attributes Name
class(dvode_t), intent(inout) :: me
integer :: lun