*DECK DSTODE
SUBROUTINE DSTODE (NEQ, Y, YH, NYH, YH1, EWT, SAVF, ACOR,
1 WM, IWM, F, JAC, PJAC, SLVS)
C***BEGIN PROLOGUE DSTODE
C***SUBSIDIARY
C***PURPOSE Performs one step of an ODEPACK integration.
C***TYPE DOUBLE PRECISION (SSTODE-S, DSTODE-D)
C***AUTHOR Hindmarsh, Alan C., (LLNL)
C***DESCRIPTION
C
C DSTODE performs one step of the integration of an initial value
C problem for a system of ordinary differential equations.
C Note: DSTODE is independent of the value of the iteration method
C indicator MITER, when this is .ne. 0, and hence is independent
C of the type of chord method used, or the Jacobian structure.
C Communication with DSTODE is done with the following variables:
C
C NEQ = integer array containing problem size in NEQ(1), and
C passed as the NEQ argument in all calls to F and JAC.
C Y = an array of length .ge. N used as the Y argument in
C all calls to F and JAC.
C YH = an NYH by LMAX array containing the dependent variables
C and their approximate scaled derivatives, where
C LMAX = MAXORD + 1. YH(i,j+1) contains the approximate
C j-th derivative of y(i), scaled by h**j/factorial(j)
C (j = 0,1,...,NQ). on entry for the first step, the first
C two columns of YH must be set from the initial values.
C NYH = a constant integer .ge. N, the first dimension of YH.
C YH1 = a one-dimensional array occupying the same space as YH.
C EWT = an array of length N containing multiplicative weights
C for local error measurements. Local errors in Y(i) are
C compared to 1.0/EWT(i) in various error tests.
C SAVF = an array of working storage, of length N.
C Also used for input of YH(*,MAXORD+2) when JSTART = -1
C and MAXORD .lt. the current order NQ.
C ACOR = a work array of length N, used for the accumulated
C corrections. On a successful return, ACOR(i) contains
C the estimated one-step local error in Y(i).
C WM,IWM = real and integer work arrays associated with matrix
C operations in chord iteration (MITER .ne. 0).
C PJAC = name of routine to evaluate and preprocess Jacobian matrix
C and P = I - h*el0*JAC, if a chord method is being used.
C SLVS = name of routine to solve linear system in chord iteration.
C CCMAX = maximum relative change in h*el0 before PJAC is called.
C H = the step size to be attempted on the next step.
C H is altered by the error control algorithm during the
C problem. H can be either positive or negative, but its
C sign must remain constant throughout the problem.
C HMIN = the minimum absolute value of the step size h to be used.
C HMXI = inverse of the maximum absolute value of h to be used.
C HMXI = 0.0 is allowed and corresponds to an infinite hmax.
C HMIN and HMXI may be changed at any time, but will not
C take effect until the next change of h is considered.
C TN = the independent variable. TN is updated on each step taken.
C JSTART = an integer used for input only, with the following
C values and meanings:
C 0 perform the first step.
C .gt.0 take a new step continuing from the last.
C -1 take the next step with a new value of H, MAXORD,
C N, METH, MITER, and/or matrix parameters.
C -2 take the next step with a new value of H,
C but with other inputs unchanged.
C On return, JSTART is set to 1 to facilitate continuation.
C KFLAG = a completion code with the following meanings:
C 0 the step was succesful.
C -1 the requested error could not be achieved.
C -2 corrector convergence could not be achieved.
C -3 fatal error in PJAC or SLVS.
C A return with KFLAG = -1 or -2 means either
C abs(H) = HMIN or 10 consecutive failures occurred.
C On a return with KFLAG negative, the values of TN and
C the YH array are as of the beginning of the last
C step, and H is the last step size attempted.
C MAXORD = the maximum order of integration method to be allowed.
C MAXCOR = the maximum number of corrector iterations allowed.
C MSBP = maximum number of steps between PJAC calls (MITER .gt. 0).
C MXNCF = maximum number of convergence failures allowed.
C METH/MITER = the method flags. See description in driver.
C N = the number of first-order differential equations.
C The values of CCMAX, H, HMIN, HMXI, TN, JSTART, KFLAG, MAXORD,
C MAXCOR, MSBP, MXNCF, METH, MITER, and N are communicated via COMMON.
C
C***SEE ALSO DLSODE
C***ROUTINES CALLED DCFODE, DVNORM
C***COMMON BLOCKS DLS001
C***REVISION HISTORY (YYMMDD)
C 791129 DATE WRITTEN
C 890501 Modified prologue to SLATEC/LDOC format. (FNF)
C 890503 Minor cosmetic changes. (FNF)
C 930809 Renamed to allow single/double precision versions. (ACH)
C 010418 Reduced size of Common block /DLS001/. (ACH)
C 031105 Restored 'own' variables to Common block /DLS001/, to
C enable interrupt/restart feature. (ACH)
C***END PROLOGUE DSTODE
C**End
EXTERNAL F, JAC, PJAC, SLVS
INTEGER NEQ, NYH, IWM
DOUBLE PRECISION Y, YH, YH1, EWT, SAVF, ACOR, WM
DIMENSION NEQ(*), Y(*), YH(NYH,*), YH1(*), EWT(*), SAVF(*),
1 ACOR(*), WM(*), IWM(*)
INTEGER IOWND, IALTH, IPUP, LMAX, MEO, NQNYH, NSLP,
1 ICF, IERPJ, IERSL, JCUR, JSTART, KFLAG, L,
2 LYH, LEWT, LACOR, LSAVF, LWM, LIWM, METH, MITER,
3 MAXORD, MAXCOR, MSBP, MXNCF, N, NQ, NST, NFE, NJE, NQU
INTEGER I, I1, IREDO, IRET, J, JB, M, NCF, NEWQ
DOUBLE PRECISION CONIT, CRATE, EL, ELCO, HOLD, RMAX, TESCO,
2 CCMAX, EL0, H, HMIN, HMXI, HU, RC, TN, UROUND
DOUBLE PRECISION DCON, DDN, DEL, DELP, DSM, DUP, EXDN, EXSM, EXUP,
1 R, RH, RHDN, RHSM, RHUP, TOLD, DVNORM
COMMON /DLS001/ CONIT, CRATE, EL(13), ELCO(13,12),
1 HOLD, RMAX, TESCO(3,12),
2 CCMAX, EL0, H, HMIN, HMXI, HU, RC, TN, UROUND,
3 IOWND(6), IALTH, IPUP, LMAX, MEO, NQNYH, NSLP,
3 ICF, IERPJ, IERSL, JCUR, JSTART, KFLAG, L,
4 LYH, LEWT, LACOR, LSAVF, LWM, LIWM, METH, MITER,
5 MAXORD, MAXCOR, MSBP, MXNCF, N, NQ, NST, NFE, NJE, NQU
C
C***FIRST EXECUTABLE STATEMENT DSTODE
KFLAG = 0
TOLD = TN
NCF = 0
IERPJ = 0
IERSL = 0
JCUR = 0
ICF = 0
DELP = 0.0D0
IF (JSTART .GT. 0) GO TO 200
IF (JSTART .EQ. -1) GO TO 100
IF (JSTART .EQ. -2) GO TO 160
C-----------------------------------------------------------------------
C On the first call, the order is set to 1, and other variables are
C initialized. RMAX is the maximum ratio by which H can be increased
C in a single step. It is initially 1.E4 to compensate for the small
C initial H, but then is normally equal to 10. If a failure
C occurs (in corrector convergence or error test), RMAX is set to 2
C for the next increase.
C-----------------------------------------------------------------------
LMAX = MAXORD + 1
NQ = 1
L = 2
IALTH = 2
RMAX = 10000.0D0
RC = 0.0D0
EL0 = 1.0D0
CRATE = 0.7D0
HOLD = H
MEO = METH
NSLP = 0
IPUP = MITER
IRET = 3
GO TO 140
C-----------------------------------------------------------------------
C The following block handles preliminaries needed when JSTART = -1.
C IPUP is set to MITER to force a matrix update.
C If an order increase is about to be considered (IALTH = 1),
C IALTH is reset to 2 to postpone consideration one more step.
C If the caller has changed METH, DCFODE is called to reset
C the coefficients of the method.
C If the caller has changed MAXORD to a value less than the current
C order NQ, NQ is reduced to MAXORD, and a new H chosen accordingly.
C If H is to be changed, YH must be rescaled.
C If H or METH is being changed, IALTH is reset to L = NQ + 1
C to prevent further changes in H for that many steps.
C-----------------------------------------------------------------------
100 IPUP = MITER
LMAX = MAXORD + 1
IF (IALTH .EQ. 1) IALTH = 2
IF (METH .EQ. MEO) GO TO 110
CALL DCFODE (METH, ELCO, TESCO)
MEO = METH
IF (NQ .GT. MAXORD) GO TO 120
IALTH = L
IRET = 1
GO TO 150
110 IF (NQ .LE. MAXORD) GO TO 160
120 NQ = MAXORD
L = LMAX
DO 125 I = 1,L
125 EL(I) = ELCO(I,NQ)
NQNYH = NQ*NYH
RC = RC*EL(1)/EL0
EL0 = EL(1)
CONIT = 0.5D0/(NQ+2)
DDN = DVNORM (N, SAVF, EWT)/TESCO(1,L)
EXDN = 1.0D0/L
RHDN = 1.0D0/(1.3D0*DDN**EXDN + 0.0000013D0)
RH = MIN(RHDN,1.0D0)
IREDO = 3
IF (H .EQ. HOLD) GO TO 170
RH = MIN(RH,ABS(H/HOLD))
H = HOLD
GO TO 175
C-----------------------------------------------------------------------
C DCFODE is called to get all the integration coefficients for the
C current METH. Then the EL vector and related constants are reset
C whenever the order NQ is changed, or at the start of the problem.
C-----------------------------------------------------------------------
140 CALL DCFODE (METH, ELCO, TESCO)
150 DO 155 I = 1,L
155 EL(I) = ELCO(I,NQ)
NQNYH = NQ*NYH
RC = RC*EL(1)/EL0
EL0 = EL(1)
CONIT = 0.5D0/(NQ+2)
GO TO (160, 170, 200), IRET
C-----------------------------------------------------------------------
C If H is being changed, the H ratio RH is checked against
C RMAX, HMIN, and HMXI, and the YH array rescaled. IALTH is set to
C L = NQ + 1 to prevent a change of H for that many steps, unless
C forced by a convergence or error test failure.
C-----------------------------------------------------------------------
160 IF (H .EQ. HOLD) GO TO 200
RH = H/HOLD
H = HOLD
IREDO = 3
GO TO 175
170 RH = MAX(RH,HMIN/ABS(H))
175 RH = MIN(RH,RMAX)
RH = RH/MAX(1.0D0,ABS(H)*HMXI*RH)
R = 1.0D0
DO 180 J = 2,L
R = R*RH
DO 180 I = 1,N
180 YH(I,J) = YH(I,J)*R
H = H*RH
RC = RC*RH
IALTH = L
IF (IREDO .EQ. 0) GO TO 690
C-----------------------------------------------------------------------
C This section computes the predicted values by effectively
C multiplying the YH array by the Pascal Triangle matrix.
C RC is the ratio of new to old values of the coefficient H*EL(1).
C When RC differs from 1 by more than CCMAX, IPUP is set to MITER
C to force PJAC to be called, if a Jacobian is involved.
C In any case, PJAC is called at least every MSBP steps.
C-----------------------------------------------------------------------
200 IF (ABS(RC-1.0D0) .GT. CCMAX) IPUP = MITER
IF (NST .GE. NSLP+MSBP) IPUP = MITER
TN = TN + H
I1 = NQNYH + 1
DO 215 JB = 1,NQ
I1 = I1 - NYH
Cdir$ ivdep
DO 210 I = I1,NQNYH
210 YH1(I) = YH1(I) + YH1(I+NYH)
215 CONTINUE
C-----------------------------------------------------------------------
C Up to MAXCOR corrector iterations are taken. A convergence test is
C made on the R.M.S. norm of each correction, weighted by the error
C weight vector EWT. The sum of the corrections is accumulated in the
C vector ACOR(i). The YH array is not altered in the corrector loop.
C-----------------------------------------------------------------------
220 M = 0
DO 230 I = 1,N
230 Y(I) = YH(I,1)
CALL F (NEQ, TN, Y, SAVF)
NFE = NFE + 1
IF (IPUP .LE. 0) GO TO 250
C-----------------------------------------------------------------------
C If indicated, the matrix P = I - h*el(1)*J is reevaluated and
C preprocessed before starting the corrector iteration. IPUP is set
C to 0 as an indicator that this has been done.
C-----------------------------------------------------------------------
CALL PJAC (NEQ, Y, YH, NYH, EWT, ACOR, SAVF, WM, IWM, F, JAC)
IPUP = 0
RC = 1.0D0
NSLP = NST
CRATE = 0.7D0
IF (IERPJ .NE. 0) GO TO 430
250 DO 260 I = 1,N
260 ACOR(I) = 0.0D0
270 IF (MITER .NE. 0) GO TO 350
C-----------------------------------------------------------------------
C In the case of functional iteration, update Y directly from
C the result of the last function evaluation.
C-----------------------------------------------------------------------
DO 290 I = 1,N
SAVF(I) = H*SAVF(I) - YH(I,2)
290 Y(I) = SAVF(I) - ACOR(I)
DEL = DVNORM (N, Y, EWT)
DO 300 I = 1,N
Y(I) = YH(I,1) + EL(1)*SAVF(I)
300 ACOR(I) = SAVF(I)
GO TO 400
C-----------------------------------------------------------------------
C In the case of the chord method, compute the corrector error,
C and solve the linear system with that as right-hand side and
C P as coefficient matrix.
C-----------------------------------------------------------------------
350 DO 360 I = 1,N
360 Y(I) = H*SAVF(I) - (YH(I,2) + ACOR(I))
CALL SLVS (WM, IWM, Y, SAVF)
IF (IERSL .LT. 0) GO TO 430
IF (IERSL .GT. 0) GO TO 410
DEL = DVNORM (N, Y, EWT)
DO 380 I = 1,N
ACOR(I) = ACOR(I) + Y(I)
380 Y(I) = YH(I,1) + EL(1)*ACOR(I)
C-----------------------------------------------------------------------
C Test for convergence. If M.gt.0, an estimate of the convergence
C rate constant is stored in CRATE, and this is used in the test.
C-----------------------------------------------------------------------
400 IF (M .NE. 0) CRATE = MAX(0.2D0*CRATE,DEL/DELP)
DCON = DEL*MIN(1.0D0,1.5D0*CRATE)/(TESCO(2,NQ)*CONIT)
IF (DCON .LE. 1.0D0) GO TO 450
M = M + 1
IF (M .EQ. MAXCOR) GO TO 410
IF (M .GE. 2 .AND. DEL .GT. 2.0D0*DELP) GO TO 410
DELP = DEL
CALL F (NEQ, TN, Y, SAVF)
NFE = NFE + 1
GO TO 270
C-----------------------------------------------------------------------
C The corrector iteration failed to converge.
C If MITER .ne. 0 and the Jacobian is out of date, PJAC is called for
C the next try. Otherwise the YH array is retracted to its values
C before prediction, and H is reduced, if possible. If H cannot be
C reduced or MXNCF failures have occurred, exit with KFLAG = -2.
C-----------------------------------------------------------------------
410 IF (MITER .EQ. 0 .OR. JCUR .EQ. 1) GO TO 430
ICF = 1
IPUP = MITER
GO TO 220
430 ICF = 2
NCF = NCF + 1
RMAX = 2.0D0
TN = TOLD
I1 = NQNYH + 1
DO 445 JB = 1,NQ
I1 = I1 - NYH
Cdir$ ivdep
DO 440 I = I1,NQNYH
440 YH1(I) = YH1(I) - YH1(I+NYH)
445 CONTINUE
IF (IERPJ .LT. 0 .OR. IERSL .LT. 0) GO TO 680
IF (ABS(H) .LE. HMIN*1.00001D0) GO TO 670
IF (NCF .EQ. MXNCF) GO TO 670
RH = 0.25D0
IPUP = MITER
IREDO = 1
GO TO 170
C-----------------------------------------------------------------------
C The corrector has converged. JCUR is set to 0
C to signal that the Jacobian involved may need updating later.
C The local error test is made and control passes to statement 500
C if it fails.
C-----------------------------------------------------------------------
450 JCUR = 0
IF (M .EQ. 0) DSM = DEL/TESCO(2,NQ)
IF (M .GT. 0) DSM = DVNORM (N, ACOR, EWT)/TESCO(2,NQ)
IF (DSM .GT. 1.0D0) GO TO 500
C-----------------------------------------------------------------------
C After a successful step, update the YH array.
C Consider changing H if IALTH = 1. Otherwise decrease IALTH by 1.
C If IALTH is then 1 and NQ .lt. MAXORD, then ACOR is saved for
C use in a possible order increase on the next step.
C If a change in H is considered, an increase or decrease in order
C by one is considered also. A change in H is made only if it is by a
C factor of at least 1.1. If not, IALTH is set to 3 to prevent
C testing for that many steps.
C-----------------------------------------------------------------------
KFLAG = 0
IREDO = 0
NST = NST + 1
HU = H
NQU = NQ
DO 470 J = 1,L
DO 470 I = 1,N
470 YH(I,J) = YH(I,J) + EL(J)*ACOR(I)
IALTH = IALTH - 1
IF (IALTH .EQ. 0) GO TO 520
IF (IALTH .GT. 1) GO TO 700
IF (L .EQ. LMAX) GO TO 700
DO 490 I = 1,N
490 YH(I,LMAX) = ACOR(I)
GO TO 700
C-----------------------------------------------------------------------
C The error test failed. KFLAG keeps track of multiple failures.
C Restore TN and the YH array to their previous values, and prepare
C to try the step again. Compute the optimum step size for this or
C one lower order. After 2 or more failures, H is forced to decrease
C by a factor of 0.2 or less.
C-----------------------------------------------------------------------
500 KFLAG = KFLAG - 1
TN = TOLD
I1 = NQNYH + 1
DO 515 JB = 1,NQ
I1 = I1 - NYH
Cdir$ ivdep
DO 510 I = I1,NQNYH
510 YH1(I) = YH1(I) - YH1(I+NYH)
515 CONTINUE
RMAX = 2.0D0
IF (ABS(H) .LE. HMIN*1.00001D0) GO TO 660
IF (KFLAG .LE. -3) GO TO 640
IREDO = 2
RHUP = 0.0D0
GO TO 540
C-----------------------------------------------------------------------
C Regardless of the success or failure of the step, factors
C RHDN, RHSM, and RHUP are computed, by which H could be multiplied
C at order NQ - 1, order NQ, or order NQ + 1, respectively.
C In the case of failure, RHUP = 0.0 to avoid an order increase.
C The largest of these is determined and the new order chosen
C accordingly. If the order is to be increased, we compute one
C additional scaled derivative.
C-----------------------------------------------------------------------
520 RHUP = 0.0D0
IF (L .EQ. LMAX) GO TO 540
DO 530 I = 1,N
530 SAVF(I) = ACOR(I) - YH(I,LMAX)
DUP = DVNORM (N, SAVF, EWT)/TESCO(3,NQ)
EXUP = 1.0D0/(L+1)
RHUP = 1.0D0/(1.4D0*DUP**EXUP + 0.0000014D0)
540 EXSM = 1.0D0/L
RHSM = 1.0D0/(1.2D0*DSM**EXSM + 0.0000012D0)
RHDN = 0.0D0
IF (NQ .EQ. 1) GO TO 560
DDN = DVNORM (N, YH(1,L), EWT)/TESCO(1,NQ)
EXDN = 1.0D0/NQ
RHDN = 1.0D0/(1.3D0*DDN**EXDN + 0.0000013D0)
560 IF (RHSM .GE. RHUP) GO TO 570
IF (RHUP .GT. RHDN) GO TO 590
GO TO 580
570 IF (RHSM .LT. RHDN) GO TO 580
NEWQ = NQ
RH = RHSM
GO TO 620
580 NEWQ = NQ - 1
RH = RHDN
IF (KFLAG .LT. 0 .AND. RH .GT. 1.0D0) RH = 1.0D0
GO TO 620
590 NEWQ = L
RH = RHUP
IF (RH .LT. 1.1D0) GO TO 610
R = EL(L)/L
DO 600 I = 1,N
600 YH(I,NEWQ+1) = ACOR(I)*R
GO TO 630
610 IALTH = 3
GO TO 700
620 IF ((KFLAG .EQ. 0) .AND. (RH .LT. 1.1D0)) GO TO 610
IF (KFLAG .LE. -2) RH = MIN(RH,0.2D0)
C-----------------------------------------------------------------------
C If there is a change of order, reset NQ, l, and the coefficients.
C In any case H is reset according to RH and the YH array is rescaled.
C Then exit from 690 if the step was OK, or redo the step otherwise.
C-----------------------------------------------------------------------
IF (NEWQ .EQ. NQ) GO TO 170
630 NQ = NEWQ
L = NQ + 1
IRET = 2
GO TO 150
C-----------------------------------------------------------------------
C Control reaches this section if 3 or more failures have occured.
C If 10 failures have occurred, exit with KFLAG = -1.
C It is assumed that the derivatives that have accumulated in the
C YH array have errors of the wrong order. Hence the first
C derivative is recomputed, and the order is set to 1. Then
C H is reduced by a factor of 10, and the step is retried,
C until it succeeds or H reaches HMIN.
C-----------------------------------------------------------------------
640 IF (KFLAG .EQ. -10) GO TO 660
RH = 0.1D0
RH = MAX(HMIN/ABS(H),RH)
H = H*RH
DO 645 I = 1,N
645 Y(I) = YH(I,1)
CALL F (NEQ, TN, Y, SAVF)
NFE = NFE + 1
DO 650 I = 1,N
650 YH(I,2) = H*SAVF(I)
IPUP = MITER
IALTH = 5
IF (NQ .EQ. 1) GO TO 200
NQ = 1
L = 2
IRET = 3
GO TO 150
C-----------------------------------------------------------------------
C All returns are made through this section. H is saved in HOLD
C to allow the caller to change H on the next step.
C-----------------------------------------------------------------------
660 KFLAG = -1
GO TO 720
670 KFLAG = -2
GO TO 720
680 KFLAG = -3
GO TO 720
690 RMAX = 10.0D0
700 R = 1.0D0/TESCO(2,NQU)
DO 710 I = 1,N
710 ACOR(I) = ACOR(I)*R
720 HOLD = H
JSTART = 1
RETURN
C----------------------- END OF SUBROUTINE DSTODE ----------------------
END