dprjis Subroutine

public subroutine dprjis(Neq, Y, Yh, Nyh, Ewt, Rtem, Savr, S, Wk, Iwk, res, jac, adda)

DPRJIS is called to compute and process the matrix P = A - HEL(1)J, where J is an approximation to the Jacobian dr/dy, where r = g(t,y) - A(t,y)*s.

J is computed by columns, either by the user-supplied routine JAC if MITER = 1, or by finite differencing if MITER = 2.

J is stored in WK, rescaled, and ADDA is called to generate P.

The matrix P is subjected to LU decomposition in CDRV. P and its LU decomposition are stored separately in WK.

In addition to variables described previously, communication with DPRJIS uses the following:

Y

array containing predicted values on entry.

RTEM

work array of length N (ACOR in DSTODI).

SAVR

array containing r evaluated at predicted y. On output it contains the residual evaluated at current values of t and y.

S

array containing predicted values of dy/dt (SAVF in DSTODI).

WK

real work space for matrices. On output it contains P and its sparse LU decomposition. Storage of matrix elements starts at WK(3). WK also contains the following matrix-related data. WK(1) = SQRT(UROUND), used in numerical Jacobian increments.

IWK

integer work space for matrix-related data, assumed to be equivalenced to WK. In addition, WK(IPRSP) and IWK(IPISP) are assumed to have identical locations.

EL0

EL(1) (input).

IERPJ

output error flag (in COMMON). = 0 if no error. = 1 if zero pivot found in CDRV. = IRES (= 2 or 3) if RES returned IRES = 2 or 3. = -1 if insufficient storage for CDRV (should not occur). = -2 if other error found in CDRV (should not occur here).

JCUR

output flag = 1 to indicate that the Jacobian matrix (or approximation) is now current.

This routine also uses other variables in global structures.

Arguments

Type IntentOptional Attributes Name
integer, dimension(*) :: Neq
real(kind=dp), intent(inout), dimension(*) :: Y
real(kind=dp), intent(in), dimension(Nyh,*) :: Yh
integer, intent(in) :: Nyh
real(kind=dp), intent(in), dimension(*) :: Ewt
real(kind=dp), intent(inout), dimension(*) :: Rtem
real(kind=dp), dimension(*) :: Savr
real(kind=dp), dimension(*) :: S
real(kind=dp), intent(inout), dimension(*) :: Wk
integer, dimension(*) :: Iwk
real :: res
integer :: jac
real :: adda

Calls

proc~~dprjis~2~~CallsGraph proc~dprjis~2 M_odepack::dprjis proc~cdrv M_odepack::cdrv proc~dprjis~2->proc~cdrv