GNU Octave  9.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
dlsode.f
Go to the documentation of this file.
1  SUBROUTINE dlsode (F, NEQ, Y, T, TOUT, ITOL, RTOL, ATOL, ITASK,
2  1 ISTATE, IOPT, RWORK, LRW, IWORK, LIW, JAC, MF)
3  EXTERNAL f, jac
4  INTEGER NEQ, ITOL, ITASK, ISTATE, IOPT, LRW, IWORK, LIW, MF
5  DOUBLE PRECISION Y, T, TOUT, RTOL, ATOL, RWORK
6  dimension neq(*), y(*), rtol(*), atol(*), rwork(lrw), iwork(liw)
7 C***BEGIN PROLOGUE DLSODE
8 C***PURPOSE Livermore Solver for Ordinary Differential Equations.
9 C DLSODE solves the initial-value problem for stiff or
10 C nonstiff systems of first-order ODE's,
11 C dy/dt = f(t,y), or, in component form,
12 C dy(i)/dt = f(i) = f(i,t,y(1),y(2),...,y(N)), i=1,...,N.
13 C***CATEGORY I1A
14 C***TYPE DOUBLE PRECISION (SLSODE-S, DLSODE-D)
15 C***KEYWORDS ORDINARY DIFFERENTIAL EQUATIONS, INITIAL VALUE PROBLEM,
16 C STIFF, NONSTIFF
17 C***AUTHOR Hindmarsh, Alan C., (LLNL)
18 C Center for Applied Scientific Computing, L-561
19 C Lawrence Livermore National Laboratory
20 C Livermore, CA 94551.
21 C***DESCRIPTION
22 C
23 C NOTE: The "Usage" and "Arguments" sections treat only a subset of
24 C available options, in condensed fashion. The options
25 C covered and the information supplied will support most
26 C standard uses of DLSODE.
27 C
28 C For more sophisticated uses, full details on all options are
29 C given in the concluding section, headed "Long Description."
30 C A synopsis of the DLSODE Long Description is provided at the
31 C beginning of that section; general topics covered are:
32 C - Elements of the call sequence; optional input and output
33 C - Optional supplemental routines in the DLSODE package
34 C - internal COMMON block
35 C
36 C *Usage:
37 C Communication between the user and the DLSODE package, for normal
38 C situations, is summarized here. This summary describes a subset
39 C of the available options. See "Long Description" for complete
40 C details, including optional communication, nonstandard options,
41 C and instructions for special situations.
42 C
43 C A sample program is given in the "Examples" section.
44 C
45 C Refer to the argument descriptions for the definitions of the
46 C quantities that appear in the following sample declarations.
47 C
48 C For MF = 10,
49 C PARAMETER (LRW = 20 + 16*NEQ, LIW = 20)
50 C For MF = 21 or 22,
51 C PARAMETER (LRW = 22 + 9*NEQ + NEQ**2, LIW = 20 + NEQ)
52 C For MF = 24 or 25,
53 C PARAMETER (LRW = 22 + 10*NEQ + (2*ML+MU)*NEQ,
54 C * LIW = 20 + NEQ)
55 C
56 C EXTERNAL F, JAC
57 C INTEGER NEQ, ITOL, ITASK, ISTATE, IOPT, LRW, IWORK(LIW),
58 C * LIW, MF
59 C DOUBLE PRECISION Y(NEQ), T, TOUT, RTOL, ATOL(ntol), RWORK(LRW)
60 C
61 C CALL DLSODE (F, NEQ, Y, T, TOUT, ITOL, RTOL, ATOL, ITASK,
62 C * ISTATE, IOPT, RWORK, LRW, IWORK, LIW, JAC, MF)
63 C
64 C *Arguments:
65 C F :EXT Name of subroutine for right-hand-side vector f.
66 C This name must be declared EXTERNAL in calling
67 C program. The form of F must be:
68 C
69 C SUBROUTINE F (NEQ, T, Y, YDOT)
70 C INTEGER NEQ
71 C DOUBLE PRECISION T, Y(*), YDOT(*)
72 C
73 C The inputs are NEQ, T, Y. F is to set
74 C
75 C YDOT(i) = f(i,T,Y(1),Y(2),...,Y(NEQ)),
76 C i = 1, ..., NEQ .
77 C
78 C NEQ :IN Number of first-order ODE's.
79 C
80 C Y :INOUT Array of values of the y(t) vector, of length NEQ.
81 C Input: For the first call, Y should contain the
82 C values of y(t) at t = T. (Y is an input
83 C variable only if ISTATE = 1.)
84 C Output: On return, Y will contain the values at the
85 C new t-value.
86 C
87 C T :INOUT Value of the independent variable. On return it
88 C will be the current value of t (normally TOUT).
89 C
90 C TOUT :IN Next point where output is desired (.NE. T).
91 C
92 C ITOL :IN 1 or 2 according as ATOL (below) is a scalar or
93 C an array.
94 C
95 C RTOL :IN Relative tolerance parameter (scalar).
96 C
97 C ATOL :IN Absolute tolerance parameter (scalar or array).
98 C If ITOL = 1, ATOL need not be dimensioned.
99 C If ITOL = 2, ATOL must be dimensioned at least NEQ.
100 C
101 C The estimated local error in Y(i) will be controlled
102 C so as to be roughly less (in magnitude) than
103 C
104 C EWT(i) = RTOL*ABS(Y(i)) + ATOL if ITOL = 1, or
105 C EWT(i) = RTOL*ABS(Y(i)) + ATOL(i) if ITOL = 2.
106 C
107 C Thus the local error test passes if, in each
108 C component, either the absolute error is less than
109 C ATOL (or ATOL(i)), or the relative error is less
110 C than RTOL.
111 C
112 C Use RTOL = 0.0 for pure absolute error control, and
113 C use ATOL = 0.0 (or ATOL(i) = 0.0) for pure relative
114 C error control. Caution: Actual (global) errors may
115 C exceed these local tolerances, so choose them
116 C conservatively.
117 C
118 C ITASK :IN Flag indicating the task DLSODE is to perform.
119 C Use ITASK = 1 for normal computation of output
120 C values of y at t = TOUT.
121 C
122 C ISTATE:INOUT Index used for input and output to specify the state
123 C of the calculation.
124 C Input:
125 C 1 This is the first call for a problem.
126 C 2 This is a subsequent call.
127 C Output:
128 C 1 Nothing was done, because TOUT was equal to T.
129 C 2 DLSODE was successful (otherwise, negative).
130 C Note that ISTATE need not be modified after a
131 C successful return.
132 C -1 Excess work done on this call (perhaps wrong
133 C MF).
134 C -2 Excess accuracy requested (tolerances too
135 C small).
136 C -3 Illegal input detected (see printed message).
137 C -4 Repeated error test failures (check all
138 C inputs).
139 C -5 Repeated convergence failures (perhaps bad
140 C Jacobian supplied or wrong choice of MF or
141 C tolerances).
142 C -6 Error weight became zero during problem
143 C (solution component i vanished, and ATOL or
144 C ATOL(i) = 0.).
145 C
146 C IOPT :IN Flag indicating whether optional inputs are used:
147 C 0 No.
148 C 1 Yes. (See "Optional inputs" under "Long
149 C Description," Part 1.)
150 C
151 C RWORK :WORK Real work array of length at least:
152 C 20 + 16*NEQ for MF = 10,
153 C 22 + 9*NEQ + NEQ**2 for MF = 21 or 22,
154 C 22 + 10*NEQ + (2*ML + MU)*NEQ for MF = 24 or 25.
155 C
156 C LRW :IN Declared length of RWORK (in user's DIMENSION
157 C statement).
158 C
159 C IWORK :WORK Integer work array of length at least:
160 C 20 for MF = 10,
161 C 20 + NEQ for MF = 21, 22, 24, or 25.
162 C
163 C If MF = 24 or 25, input in IWORK(1),IWORK(2) the
164 C lower and upper Jacobian half-bandwidths ML,MU.
165 C
166 C On return, IWORK contains information that may be
167 C of interest to the user:
168 C
169 C Name Location Meaning
170 C ----- --------- -----------------------------------------
171 C NST IWORK(11) Number of steps taken for the problem so
172 C far.
173 C NFE IWORK(12) Number of f evaluations for the problem
174 C so far.
175 C NJE IWORK(13) Number of Jacobian evaluations (and of
176 C matrix LU decompositions) for the problem
177 C so far.
178 C NQU IWORK(14) Method order last used (successfully).
179 C LENRW IWORK(17) Length of RWORK actually required. This
180 C is defined on normal returns and on an
181 C illegal input return for insufficient
182 C storage.
183 C LENIW IWORK(18) Length of IWORK actually required. This
184 C is defined on normal returns and on an
185 C illegal input return for insufficient
186 C storage.
187 C
188 C LIW :IN Declared length of IWORK (in user's DIMENSION
189 C statement).
190 C
191 C JAC :EXT Name of subroutine for Jacobian matrix (MF =
192 C 21 or 24). If used, this name must be declared
193 C EXTERNAL in calling program. If not used, pass a
194 C dummy name. The form of JAC must be:
195 C
196 C SUBROUTINE JAC (NEQ, T, Y, ML, MU, PD, NROWPD)
197 C INTEGER NEQ, ML, MU, NROWPD
198 C DOUBLE PRECISION T, Y(*), PD(NROWPD,*)
199 C
200 C See item c, under "Description" below for more
201 C information about JAC.
202 C
203 C MF :IN Method flag. Standard values are:
204 C 10 Nonstiff (Adams) method, no Jacobian used.
205 C 21 Stiff (BDF) method, user-supplied full Jacobian.
206 C 22 Stiff method, internally generated full
207 C Jacobian.
208 C 24 Stiff method, user-supplied banded Jacobian.
209 C 25 Stiff method, internally generated banded
210 C Jacobian.
211 C
212 C *Description:
213 C DLSODE solves the initial value problem for stiff or nonstiff
214 C systems of first-order ODE's,
215 C
216 C dy/dt = f(t,y) ,
217 C
218 C or, in component form,
219 C
220 C dy(i)/dt = f(i) = f(i,t,y(1),y(2),...,y(NEQ))
221 C (i = 1, ..., NEQ) .
222 C
223 C DLSODE is a package based on the GEAR and GEARB packages, and on
224 C the October 23, 1978, version of the tentative ODEPACK user
225 C interface standard, with minor modifications.
226 C
227 C The steps in solving such a problem are as follows.
228 C
229 C a. First write a subroutine of the form
230 C
231 C SUBROUTINE F (NEQ, T, Y, YDOT)
232 C INTEGER NEQ
233 C DOUBLE PRECISION T, Y(*), YDOT(*)
234 C
235 C which supplies the vector function f by loading YDOT(i) with
236 C f(i).
237 C
238 C b. Next determine (or guess) whether or not the problem is stiff.
239 C Stiffness occurs when the Jacobian matrix df/dy has an
240 C eigenvalue whose real part is negative and large in magnitude
241 C compared to the reciprocal of the t span of interest. If the
242 C problem is nonstiff, use method flag MF = 10. If it is stiff,
243 C there are four standard choices for MF, and DLSODE requires the
244 C Jacobian matrix in some form. This matrix is regarded either
245 C as full (MF = 21 or 22), or banded (MF = 24 or 25). In the
246 C banded case, DLSODE requires two half-bandwidth parameters ML
247 C and MU. These are, respectively, the widths of the lower and
248 C upper parts of the band, excluding the main diagonal. Thus the
249 C band consists of the locations (i,j) with
250 C
251 C i - ML <= j <= i + MU ,
252 C
253 C and the full bandwidth is ML + MU + 1 .
254 C
255 C c. If the problem is stiff, you are encouraged to supply the
256 C Jacobian directly (MF = 21 or 24), but if this is not feasible,
257 C DLSODE will compute it internally by difference quotients (MF =
258 C 22 or 25). If you are supplying the Jacobian, write a
259 C subroutine of the form
260 C
261 C SUBROUTINE JAC (NEQ, T, Y, ML, MU, PD, NROWPD)
262 C INTEGER NEQ, ML, MU, NRWOPD
263 C DOUBLE PRECISION T, Y(*), PD(NROWPD,*)
264 C
265 C which provides df/dy by loading PD as follows:
266 C - For a full Jacobian (MF = 21), load PD(i,j) with df(i)/dy(j),
267 C the partial derivative of f(i) with respect to y(j). (Ignore
268 C the ML and MU arguments in this case.)
269 C - For a banded Jacobian (MF = 24), load PD(i-j+MU+1,j) with
270 C df(i)/dy(j); i.e., load the diagonal lines of df/dy into the
271 C rows of PD from the top down.
272 C - In either case, only nonzero elements need be loaded.
273 C
274 C d. Write a main program that calls subroutine DLSODE once for each
275 C point at which answers are desired. This should also provide
276 C for possible use of logical unit 6 for output of error messages
277 C by DLSODE.
278 C
279 C Before the first call to DLSODE, set ISTATE = 1, set Y and T to
280 C the initial values, and set TOUT to the first output point. To
281 C continue the integration after a successful return, simply
282 C reset TOUT and call DLSODE again. No other parameters need be
283 C reset.
284 C
285 C *Examples:
286 C The following is a simple example problem, with the coding needed
287 C for its solution by DLSODE. The problem is from chemical kinetics,
288 C and consists of the following three rate equations:
289 C
290 C dy1/dt = -.04*y1 + 1.E4*y2*y3
291 C dy2/dt = .04*y1 - 1.E4*y2*y3 - 3.E7*y2**2
292 C dy3/dt = 3.E7*y2**2
293 C
294 C on the interval from t = 0.0 to t = 4.E10, with initial conditions
295 C y1 = 1.0, y2 = y3 = 0. The problem is stiff.
296 C
297 C The following coding solves this problem with DLSODE, using
298 C MF = 21 and printing results at t = .4, 4., ..., 4.E10. It uses
299 C ITOL = 2 and ATOL much smaller for y2 than for y1 or y3 because y2
300 C has much smaller values. At the end of the run, statistical
301 C quantities of interest are printed.
302 C
303 C EXTERNAL FEX, JEX
304 C INTEGER IOPT, IOUT, ISTATE, ITASK, ITOL, IWORK(23), LIW, LRW,
305 C * MF, NEQ
306 C DOUBLE PRECISION ATOL(3), RTOL, RWORK(58), T, TOUT, Y(3)
307 C NEQ = 3
308 C Y(1) = 1.D0
309 C Y(2) = 0.D0
310 C Y(3) = 0.D0
311 C T = 0.D0
312 C TOUT = .4D0
313 C ITOL = 2
314 C RTOL = 1.D-4
315 C ATOL(1) = 1.D-6
316 C ATOL(2) = 1.D-10
317 C ATOL(3) = 1.D-6
318 C ITASK = 1
319 C ISTATE = 1
320 C IOPT = 0
321 C LRW = 58
322 C LIW = 23
323 C MF = 21
324 C DO 40 IOUT = 1,12
325 C CALL DLSODE (FEX, NEQ, Y, T, TOUT, ITOL, RTOL, ATOL, ITASK,
326 C * ISTATE, IOPT, RWORK, LRW, IWORK, LIW, JEX, MF)
327 C WRITE(6,20) T, Y(1), Y(2), Y(3)
328 C 20 FORMAT(' At t =',D12.4,' y =',3D14.6)
329 C IF (ISTATE .LT. 0) GO TO 80
330 C 40 TOUT = TOUT*10.D0
331 C WRITE(6,60) IWORK(11), IWORK(12), IWORK(13)
332 C 60 FORMAT(/' No. steps =',i4,', No. f-s =',i4,', No. J-s =',i4)
333 C STOP
334 C 80 WRITE(6,90) ISTATE
335 C 90 FORMAT(///' Error halt.. ISTATE =',I3)
336 C STOP
337 C END
338 C
339 C SUBROUTINE FEX (NEQ, T, Y, YDOT)
340 C INTEGER NEQ
341 C DOUBLE PRECISION T, Y(3), YDOT(3)
342 C YDOT(1) = -.04D0*Y(1) + 1.D4*Y(2)*Y(3)
343 C YDOT(3) = 3.D7*Y(2)*Y(2)
344 C YDOT(2) = -YDOT(1) - YDOT(3)
345 C RETURN
346 C END
347 C
348 C SUBROUTINE JEX (NEQ, T, Y, ML, MU, PD, NRPD)
349 C INTEGER NEQ, ML, MU, NRPD
350 C DOUBLE PRECISION T, Y(3), PD(NRPD,3)
351 C PD(1,1) = -.04D0
352 C PD(1,2) = 1.D4*Y(3)
353 C PD(1,3) = 1.D4*Y(2)
354 C PD(2,1) = .04D0
355 C PD(2,3) = -PD(1,3)
356 C PD(3,2) = 6.D7*Y(2)
357 C PD(2,2) = -PD(1,2) - PD(3,2)
358 C RETURN
359 C END
360 C
361 C The output from this program (on a Cray-1 in single precision)
362 C is as follows.
363 C
364 C At t = 4.0000e-01 y = 9.851726e-01 3.386406e-05 1.479357e-02
365 C At t = 4.0000e+00 y = 9.055142e-01 2.240418e-05 9.446344e-02
366 C At t = 4.0000e+01 y = 7.158050e-01 9.184616e-06 2.841858e-01
367 C At t = 4.0000e+02 y = 4.504846e-01 3.222434e-06 5.495122e-01
368 C At t = 4.0000e+03 y = 1.831701e-01 8.940379e-07 8.168290e-01
369 C At t = 4.0000e+04 y = 3.897016e-02 1.621193e-07 9.610297e-01
370 C At t = 4.0000e+05 y = 4.935213e-03 1.983756e-08 9.950648e-01
371 C At t = 4.0000e+06 y = 5.159269e-04 2.064759e-09 9.994841e-01
372 C At t = 4.0000e+07 y = 5.306413e-05 2.122677e-10 9.999469e-01
373 C At t = 4.0000e+08 y = 5.494530e-06 2.197825e-11 9.999945e-01
374 C At t = 4.0000e+09 y = 5.129458e-07 2.051784e-12 9.999995e-01
375 C At t = 4.0000e+10 y = -7.170603e-08 -2.868241e-13 1.000000e+00
376 C
377 C No. steps = 330, No. f-s = 405, No. J-s = 69
378 C
379 C *Accuracy:
380 C The accuracy of the solution depends on the choice of tolerances
381 C RTOL and ATOL. Actual (global) errors may exceed these local
382 C tolerances, so choose them conservatively.
383 C
384 C *Cautions:
385 C The work arrays should not be altered between calls to DLSODE for
386 C the same problem, except possibly for the conditional and optional
387 C inputs.
388 C
389 C *Portability:
390 C Since NEQ is dimensioned inside DLSODE, some compilers may object
391 C to a call to DLSODE with NEQ a scalar variable. In this event,
392 C use DIMENSION NEQ(1). Similar remarks apply to RTOL and ATOL.
393 C
394 C Note to Cray users:
395 C For maximum efficiency, use the CFT77 compiler. Appropriate
396 C compiler optimization directives have been inserted for CFT77.
397 C
398 C *Reference:
399 C Alan C. Hindmarsh, "ODEPACK, A Systematized Collection of ODE
400 C Solvers," in Scientific Computing, R. S. Stepleman, et al., Eds.
401 C (North-Holland, Amsterdam, 1983), pp. 55-64.
402 C
403 C *Long Description:
404 C The following complete description of the user interface to
405 C DLSODE consists of four parts:
406 C
407 C 1. The call sequence to subroutine DLSODE, which is a driver
408 C routine for the solver. This includes descriptions of both
409 C the call sequence arguments and user-supplied routines.
410 C Following these descriptions is a description of optional
411 C inputs available through the call sequence, and then a
412 C description of optional outputs in the work arrays.
413 C
414 C 2. Descriptions of other routines in the DLSODE package that may
415 C be (optionally) called by the user. These provide the ability
416 C to alter error message handling, save and restore the internal
417 C COMMON, and obtain specified derivatives of the solution y(t).
418 C
419 C 3. Descriptions of COMMON block to be declared in overlay or
420 C similar environments, or to be saved when doing an interrupt
421 C of the problem and continued solution later.
422 C
423 C 4. Description of two routines in the DLSODE package, either of
424 C which the user may replace with his own version, if desired.
425 C These relate to the measurement of errors.
426 C
427 C
428 C Part 1. Call Sequence
429 C ----------------------
430 C
431 C Arguments
432 C ---------
433 C The call sequence parameters used for input only are
434 C
435 C F, NEQ, TOUT, ITOL, RTOL, ATOL, ITASK, IOPT, LRW, LIW, JAC, MF,
436 C
437 C and those used for both input and output are
438 C
439 C Y, T, ISTATE.
440 C
441 C The work arrays RWORK and IWORK are also used for conditional and
442 C optional inputs and optional outputs. (The term output here
443 C refers to the return from subroutine DLSODE to the user's calling
444 C program.)
445 C
446 C The legality of input parameters will be thoroughly checked on the
447 C initial call for the problem, but not checked thereafter unless a
448 C change in input parameters is flagged by ISTATE = 3 on input.
449 C
450 C The descriptions of the call arguments are as follows.
451 C
452 C F The name of the user-supplied subroutine defining the ODE
453 C system. The system must be put in the first-order form
454 C dy/dt = f(t,y), where f is a vector-valued function of
455 C the scalar t and the vector y. Subroutine F is to compute
456 C the function f. It is to have the form
457 C
458 C SUBROUTINE F (NEQ, T, Y, YDOT)
459 C DOUBLE PRECISION T, Y(*), YDOT(*)
460 C
461 C where NEQ, T, and Y are input, and the array YDOT =
462 C f(T,Y) is output. Y and YDOT are arrays of length NEQ.
463 C Subroutine F should not alter Y(1),...,Y(NEQ). F must be
464 C declared EXTERNAL in the calling program.
465 C
466 C Subroutine F may access user-defined quantities in
467 C NEQ(2),... and/or in Y(NEQ(1)+1),..., if NEQ is an array
468 C (dimensioned in F) and/or Y has length exceeding NEQ(1).
469 C See the descriptions of NEQ and Y below.
470 C
471 C If quantities computed in the F routine are needed
472 C externally to DLSODE, an extra call to F should be made
473 C for this purpose, for consistent and accurate results.
474 C If only the derivative dy/dt is needed, use DINTDY
475 C instead.
476 C
477 C NEQ The size of the ODE system (number of first-order
478 C ordinary differential equations). Used only for input.
479 C NEQ may be decreased, but not increased, during the
480 C problem. If NEQ is decreased (with ISTATE = 3 on input),
481 C the remaining components of Y should be left undisturbed,
482 C if these are to be accessed in F and/or JAC.
483 C
484 C Normally, NEQ is a scalar, and it is generally referred
485 C to as a scalar in this user interface description.
486 C However, NEQ may be an array, with NEQ(1) set to the
487 C system size. (The DLSODE package accesses only NEQ(1).)
488 C In either case, this parameter is passed as the NEQ
489 C argument in all calls to F and JAC. Hence, if it is an
490 C array, locations NEQ(2),... may be used to store other
491 C integer data and pass it to F and/or JAC. Subroutines
492 C F and/or JAC must include NEQ in a DIMENSION statement
493 C in that case.
494 C
495 C Y A real array for the vector of dependent variables, of
496 C length NEQ or more. Used for both input and output on
497 C the first call (ISTATE = 1), and only for output on
498 C other calls. On the first call, Y must contain the
499 C vector of initial values. On output, Y contains the
500 C computed solution vector, evaluated at T. If desired,
501 C the Y array may be used for other purposes between
502 C calls to the solver.
503 C
504 C This array is passed as the Y argument in all calls to F
505 C and JAC. Hence its length may exceed NEQ, and locations
506 C Y(NEQ+1),... may be used to store other real data and
507 C pass it to F and/or JAC. (The DLSODE package accesses
508 C only Y(1),...,Y(NEQ).)
509 C
510 C T The independent variable. On input, T is used only on
511 C the first call, as the initial point of the integration.
512 C On output, after each call, T is the value at which a
513 C computed solution Y is evaluated (usually the same as
514 C TOUT). On an error return, T is the farthest point
515 C reached.
516 C
517 C TOUT The next value of T at which a computed solution is
518 C desired. Used only for input.
519 C
520 C When starting the problem (ISTATE = 1), TOUT may be equal
521 C to T for one call, then should not equal T for the next
522 C call. For the initial T, an input value of TOUT .NE. T
523 C is used in order to determine the direction of the
524 C integration (i.e., the algebraic sign of the step sizes)
525 C and the rough scale of the problem. Integration in
526 C either direction (forward or backward in T) is permitted.
527 C
528 C If ITASK = 2 or 5 (one-step modes), TOUT is ignored
529 C after the first call (i.e., the first call with
530 C TOUT .NE. T). Otherwise, TOUT is required on every call.
531 C
532 C If ITASK = 1, 3, or 4, the values of TOUT need not be
533 C monotone, but a value of TOUT which backs up is limited
534 C to the current internal T interval, whose endpoints are
535 C TCUR - HU and TCUR. (See "Optional Outputs" below for
536 C TCUR and HU.)
537 C
538 C
539 C ITOL An indicator for the type of error control. See
540 C description below under ATOL. Used only for input.
541 C
542 C RTOL A relative error tolerance parameter, either a scalar or
543 C an array of length NEQ. See description below under
544 C ATOL. Input only.
545 C
546 C ATOL An absolute error tolerance parameter, either a scalar or
547 C an array of length NEQ. Input only.
548 C
549 C The input parameters ITOL, RTOL, and ATOL determine the
550 C error control performed by the solver. The solver will
551 C control the vector e = (e(i)) of estimated local errors
552 C in Y, according to an inequality of the form
553 C
554 C rms-norm of ( e(i)/EWT(i) ) <= 1,
555 C
556 C where
557 C
558 C EWT(i) = RTOL(i)*ABS(Y(i)) + ATOL(i),
559 C
560 C and the rms-norm (root-mean-square norm) here is
561 C
562 C rms-norm(v) = SQRT(sum v(i)**2 / NEQ).
563 C
564 C Here EWT = (EWT(i)) is a vector of weights which must
565 C always be positive, and the values of RTOL and ATOL
566 C should all be nonnegative. The following table gives the
567 C types (scalar/array) of RTOL and ATOL, and the
568 C corresponding form of EWT(i).
569 C
570 C ITOL RTOL ATOL EWT(i)
571 C ---- ------ ------ -----------------------------
572 C 1 scalar scalar RTOL*ABS(Y(i)) + ATOL
573 C 2 scalar array RTOL*ABS(Y(i)) + ATOL(i)
574 C 3 array scalar RTOL(i)*ABS(Y(i)) + ATOL
575 C 4 array array RTOL(i)*ABS(Y(i)) + ATOL(i)
576 C
577 C When either of these parameters is a scalar, it need not
578 C be dimensioned in the user's calling program.
579 C
580 C If none of the above choices (with ITOL, RTOL, and ATOL
581 C fixed throughout the problem) is suitable, more general
582 C error controls can be obtained by substituting
583 C user-supplied routines for the setting of EWT and/or for
584 C the norm calculation. See Part 4 below.
585 C
586 C If global errors are to be estimated by making a repeated
587 C run on the same problem with smaller tolerances, then all
588 C components of RTOL and ATOL (i.e., of EWT) should be
589 C scaled down uniformly.
590 C
591 C ITASK An index specifying the task to be performed. Input
592 C only. ITASK has the following values and meanings:
593 C 1 Normal computation of output values of y(t) at
594 C t = TOUT (by overshooting and interpolating).
595 C 2 Take one step only and return.
596 C 3 Stop at the first internal mesh point at or beyond
597 C t = TOUT and return.
598 C 4 Normal computation of output values of y(t) at
599 C t = TOUT but without overshooting t = TCRIT. TCRIT
600 C must be input as RWORK(1). TCRIT may be equal to or
601 C beyond TOUT, but not behind it in the direction of
602 C integration. This option is useful if the problem
603 C has a singularity at or beyond t = TCRIT.
604 C 5 Take one step, without passing TCRIT, and return.
605 C TCRIT must be input as RWORK(1).
606 C
607 C Note: If ITASK = 4 or 5 and the solver reaches TCRIT
608 C (within roundoff), it will return T = TCRIT (exactly) to
609 C indicate this (unless ITASK = 4 and TOUT comes before
610 C TCRIT, in which case answers at T = TOUT are returned
611 C first).
612 C
613 C ISTATE An index used for input and output to specify the state
614 C of the calculation.
615 C
616 C On input, the values of ISTATE are as follows:
617 C 1 This is the first call for the problem
618 C (initializations will be done). See "Note" below.
619 C 2 This is not the first call, and the calculation is to
620 C continue normally, with no change in any input
621 C parameters except possibly TOUT and ITASK. (If ITOL,
622 C RTOL, and/or ATOL are changed between calls with
623 C ISTATE = 2, the new values will be used but not
624 C tested for legality.)
625 C 3 This is not the first call, and the calculation is to
626 C continue normally, but with a change in input
627 C parameters other than TOUT and ITASK. Changes are
628 C allowed in NEQ, ITOL, RTOL, ATOL, IOPT, LRW, LIW, MF,
629 C ML, MU, and any of the optional inputs except H0.
630 C (See IWORK description for ML and MU.)
631 C
632 C Note: A preliminary call with TOUT = T is not counted as
633 C a first call here, as no initialization or checking of
634 C input is done. (Such a call is sometimes useful for the
635 C purpose of outputting the initial conditions.) Thus the
636 C first call for which TOUT .NE. T requires ISTATE = 1 on
637 C input.
638 C
639 C On output, ISTATE has the following values and meanings:
640 C 1 Nothing was done, as TOUT was equal to T with
641 C ISTATE = 1 on input.
642 C 2 The integration was performed successfully.
643 C -1 An excessive amount of work (more than MXSTEP steps)
644 C was done on this call, before completing the
645 C requested task, but the integration was otherwise
646 C successful as far as T. (MXSTEP is an optional input
647 C and is normally 500.) To continue, the user may
648 C simply reset ISTATE to a value >1 and call again (the
649 C excess work step counter will be reset to 0). In
650 C addition, the user may increase MXSTEP to avoid this
651 C error return; see "Optional Inputs" below.
652 C -2 Too much accuracy was requested for the precision of
653 C the machine being used. This was detected before
654 C completing the requested task, but the integration
655 C was successful as far as T. To continue, the
656 C tolerance parameters must be reset, and ISTATE must
657 C be set to 3. The optional output TOLSF may be used
658 C for this purpose. (Note: If this condition is
659 C detected before taking any steps, then an illegal
660 C input return (ISTATE = -3) occurs instead.)
661 C -3 Illegal input was detected, before taking any
662 C integration steps. See written message for details.
663 C (Note: If the solver detects an infinite loop of
664 C calls to the solver with illegal input, it will cause
665 C the run to stop.)
666 C -4 There were repeated error-test failures on one
667 C attempted step, before completing the requested task,
668 C but the integration was successful as far as T. The
669 C problem may have a singularity, or the input may be
670 C inappropriate.
671 C -5 There were repeated convergence-test failures on one
672 C attempted step, before completing the requested task,
673 C but the integration was successful as far as T. This
674 C may be caused by an inaccurate Jacobian matrix, if
675 C one is being used.
676 C -6 EWT(i) became zero for some i during the integration.
677 C Pure relative error control (ATOL(i)=0.0) was
678 C requested on a variable which has now vanished. The
679 C integration was successful as far as T.
680 C
681 C Note: Since the normal output value of ISTATE is 2, it
682 C does not need to be reset for normal continuation. Also,
683 C since a negative input value of ISTATE will be regarded
684 C as illegal, a negative output value requires the user to
685 C change it, and possibly other inputs, before calling the
686 C solver again.
687 C
688 C IOPT An integer flag to specify whether any optional inputs
689 C are being used on this call. Input only. The optional
690 C inputs are listed under a separate heading below.
691 C 0 No optional inputs are being used. Default values
692 C will be used in all cases.
693 C 1 One or more optional inputs are being used.
694 C
695 C RWORK A real working array (double precision). The length of
696 C RWORK must be at least
697 C
698 C 20 + NYH*(MAXORD + 1) + 3*NEQ + LWM
699 C
700 C where
701 C NYH = the initial value of NEQ,
702 C MAXORD = 12 (if METH = 1) or 5 (if METH = 2) (unless a
703 C smaller value is given as an optional input),
704 C LWM = 0 if MITER = 0,
705 C LWM = NEQ**2 + 2 if MITER = 1 or 2,
706 C LWM = NEQ + 2 if MITER = 3, and
707 C LWM = (2*ML + MU + 1)*NEQ + 2
708 C if MITER = 4 or 5.
709 C (See the MF description below for METH and MITER.)
710 C
711 C Thus if MAXORD has its default value and NEQ is constant,
712 C this length is:
713 C 20 + 16*NEQ for MF = 10,
714 C 22 + 16*NEQ + NEQ**2 for MF = 11 or 12,
715 C 22 + 17*NEQ for MF = 13,
716 C 22 + 17*NEQ + (2*ML + MU)*NEQ for MF = 14 or 15,
717 C 20 + 9*NEQ for MF = 20,
718 C 22 + 9*NEQ + NEQ**2 for MF = 21 or 22,
719 C 22 + 10*NEQ for MF = 23,
720 C 22 + 10*NEQ + (2*ML + MU)*NEQ for MF = 24 or 25.
721 C
722 C The first 20 words of RWORK are reserved for conditional
723 C and optional inputs and optional outputs.
724 C
725 C The following word in RWORK is a conditional input:
726 C RWORK(1) = TCRIT, the critical value of t which the
727 C solver is not to overshoot. Required if ITASK
728 C is 4 or 5, and ignored otherwise. See ITASK.
729 C
730 C LRW The length of the array RWORK, as declared by the user.
731 C (This will be checked by the solver.)
732 C
733 C IWORK An integer work array. Its length must be at least
734 C 20 if MITER = 0 or 3 (MF = 10, 13, 20, 23), or
735 C 20 + NEQ otherwise (MF = 11, 12, 14, 15, 21, 22, 24, 25).
736 C (See the MF description below for MITER.) The first few
737 C words of IWORK are used for conditional and optional
738 C inputs and optional outputs.
739 C
740 C The following two words in IWORK are conditional inputs:
741 C IWORK(1) = ML These are the lower and upper half-
742 C IWORK(2) = MU bandwidths, respectively, of the banded
743 C Jacobian, excluding the main diagonal.
744 C The band is defined by the matrix locations
745 C (i,j) with i - ML <= j <= i + MU. ML and MU
746 C must satisfy 0 <= ML,MU <= NEQ - 1. These are
747 C required if MITER is 4 or 5, and ignored
748 C otherwise. ML and MU may in fact be the band
749 C parameters for a matrix to which df/dy is only
750 C approximately equal.
751 C
752 C LIW The length of the array IWORK, as declared by the user.
753 C (This will be checked by the solver.)
754 C
755 C Note: The work arrays must not be altered between calls to DLSODE
756 C for the same problem, except possibly for the conditional and
757 C optional inputs, and except for the last 3*NEQ words of RWORK.
758 C The latter space is used for internal scratch space, and so is
759 C available for use by the user outside DLSODE between calls, if
760 C desired (but not for use by F or JAC).
761 C
762 C JAC The name of the user-supplied routine (MITER = 1 or 4) to
763 C compute the Jacobian matrix, df/dy, as a function of the
764 C scalar t and the vector y. (See the MF description below
765 C for MITER.) It is to have the form
766 C
767 C SUBROUTINE JAC (NEQ, T, Y, ML, MU, PD, NROWPD)
768 C DOUBLE PRECISION T, Y(*), PD(NROWPD,*)
769 C
770 C where NEQ, T, Y, ML, MU, and NROWPD are input and the
771 C array PD is to be loaded with partial derivatives
772 C (elements of the Jacobian matrix) on output. PD must be
773 C given a first dimension of NROWPD. T and Y have the same
774 C meaning as in subroutine F.
775 C
776 C In the full matrix case (MITER = 1), ML and MU are
777 C ignored, and the Jacobian is to be loaded into PD in
778 C columnwise manner, with df(i)/dy(j) loaded into PD(i,j).
779 C
780 C In the band matrix case (MITER = 4), the elements within
781 C the band are to be loaded into PD in columnwise manner,
782 C with diagonal lines of df/dy loaded into the rows of PD.
783 C Thus df(i)/dy(j) is to be loaded into PD(i-j+MU+1,j). ML
784 C and MU are the half-bandwidth parameters (see IWORK).
785 C The locations in PD in the two triangular areas which
786 C correspond to nonexistent matrix elements can be ignored
787 C or loaded arbitrarily, as they are overwritten by DLSODE.
788 C
789 C JAC need not provide df/dy exactly. A crude approximation
790 C (possibly with a smaller bandwidth) will do.
791 C
792 C In either case, PD is preset to zero by the solver, so
793 C that only the nonzero elements need be loaded by JAC.
794 C Each call to JAC is preceded by a call to F with the same
795 C arguments NEQ, T, and Y. Thus to gain some efficiency,
796 C intermediate quantities shared by both calculations may
797 C be saved in a user COMMON block by F and not recomputed
798 C by JAC, if desired. Also, JAC may alter the Y array, if
799 C desired. JAC must be declared EXTERNAL in the calling
800 C program.
801 C
802 C Subroutine JAC may access user-defined quantities in
803 C NEQ(2),... and/or in Y(NEQ(1)+1),... if NEQ is an array
804 C (dimensioned in JAC) and/or Y has length exceeding
805 C NEQ(1). See the descriptions of NEQ and Y above.
806 C
807 C MF The method flag. Used only for input. The legal values
808 C of MF are 10, 11, 12, 13, 14, 15, 20, 21, 22, 23, 24,
809 C and 25. MF has decimal digits METH and MITER:
810 C MF = 10*METH + MITER .
811 C
812 C METH indicates the basic linear multistep method:
813 C 1 Implicit Adams method.
814 C 2 Method based on backward differentiation formulas
815 C (BDF's).
816 C
817 C MITER indicates the corrector iteration method:
818 C 0 Functional iteration (no Jacobian matrix is
819 C involved).
820 C 1 Chord iteration with a user-supplied full (NEQ by
821 C NEQ) Jacobian.
822 C 2 Chord iteration with an internally generated
823 C (difference quotient) full Jacobian (using NEQ
824 C extra calls to F per df/dy value).
825 C 3 Chord iteration with an internally generated
826 C diagonal Jacobian approximation (using one extra call
827 C to F per df/dy evaluation).
828 C 4 Chord iteration with a user-supplied banded Jacobian.
829 C 5 Chord iteration with an internally generated banded
830 C Jacobian (using ML + MU + 1 extra calls to F per
831 C df/dy evaluation).
832 C
833 C If MITER = 1 or 4, the user must supply a subroutine JAC
834 C (the name is arbitrary) as described above under JAC.
835 C For other values of MITER, a dummy argument can be used.
836 C
837 C Optional Inputs
838 C ---------------
839 C The following is a list of the optional inputs provided for in the
840 C call sequence. (See also Part 2.) For each such input variable,
841 C this table lists its name as used in this documentation, its
842 C location in the call sequence, its meaning, and the default value.
843 C The use of any of these inputs requires IOPT = 1, and in that case
844 C all of these inputs are examined. A value of zero for any of
845 C these optional inputs will cause the default value to be used.
846 C Thus to use a subset of the optional inputs, simply preload
847 C locations 5 to 10 in RWORK and IWORK to 0.0 and 0 respectively,
848 C and then set those of interest to nonzero values.
849 C
850 C Name Location Meaning and default value
851 C ------ --------- -----------------------------------------------
852 C H0 RWORK(5) Step size to be attempted on the first step.
853 C The default value is determined by the solver.
854 C HMAX RWORK(6) Maximum absolute step size allowed. The
855 C default value is infinite.
856 C HMIN RWORK(7) Minimum absolute step size allowed. The
857 C default value is 0. (This lower bound is not
858 C enforced on the final step before reaching
859 C TCRIT when ITASK = 4 or 5.)
860 C MAXORD IWORK(5) Maximum order to be allowed. The default value
861 C is 12 if METH = 1, and 5 if METH = 2. (See the
862 C MF description above for METH.) If MAXORD
863 C exceeds the default value, it will be reduced
864 C to the default value. If MAXORD is changed
865 C during the problem, it may cause the current
866 C order to be reduced.
867 C MXSTEP IWORK(6) Maximum number of (internally defined) steps
868 C allowed during one call to the solver. The
869 C default value is 500.
870 C MXHNIL IWORK(7) Maximum number of messages printed (per
871 C problem) warning that T + H = T on a step
872 C (H = step size). This must be positive to
873 C result in a nondefault value. The default
874 C value is 10.
875 C
876 C Optional Outputs
877 C ----------------
878 C As optional additional output from DLSODE, the variables listed
879 C below are quantities related to the performance of DLSODE which
880 C are available to the user. These are communicated by way of the
881 C work arrays, but also have internal mnemonic names as shown.
882 C Except where stated otherwise, all of these outputs are defined on
883 C any successful return from DLSODE, and on any return with ISTATE =
884 C -1, -2, -4, -5, or -6. On an illegal input return (ISTATE = -3),
885 C they will be unchanged from their existing values (if any), except
886 C possibly for TOLSF, LENRW, and LENIW. On any error return,
887 C outputs relevant to the error will be defined, as noted below.
888 C
889 C Name Location Meaning
890 C ----- --------- ------------------------------------------------
891 C HU RWORK(11) Step size in t last used (successfully).
892 C HCUR RWORK(12) Step size to be attempted on the next step.
893 C TCUR RWORK(13) Current value of the independent variable which
894 C the solver has actually reached, i.e., the
895 C current internal mesh point in t. On output,
896 C TCUR will always be at least as far as the
897 C argument T, but may be farther (if interpolation
898 C was done).
899 C TOLSF RWORK(14) Tolerance scale factor, greater than 1.0,
900 C computed when a request for too much accuracy
901 C was detected (ISTATE = -3 if detected at the
902 C start of the problem, ISTATE = -2 otherwise).
903 C If ITOL is left unaltered but RTOL and ATOL are
904 C uniformly scaled up by a factor of TOLSF for the
905 C next call, then the solver is deemed likely to
906 C succeed. (The user may also ignore TOLSF and
907 C alter the tolerance parameters in any other way
908 C appropriate.)
909 C NST IWORK(11) Number of steps taken for the problem so far.
910 C NFE IWORK(12) Number of F evaluations for the problem so far.
911 C NJE IWORK(13) Number of Jacobian evaluations (and of matrix LU
912 C decompositions) for the problem so far.
913 C NQU IWORK(14) Method order last used (successfully).
914 C NQCUR IWORK(15) Order to be attempted on the next step.
915 C IMXER IWORK(16) Index of the component of largest magnitude in
916 C the weighted local error vector ( e(i)/EWT(i) ),
917 C on an error return with ISTATE = -4 or -5.
918 C LENRW IWORK(17) Length of RWORK actually required. This is
919 C defined on normal returns and on an illegal
920 C input return for insufficient storage.
921 C LENIW IWORK(18) Length of IWORK actually required. This is
922 C defined on normal returns and on an illegal
923 C input return for insufficient storage.
924 C
925 C The following two arrays are segments of the RWORK array which may
926 C also be of interest to the user as optional outputs. For each
927 C array, the table below gives its internal name, its base address
928 C in RWORK, and its description.
929 C
930 C Name Base address Description
931 C ---- ------------ ----------------------------------------------
932 C YH 21 The Nordsieck history array, of size NYH by
933 C (NQCUR + 1), where NYH is the initial value of
934 C NEQ. For j = 0,1,...,NQCUR, column j + 1 of
935 C YH contains HCUR**j/factorial(j) times the jth
936 C derivative of the interpolating polynomial
937 C currently representing the solution, evaluated
938 C at t = TCUR.
939 C ACOR LENRW-NEQ+1 Array of size NEQ used for the accumulated
940 C corrections on each step, scaled on output to
941 C represent the estimated local error in Y on
942 C the last step. This is the vector e in the
943 C description of the error control. It is
944 C defined only on successful return from DLSODE.
945 C
946 C
947 C Part 2. Other Callable Routines
948 C --------------------------------
949 C
950 C The following are optional calls which the user may make to gain
951 C additional capabilities in conjunction with DLSODE.
952 C
953 C Form of call Function
954 C ------------------------ ----------------------------------------
955 C CALL XSETUN(LUN) Set the logical unit number, LUN, for
956 C output of messages from DLSODE, if the
957 C default is not desired. The default
958 C value of LUN is 6. This call may be made
959 C at any time and will take effect
960 C immediately.
961 C CALL XSETF(MFLAG) Set a flag to control the printing of
962 C messages by DLSODE. MFLAG = 0 means do
963 C not print. (Danger: this risks losing
964 C valuable information.) MFLAG = 1 means
965 C print (the default). This call may be
966 C made at any time and will take effect
967 C immediately.
968 C CALL DSRCOM(RSAV,ISAV,JOB) Saves and restores the contents of the
969 C internal COMMON blocks used by DLSODE
970 C (see Part 3 below). RSAV must be a
971 C real array of length 218 or more, and
972 C ISAV must be an integer array of length
973 C 37 or more. JOB = 1 means save COMMON
974 C into RSAV/ISAV. JOB = 2 means restore
975 C COMMON from same. DSRCOM is useful if
976 C one is interrupting a run and restarting
977 C later, or alternating between two or
978 C more problems solved with DLSODE.
979 C CALL DINTDY(,,,,,) Provide derivatives of y, of various
980 C (see below) orders, at a specified point t, if
981 C desired. It may be called only after a
982 C successful return from DLSODE. Detailed
983 C instructions follow.
984 C
985 C Detailed instructions for using DINTDY
986 C --------------------------------------
987 C The form of the CALL is:
988 C
989 C CALL DINTDY (T, K, RWORK(21), NYH, DKY, IFLAG)
990 C
991 C The input parameters are:
992 C
993 C T Value of independent variable where answers are
994 C desired (normally the same as the T last returned by
995 C DLSODE). For valid results, T must lie between
996 C TCUR - HU and TCUR. (See "Optional Outputs" above
997 C for TCUR and HU.)
998 C K Integer order of the derivative desired. K must
999 C satisfy 0 <= K <= NQCUR, where NQCUR is the current
1000 C order (see "Optional Outputs"). The capability
1001 C corresponding to K = 0, i.e., computing y(t), is
1002 C already provided by DLSODE directly. Since
1003 C NQCUR >= 1, the first derivative dy/dt is always
1004 C available with DINTDY.
1005 C RWORK(21) The base address of the history array YH.
1006 C NYH Column length of YH, equal to the initial value of NEQ.
1007 C
1008 C The output parameters are:
1009 C
1010 C DKY Real array of length NEQ containing the computed value
1011 C of the Kth derivative of y(t).
1012 C IFLAG Integer flag, returned as 0 if K and T were legal,
1013 C -1 if K was illegal, and -2 if T was illegal.
1014 C On an error return, a message is also written.
1015 C
1016 C
1017 C Part 3. Common Blocks
1018 C ----------------------
1019 C
1020 C If DLSODE is to be used in an overlay situation, the user must
1021 C declare, in the primary overlay, the variables in:
1022 C (1) the call sequence to DLSODE,
1023 C (2) the internal COMMON block /DLS001/, of length 255
1024 C (218 double precision words followed by 37 integer words).
1025 C
1026 C If DLSODE is used on a system in which the contents of internal
1027 C COMMON blocks are not preserved between calls, the user should
1028 C declare the above COMMON block in his main program to insure that
1029 C its contents are preserved.
1030 C
1031 C If the solution of a given problem by DLSODE is to be interrupted
1032 C and then later continued, as when restarting an interrupted run or
1033 C alternating between two or more problems, the user should save,
1034 C following the return from the last DLSODE call prior to the
1035 C interruption, the contents of the call sequence variables and the
1036 C internal COMMON block, and later restore these values before the
1037 C next DLSODE call for that problem. In addition, if XSETUN and/or
1038 C XSETF was called for non-default handling of error messages, then
1039 C these calls must be repeated. To save and restore the COMMON
1040 C block, use subroutine DSRCOM (see Part 2 above).
1041 C
1042 C
1043 C Part 4. Optionally Replaceable Solver Routines
1044 C -----------------------------------------------
1045 C
1046 C Below are descriptions of two routines in the DLSODE package which
1047 C relate to the measurement of errors. Either routine can be
1048 C replaced by a user-supplied version, if desired. However, since
1049 C such a replacement may have a major impact on performance, it
1050 C should be done only when absolutely necessary, and only with great
1051 C caution. (Note: The means by which the package version of a
1052 C routine is superseded by the user's version may be system-
1053 C dependent.)
1054 C
1055 C DEWSET
1056 C ------
1057 C The following subroutine is called just before each internal
1058 C integration step, and sets the array of error weights, EWT, as
1059 C described under ITOL/RTOL/ATOL above:
1060 C
1061 C SUBROUTINE DEWSET (NEQ, ITOL, RTOL, ATOL, YCUR, EWT)
1062 C
1063 C where NEQ, ITOL, RTOL, and ATOL are as in the DLSODE call
1064 C sequence, YCUR contains the current dependent variable vector,
1065 C and EWT is the array of weights set by DEWSET.
1066 C
1067 C If the user supplies this subroutine, it must return in EWT(i)
1068 C (i = 1,...,NEQ) a positive quantity suitable for comparing errors
1069 C in Y(i) to. The EWT array returned by DEWSET is passed to the
1070 C DVNORM routine (see below), and also used by DLSODE in the
1071 C computation of the optional output IMXER, the diagonal Jacobian
1072 C approximation, and the increments for difference quotient
1073 C Jacobians.
1074 C
1075 C In the user-supplied version of DEWSET, it may be desirable to use
1076 C the current values of derivatives of y. Derivatives up to order NQ
1077 C are available from the history array YH, described above under
1078 C optional outputs. In DEWSET, YH is identical to the YCUR array,
1079 C extended to NQ + 1 columns with a column length of NYH and scale
1080 C factors of H**j/factorial(j). On the first call for the problem,
1081 C given by NST = 0, NQ is 1 and H is temporarily set to 1.0.
1082 C NYH is the initial value of NEQ. The quantities NQ, H, and NST
1083 C can be obtained by including in SEWSET the statements:
1084 C DOUBLE PRECISION RLS
1085 C COMMON /DLS001/ RLS(218),ILS(37)
1086 C NQ = ILS(33)
1087 C NST = ILS(34)
1088 C H = RLS(212)
1089 C Thus, for example, the current value of dy/dt can be obtained as
1090 C YCUR(NYH+i)/H (i=1,...,NEQ) (and the division by H is unnecessary
1091 C when NST = 0).
1092 C
1093 C DVNORM
1094 C ------
1095 C DVNORM is a real function routine which computes the weighted
1096 C root-mean-square norm of a vector v:
1097 C
1098 C d = DVNORM (n, v, w)
1099 C
1100 C where:
1101 C n = the length of the vector,
1102 C v = real array of length n containing the vector,
1103 C w = real array of length n containing weights,
1104 C d = SQRT( (1/n) * sum(v(i)*w(i))**2 ).
1105 C
1106 C DVNORM is called with n = NEQ and with w(i) = 1.0/EWT(i), where
1107 C EWT is as set by subroutine DEWSET.
1108 C
1109 C If the user supplies this function, it should return a nonnegative
1110 C value of DVNORM suitable for use in the error control in DLSODE.
1111 C None of the arguments should be altered by DVNORM. For example, a
1112 C user-supplied DVNORM routine might:
1113 C - Substitute a max-norm of (v(i)*w(i)) for the rms-norm, or
1114 C - Ignore some components of v in the norm, with the effect of
1115 C suppressing the error control on those components of Y.
1116 C ---------------------------------------------------------------------
1117 C***ROUTINES CALLED DEWSET, DINTDY, D1MACH, DSTODE, DVNORM, XERRWD
1118 C***COMMON BLOCKS DLS001
1119 C***REVISION HISTORY (YYYYMMDD)
1120 C 19791129 DATE WRITTEN
1121 C 19791213 Minor changes to declarations; DELP init. in STODE.
1122 C 19800118 Treat NEQ as array; integer declarations added throughout;
1123 C minor changes to prologue.
1124 C 19800306 Corrected TESCO(1,NQP1) setting in CFODE.
1125 C 19800519 Corrected access of YH on forced order reduction;
1126 C numerous corrections to prologues and other comments.
1127 C 19800617 In main driver, added loading of SQRT(UROUND) in RWORK;
1128 C minor corrections to main prologue.
1129 C 19800923 Added zero initialization of HU and NQU.
1130 C 19801218 Revised XERRWD routine; minor corrections to main prologue.
1131 C 19810401 Minor changes to comments and an error message.
1132 C 19810814 Numerous revisions: replaced EWT by 1/EWT; used flags
1133 C JCUR, ICF, IERPJ, IERSL between STODE and subordinates;
1134 C added tuning parameters CCMAX, MAXCOR, MSBP, MXNCF;
1135 C reorganized returns from STODE; reorganized type decls.;
1136 C fixed message length in XERRWD; changed default LUNIT to 6;
1137 C changed Common lengths; changed comments throughout.
1138 C 19870330 Major update by ACH: corrected comments throughout;
1139 C removed TRET from Common; rewrote EWSET with 4 loops;
1140 C fixed t test in INTDY; added Cray directives in STODE;
1141 C in STODE, fixed DELP init. and logic around PJAC call;
1142 C combined routines to save/restore Common;
1143 C passed LEVEL = 0 in error message calls (except run abort).
1144 C 19890426 Modified prologue to SLATEC/LDOC format. (FNF)
1145 C 19890501 Many improvements to prologue. (FNF)
1146 C 19890503 A few final corrections to prologue. (FNF)
1147 C 19890504 Minor cosmetic changes. (FNF)
1148 C 19890510 Corrected description of Y in Arguments section. (FNF)
1149 C 19890517 Minor corrections to prologue. (FNF)
1150 C 19920514 Updated with prologue edited 891025 by G. Shaw for manual.
1151 C 19920515 Converted source lines to upper case. (FNF)
1152 C 19920603 Revised XERRWD calls using mixed upper-lower case. (ACH)
1153 C 19920616 Revised prologue comment regarding CFT. (ACH)
1154 C 19921116 Revised prologue comments regarding Common. (ACH).
1155 C 19930326 Added comment about non-reentrancy. (FNF)
1156 C 19930723 Changed D1MACH to DUMACH. (FNF)
1157 C 19930801 Removed ILLIN and NTREP from Common (affects driver logic);
1158 C minor changes to prologue and internal comments;
1159 C changed Hollerith strings to quoted strings;
1160 C changed internal comments to mixed case;
1161 C replaced XERRWD with new version using character type;
1162 C changed dummy dimensions from 1 to *. (ACH)
1163 C 19930809 Changed to generic intrinsic names; changed names of
1164 C subprograms and Common blocks to DLSODE etc. (ACH)
1165 C 19930929 Eliminated use of REAL intrinsic; other minor changes. (ACH)
1166 C 20010412 Removed all 'own' variables from Common block /DLS001/
1167 C (affects declarations in 6 routines). (ACH)
1168 C 20010509 Minor corrections to prologue. (ACH)
1169 C 20031105 Restored 'own' variables to Common block /DLS001/, to
1170 C enable interrupt/restart feature. (ACH)
1171 C 20031112 Added SAVE statements for data-loaded constants.
1172 C
1173 C***END PROLOGUE DLSODE
1174 C
1175 C*Internal Notes:
1176 C
1177 C Other Routines in the DLSODE Package.
1178 C
1179 C In addition to Subroutine DLSODE, the DLSODE package includes the
1180 C following subroutines and function routines:
1181 C DINTDY computes an interpolated value of the y vector at t = TOUT.
1182 C DSTODE is the core integrator, which does one step of the
1183 C integration and the associated error control.
1184 C DCFODE sets all method coefficients and test constants.
1185 C DPREPJ computes and preprocesses the Jacobian matrix J = df/dy
1186 C and the Newton iteration matrix P = I - h*l0*J.
1187 C DSOLSY manages solution of linear system in chord iteration.
1188 C DEWSET sets the error weight vector EWT before each step.
1189 C DVNORM computes the weighted R.M.S. norm of a vector.
1190 C DSRCOM is a user-callable routine to save and restore
1191 C the contents of the internal Common block.
1192 C DGETRF AND DGETRS ARE ROUTINES FROM LAPACK FOR SOLVING FULL
1193 C SYSTEMS OF LINEAR ALGEBRAIC EQUATIONS.
1194 C DGBTRF AND DGBTRS ARE ROUTINES FROM LAPACK FOR SOLVING BANDED
1195 C LINEAR SYSTEMS.
1196 C D1MACH computes the unit roundoff in a machine-independent manner.
1197 C XERRWD, XSETUN, XSETF, IXSAV, IUMACH handle the printing of all
1198 C error messages and warnings. XERRWD is machine-dependent.
1199 C Note: DVNORM, DUMACH, IXSAV, and IUMACH are function routines.
1200 C All the others are subroutines.
1201 C
1202 C**End
1203 C
1204 C Declare externals.
1205  EXTERNAL dprepj, dsolsy
1206  DOUBLE PRECISION D1MACH, DVNORM
1207 C
1208 C Declare all other variables.
1209  INTEGER ILLIN, INIT, LYH, LEWT, LACOR, LSAVF, LWM, LIWM,
1210  1 mxstep, mxhnil, nhnil, ntrep, nslast, nyh,
1211  2 ialth, ipup, lmax, meo, nqnyh, nslp
1212  INTEGER ICF, IERPJ, IERSL, JCUR, JSTART, KFLAG, L, METH, MITER,
1213  1 maxord, maxcor, msbp, mxncf, n, nq, nst, nfe, nje, nqu
1214  INTEGER I, I1, I2, IFLAG, IMXER, KGO, LF0,
1215  1 leniw, lenrw, lenwm, ml, mord, mu, mxhnl0, mxstp0
1216  DOUBLE PRECISION CONIT, CRATE, EL, ELCO, HOLD, RMAX, TESCO,
1217  1 ccmax, el0, h, hmin, hmxi, hu, rc, tn, uround
1218  DOUBLE PRECISION ATOLI, AYI, BIG, EWTI, H0, HMAX, HMX, RH, RTOLI,
1219  1 tcrit, tdist, tnext, tol, tolsf, tp, SIZE, sum, w0
1220  dimension mord(2)
1221  LOGICAL IHIT
1222  CHARACTER*80 MSG
1223  SAVE mord, mxstp0, mxhnl0
1224 C-----------------------------------------------------------------------
1225 C The following internal Common block contains
1226 C (a) variables which are local to any subroutine but whose values must
1227 C be preserved between calls to the routine ("own" variables), and
1228 C (b) variables which are communicated between subroutines.
1229 C The block DLS001 is declared in subroutines DLSODE, DINTDY, DSTODE,
1230 C DPREPJ, and DSOLSY.
1231 C Groups of variables are replaced by dummy arrays in the Common
1232 C declarations in routines where those variables are not used.
1233 C-----------------------------------------------------------------------
1234  COMMON /dls001/ conit, crate, el(13), elco(13,12),
1235  1 hold, rmax, tesco(3,12),
1236  1 ccmax, el0, h, hmin, hmxi, hu, rc, tn, uround,
1237  2 illin, init, lyh, lewt, lacor, lsavf, lwm, liwm,
1238  3 mxstep, mxhnil, nhnil, ntrep, nslast, nyh,
1239  3 ialth, ipup, lmax, meo, nqnyh, nslp,
1240  4 icf, ierpj, iersl, jcur, jstart, kflag, l, meth, miter,
1241  5 maxord, maxcor, msbp, mxncf, n, nq, nst, nfe, nje, nqu
1242 C
1243  DATA mord(1),mord(2)/12,5/, mxstp0/500/, mxhnl0/10/
1244 C-----------------------------------------------------------------------
1245 C Block A.
1246 C This code block is executed on every call.
1247 C It tests ISTATE and ITASK for legality and branches appropriately.
1248 C If ISTATE .GT. 1 but the flag INIT shows that initialization has
1249 C not yet been done, an error return occurs.
1250 C If ISTATE = 1 and TOUT = T, return immediately.
1251 C-----------------------------------------------------------------------
1252 C
1253 C***FIRST EXECUTABLE STATEMENT DLSODE
1254  IF (istate .LT. 1 .OR. istate .GT. 3) GO TO 601
1255  IF (itask .LT. 1 .OR. itask .GT. 5) GO TO 602
1256  IF (istate .EQ. 1) GO TO 10
1257  IF (init .EQ. 0) GO TO 603
1258  IF (istate .EQ. 2) GO TO 200
1259  GO TO 20
1260  10 init = 0
1261  IF (tout .EQ. t) RETURN
1262 C-----------------------------------------------------------------------
1263 C Block B.
1264 C The next code block is executed for the initial call (ISTATE = 1),
1265 C or for a continuation call with parameter changes (ISTATE = 3).
1266 C It contains checking of all inputs and various initializations.
1267 C
1268 C First check legality of the non-optional inputs NEQ, ITOL, IOPT,
1269 C MF, ML, and MU.
1270 C-----------------------------------------------------------------------
1271  20 IF (neq(1) .LE. 0) GO TO 604
1272  IF (istate .EQ. 1) GO TO 25
1273  IF (neq(1) .GT. n) GO TO 605
1274  25 n = neq(1)
1275  IF (itol .LT. 1 .OR. itol .GT. 4) GO TO 606
1276  IF (iopt .LT. 0 .OR. iopt .GT. 1) GO TO 607
1277  meth = mf/10
1278  miter = mf - 10*meth
1279  IF (meth .LT. 1 .OR. meth .GT. 2) GO TO 608
1280  IF (miter .LT. 0 .OR. miter .GT. 5) GO TO 608
1281  IF (miter .LE. 3) GO TO 30
1282  ml = iwork(1)
1283  mu = iwork(2)
1284  IF (ml .LT. 0 .OR. ml .GE. n) GO TO 609
1285  IF (mu .LT. 0 .OR. mu .GE. n) GO TO 610
1286  30 CONTINUE
1287 C Next process and check the optional inputs. --------------------------
1288  IF (iopt .EQ. 1) GO TO 40
1289  maxord = mord(meth)
1290  mxstep = mxstp0
1291  mxhnil = mxhnl0
1292  IF (istate .EQ. 1) h0 = 0.0d0
1293  hmxi = 0.0d0
1294  hmin = 0.0d0
1295  GO TO 60
1296  40 maxord = iwork(5)
1297  IF (maxord .LT. 0) GO TO 611
1298  IF (maxord .EQ. 0) maxord = 100
1299  maxord = min(maxord,mord(meth))
1300  mxstep = iwork(6)
1301  IF (mxstep .LT. 0) GO TO 612
1302  IF (mxstep .EQ. 0) mxstep = mxstp0
1303  mxhnil = iwork(7)
1304  IF (mxhnil .LT. 0) GO TO 613
1305  IF (mxhnil .EQ. 0) mxhnil = mxhnl0
1306  IF (istate .NE. 1) GO TO 50
1307  h0 = rwork(5)
1308  IF ((tout - t)*h0 .LT. 0.0d0) GO TO 614
1309  50 hmax = rwork(6)
1310  IF (hmax .LT. 0.0d0) GO TO 615
1311  hmxi = 0.0d0
1312  IF (hmax .GT. 0.0d0) hmxi = 1.0d0/hmax
1313  hmin = rwork(7)
1314  IF (hmin .LT. 0.0d0) GO TO 616
1315 C-----------------------------------------------------------------------
1316 C Set work array pointers and check lengths LRW and LIW.
1317 C Pointers to segments of RWORK and IWORK are named by prefixing L to
1318 C the name of the segment. E.g., the segment YH starts at RWORK(LYH).
1319 C Segments of RWORK (in order) are denoted YH, WM, EWT, SAVF, ACOR.
1320 C-----------------------------------------------------------------------
1321  60 lyh = 21
1322  IF (istate .EQ. 1) nyh = n
1323  lwm = lyh + (maxord + 1)*nyh
1324  IF (miter .EQ. 0) lenwm = 0
1325  IF (miter .EQ. 1 .OR. miter .EQ. 2) lenwm = n*n + 2
1326  IF (miter .EQ. 3) lenwm = n + 2
1327  IF (miter .GE. 4) lenwm = (2*ml + mu + 1)*n + 2
1328  lewt = lwm + lenwm
1329  lsavf = lewt + n
1330  lacor = lsavf + n
1331  lenrw = lacor + n - 1
1332  iwork(17) = lenrw
1333  liwm = 1
1334  leniw = 20 + n
1335  IF (miter .EQ. 0 .OR. miter .EQ. 3) leniw = 20
1336  iwork(18) = leniw
1337  IF (lenrw .GT. lrw) GO TO 617
1338  IF (leniw .GT. liw) GO TO 618
1339 C Check RTOL and ATOL for legality. ------------------------------------
1340  rtoli = rtol(1)
1341  atoli = atol(1)
1342  DO 70 i = 1,n
1343  IF (itol .GE. 3) rtoli = rtol(i)
1344  IF (itol .EQ. 2 .OR. itol .EQ. 4) atoli = atol(i)
1345  IF (rtoli .LT. 0.0d0) GO TO 619
1346  IF (atoli .LT. 0.0d0) GO TO 620
1347  70 CONTINUE
1348  IF (istate .EQ. 1) GO TO 100
1349 C If ISTATE = 3, set flag to signal parameter changes to DSTODE. -------
1350  jstart = -1
1351  IF (nq .LE. maxord) GO TO 90
1352 C MAXORD was reduced below NQ. Copy YH(*,MAXORD+2) into SAVF. ---------
1353  DO 80 i = 1,n
1354  80 rwork(i+lsavf-1) = rwork(i+lwm-1)
1355 C Reload WM(1) = RWORK(LWM), since LWM may have changed. ---------------
1356  90 IF (miter .GT. 0) rwork(lwm) = sqrt(uround)
1357  IF (n .EQ. nyh) GO TO 200
1358 C NEQ was reduced. Zero part of YH to avoid undefined references. -----
1359  i1 = lyh + l*nyh
1360  i2 = lyh + (maxord + 1)*nyh - 1
1361  IF (i1 .GT. i2) GO TO 200
1362  DO 95 i = i1,i2
1363  95 rwork(i) = 0.0d0
1364  GO TO 200
1365 C-----------------------------------------------------------------------
1366 C Block C.
1367 C The next block is for the initial call only (ISTATE = 1).
1368 C It contains all remaining initializations, the initial call to F,
1369 C and the calculation of the initial step size.
1370 C The error weights in EWT are inverted after being loaded.
1371 C-----------------------------------------------------------------------
1372  100 uround = d1mach(4)
1373  tn = t
1374  IF (itask .NE. 4 .AND. itask .NE. 5) GO TO 110
1375  tcrit = rwork(1)
1376  IF ((tcrit - tout)*(tout - t) .LT. 0.0d0) GO TO 625
1377  IF (h0 .NE. 0.0d0 .AND. (t + h0 - tcrit)*h0 .GT. 0.0d0)
1378  1 h0 = tcrit - t
1379  110 jstart = 0
1380  IF (miter .GT. 0) rwork(lwm) = sqrt(uround)
1381  nhnil = 0
1382  nst = 0
1383  nje = 0
1384  nslast = 0
1385  hu = 0.0d0
1386  nqu = 0
1387  ccmax = 0.3d0
1388  maxcor = 3
1389  msbp = 20
1390  mxncf = 10
1391 C Initial call to F. (LF0 points to YH(*,2).) -------------------------
1392  lf0 = lyh + nyh
1393  CALL f (neq, t, y, rwork(lf0))
1394  nfe = 1
1395 C Load the initial value vector in YH. ---------------------------------
1396  DO 115 i = 1,n
1397  115 rwork(i+lyh-1) = y(i)
1398 C Load and invert the EWT array. (H is temporarily set to 1.0.) -------
1399  nq = 1
1400  h = 1.0d0
1401  CALL dewset (n, itol, rtol, atol, rwork(lyh), rwork(lewt))
1402  DO 120 i = 1,n
1403  IF (rwork(i+lewt-1) .LE. 0.0d0) GO TO 621
1404  120 rwork(i+lewt-1) = 1.0d0/rwork(i+lewt-1)
1405 C-----------------------------------------------------------------------
1406 C The coding below computes the step size, H0, to be attempted on the
1407 C first step, unless the user has supplied a value for this.
1408 C First check that TOUT - T differs significantly from zero.
1409 C A scalar tolerance quantity TOL is computed, as MAX(RTOL(I))
1410 C if this is positive, or MAX(ATOL(I)/ABS(Y(I))) otherwise, adjusted
1411 C so as to be between 100*UROUND and 1.0E-3.
1412 C Then the computed value H0 is given by..
1413 C NEQ
1414 C H0**2 = TOL / ( w0**-2 + (1/NEQ) * SUM ( f(i)/ywt(i) )**2 )
1415 C 1
1416 C where w0 = MAX ( ABS(T), ABS(TOUT) ),
1417 C f(i) = i-th component of initial value of f,
1418 C ywt(i) = EWT(i)/TOL (a weight for y(i)).
1419 C The sign of H0 is inferred from the initial values of TOUT and T.
1420 C-----------------------------------------------------------------------
1421  IF (h0 .NE. 0.0d0) GO TO 180
1422  tdist = abs(tout - t)
1423  w0 = max(abs(t),abs(tout))
1424  IF (tdist .LT. 2.0d0*uround*w0) GO TO 622
1425  tol = rtol(1)
1426  IF (itol .LE. 2) GO TO 140
1427  DO 130 i = 1,n
1428  130 tol = max(tol,rtol(i))
1429  140 IF (tol .GT. 0.0d0) GO TO 160
1430  atoli = atol(1)
1431  DO 150 i = 1,n
1432  IF (itol .EQ. 2 .OR. itol .EQ. 4) atoli = atol(i)
1433  ayi = abs(y(i))
1434  IF (ayi .NE. 0.0d0) tol = max(tol,atoli/ayi)
1435  150 CONTINUE
1436  160 tol = max(tol,100.0d0*uround)
1437  tol = min(tol,0.001d0)
1438  sum = dvnorm(n, rwork(lf0), rwork(lewt))
1439  sum = 1.0d0/(tol*w0*w0) + tol*sum**2
1440  h0 = 1.0d0/sqrt(sum)
1441  h0 = min(h0,tdist)
1442  h0 = sign(h0,tout-t)
1443 C Adjust H0 if necessary to meet HMAX bound. ---------------------------
1444  180 rh = abs(h0)*hmxi
1445  IF (rh .GT. 1.0d0) h0 = h0/rh
1446 C Load H with H0 and scale YH(*,2) by H0. ------------------------------
1447  h = h0
1448  DO 190 i = 1,n
1449  190 rwork(i+lf0-1) = h0*rwork(i+lf0-1)
1450  GO TO 270
1451 C-----------------------------------------------------------------------
1452 C Block D.
1453 C The next code block is for continuation calls only (ISTATE = 2 or 3)
1454 C and is to check stop conditions before taking a step.
1455 C-----------------------------------------------------------------------
1456  200 nslast = nst
1457  GO TO (210, 250, 220, 230, 240), itask
1458  210 IF ((tn - tout)*h .LT. 0.0d0) GO TO 250
1459  CALL dintdy (tout, 0, rwork(lyh), nyh, y, iflag)
1460  IF (iflag .NE. 0) GO TO 627
1461  t = tout
1462  GO TO 420
1463  220 tp = tn - hu*(1.0d0 + 100.0d0*uround)
1464  IF ((tp - tout)*h .GT. 0.0d0) GO TO 623
1465  IF ((tn - tout)*h .LT. 0.0d0) GO TO 250
1466  GO TO 400
1467  230 tcrit = rwork(1)
1468  IF ((tn - tcrit)*h .GT. 0.0d0) GO TO 624
1469  IF ((tcrit - tout)*h .LT. 0.0d0) GO TO 625
1470  IF ((tn - tout)*h .LT. 0.0d0) GO TO 245
1471  CALL dintdy (tout, 0, rwork(lyh), nyh, y, iflag)
1472  IF (iflag .NE. 0) GO TO 627
1473  t = tout
1474  GO TO 420
1475  240 tcrit = rwork(1)
1476  IF ((tn - tcrit)*h .GT. 0.0d0) GO TO 624
1477  245 hmx = abs(tn) + abs(h)
1478  ihit = abs(tn - tcrit) .LE. 100.0d0*uround*hmx
1479  IF (ihit) GO TO 400
1480  tnext = tn + h*(1.0d0 + 4.0d0*uround)
1481  IF ((tnext - tcrit)*h .LE. 0.0d0) GO TO 250
1482  h = (tcrit - tn)*(1.0d0 - 4.0d0*uround)
1483  IF (istate .EQ. 2) jstart = -2
1484 C-----------------------------------------------------------------------
1485 C Block E.
1486 C The next block is normally executed for all calls and contains
1487 C the call to the one-step core integrator DSTODE.
1488 C
1489 C This is a looping point for the integration steps.
1490 C
1491 C First check for too many steps being taken, update EWT (if not at
1492 C start of problem), check for too much accuracy being requested, and
1493 C check for H below the roundoff level in T.
1494 C-----------------------------------------------------------------------
1495  250 CONTINUE
1496  IF ((nst-nslast) .GE. mxstep) GO TO 500
1497  CALL dewset (n, itol, rtol, atol, rwork(lyh), rwork(lewt))
1498  DO 260 i = 1,n
1499  IF (rwork(i+lewt-1) .LE. 0.0d0) GO TO 510
1500  260 rwork(i+lewt-1) = 1.0d0/rwork(i+lewt-1)
1501  270 tolsf = uround*dvnorm(n, rwork(lyh), rwork(lewt))
1502  IF (tolsf .LE. 1.0d0) GO TO 280
1503  tolsf = tolsf*2.0d0
1504  IF (nst .EQ. 0) GO TO 626
1505  GO TO 520
1506  280 IF ((tn + h) .NE. tn) GO TO 290
1507  nhnil = nhnil + 1
1508  IF (nhnil .GT. mxhnil) GO TO 290
1509  msg = 'DLSODE- Warning..internal T (=R1) and H (=R2) are'
1510  CALL xerrwd (msg, 50, 101, 0, 0, 0, 0, 0, 0.0d0, 0.0d0)
1511  msg=' such that in the machine, T + H = T on the next step '
1512  CALL xerrwd (msg, 60, 101, 0, 0, 0, 0, 0, 0.0d0, 0.0d0)
1513  msg = ' (H = step size). Solver will continue anyway'
1514  CALL xerrwd (msg, 50, 101, 0, 0, 0, 0, 2, tn, h)
1515  IF (nhnil .LT. mxhnil) GO TO 290
1516  msg = 'DLSODE- Above warning has been issued I1 times. '
1517  CALL xerrwd (msg, 50, 102, 0, 0, 0, 0, 0, 0.0d0, 0.0d0)
1518  msg = ' It will not be issued again for this problem'
1519  CALL xerrwd (msg, 50, 102, 0, 1, mxhnil, 0, 0, 0.0d0, 0.0d0)
1520  290 CONTINUE
1521 C-----------------------------------------------------------------------
1522 C CALL DSTODE(NEQ,Y,YH,NYH,YH,EWT,SAVF,ACOR,WM,IWM,F,JAC,DPREPJ,DSOLSY)
1523 C-----------------------------------------------------------------------
1524  CALL dstode (neq, y, rwork(lyh), nyh, rwork(lyh), rwork(lewt),
1525  1 rwork(lsavf), rwork(lacor), rwork(lwm), iwork(liwm),
1526  2 f, jac, dprepj, dsolsy)
1527  kgo = 1 - kflag
1528  GO TO (300, 530, 540), kgo
1529 C-----------------------------------------------------------------------
1530 C Block F.
1531 C The following block handles the case of a successful return from the
1532 C core integrator (KFLAG = 0). Test for stop conditions.
1533 C-----------------------------------------------------------------------
1534  300 init = 1
1535  GO TO (310, 400, 330, 340, 350), itask
1536 C ITASK = 1. If TOUT has been reached, interpolate. -------------------
1537  310 IF ((tn - tout)*h .LT. 0.0d0) GO TO 250
1538  CALL dintdy (tout, 0, rwork(lyh), nyh, y, iflag)
1539  t = tout
1540  GO TO 420
1541 C ITASK = 3. Jump to exit if TOUT was reached. ------------------------
1542  330 IF ((tn - tout)*h .GE. 0.0d0) GO TO 400
1543  GO TO 250
1544 C ITASK = 4. See if TOUT or TCRIT was reached. Adjust H if necessary.
1545  340 IF ((tn - tout)*h .LT. 0.0d0) GO TO 345
1546  CALL dintdy (tout, 0, rwork(lyh), nyh, y, iflag)
1547  t = tout
1548  GO TO 420
1549  345 hmx = abs(tn) + abs(h)
1550  ihit = abs(tn - tcrit) .LE. 100.0d0*uround*hmx
1551  IF (ihit) GO TO 400
1552  tnext = tn + h*(1.0d0 + 4.0d0*uround)
1553  IF ((tnext - tcrit)*h .LE. 0.0d0) GO TO 250
1554  h = (tcrit - tn)*(1.0d0 - 4.0d0*uround)
1555  jstart = -2
1556  GO TO 250
1557 C ITASK = 5. See if TCRIT was reached and jump to exit. ---------------
1558  350 hmx = abs(tn) + abs(h)
1559  ihit = abs(tn - tcrit) .LE. 100.0d0*uround*hmx
1560 C-----------------------------------------------------------------------
1561 C Block G.
1562 C The following block handles all successful returns from DLSODE.
1563 C If ITASK .NE. 1, Y is loaded from YH and T is set accordingly.
1564 C ISTATE is set to 2, and the optional outputs are loaded into the
1565 C work arrays before returning.
1566 C-----------------------------------------------------------------------
1567  400 DO 410 i = 1,n
1568  410 y(i) = rwork(i+lyh-1)
1569  t = tn
1570  IF (itask .NE. 4 .AND. itask .NE. 5) GO TO 420
1571  IF (ihit) t = tcrit
1572  420 istate = 2
1573  rwork(11) = hu
1574  rwork(12) = h
1575  rwork(13) = tn
1576  iwork(11) = nst
1577  iwork(12) = nfe
1578  iwork(13) = nje
1579  iwork(14) = nqu
1580  iwork(15) = nq
1581  RETURN
1582 C-----------------------------------------------------------------------
1583 C Block H.
1584 C The following block handles all unsuccessful returns other than
1585 C those for illegal input. First the error message routine is called.
1586 C If there was an error test or convergence test failure, IMXER is set.
1587 C Then Y is loaded from YH and T is set to TN. The optional outputs
1588 C are loaded into the work arrays before returning.
1589 C-----------------------------------------------------------------------
1590 C The maximum number of steps was taken before reaching TOUT. ----------
1591  500 msg = 'DLSODE- At current T (=R1), MXSTEP (=I1) steps '
1592  CALL xerrwd (msg, 50, 201, 0, 0, 0, 0, 0, 0.0d0, 0.0d0)
1593  msg = ' taken on this call before reaching TOUT '
1594  CALL xerrwd (msg, 50, 201, 0, 1, mxstep, 0, 1, tn, 0.0d0)
1595  istate = -1
1596  GO TO 580
1597 C EWT(I) .LE. 0.0 for some I (not at start of problem). ----------------
1598  510 ewti = rwork(lewt+i-1)
1599  msg = .LE.'DLSODE- At T (=R1), EWT(I1) has become R2 0.'
1600  CALL xerrwd (msg, 50, 202, 0, 1, i, 0, 2, tn, ewti)
1601  istate = -6
1602  GO TO 580
1603 C Too much accuracy requested for machine precision. -------------------
1604  520 msg = 'DLSODE- At T (=R1), too much accuracy requested '
1605  CALL xerrwd (msg, 50, 203, 0, 0, 0, 0, 0, 0.0d0, 0.0d0)
1606  msg = ' for precision of machine.. see TOLSF (=R2) '
1607  CALL xerrwd (msg, 50, 203, 0, 0, 0, 0, 2, tn, tolsf)
1608  rwork(14) = tolsf
1609  istate = -2
1610  GO TO 580
1611 C KFLAG = -1. Error test failed repeatedly or with ABS(H) = HMIN. -----
1612  530 msg = 'DLSODE- At T(=R1) and step size H(=R2), the error'
1613  CALL xerrwd (msg, 50, 204, 0, 0, 0, 0, 0, 0.0d0, 0.0d0)
1614  msg = ' test failed repeatedly or with ABS(H) = HMIN'
1615  CALL xerrwd (msg, 50, 204, 0, 0, 0, 0, 2, tn, h)
1616  istate = -4
1617  GO TO 560
1618 C KFLAG = -2. Convergence failed repeatedly or with ABS(H) = HMIN. ----
1619  540 msg = 'DLSODE- At T (=R1) and step size H (=R2), the '
1620  CALL xerrwd (msg, 50, 205, 0, 0, 0, 0, 0, 0.0d0, 0.0d0)
1621  msg = ' corrector convergence failed repeatedly '
1622  CALL xerrwd (msg, 50, 205, 0, 0, 0, 0, 0, 0.0d0, 0.0d0)
1623  msg = ' or with ABS(H) = HMIN '
1624  CALL xerrwd (msg, 30, 205, 0, 0, 0, 0, 2, tn, h)
1625  istate = -5
1626 C Compute IMXER if relevant. -------------------------------------------
1627  560 big = 0.0d0
1628  imxer = 1
1629  DO 570 i = 1,n
1630  SIZE = abs(rwork(i+lacor-1)*rwork(i+lewt-1))
1631  IF (big .GE. size) GO TO 570
1632  big = SIZE
1633  imxer = i
1634  570 CONTINUE
1635  iwork(16) = imxer
1636 C Set Y vector, T, and optional outputs. -------------------------------
1637  580 DO 590 i = 1,n
1638  590 y(i) = rwork(i+lyh-1)
1639  t = tn
1640  rwork(11) = hu
1641  rwork(12) = h
1642  rwork(13) = tn
1643  iwork(11) = nst
1644  iwork(12) = nfe
1645  iwork(13) = nje
1646  iwork(14) = nqu
1647  iwork(15) = nq
1648  RETURN
1649 C-----------------------------------------------------------------------
1650 C Block I.
1651 C The following block handles all error returns due to illegal input
1652 C (ISTATE = -3), as detected before calling the core integrator.
1653 C First the error message routine is called. If the illegal input
1654 C is a negative ISTATE, the run is aborted (apparent infinite loop).
1655 C-----------------------------------------------------------------------
1656  601 msg = 'DLSODE- ISTATE (=I1) illegal '
1657  CALL xerrwd (msg, 30, 1, 0, 1, istate, 0, 0, 0.0d0, 0.0d0)
1658  IF (istate .LT. 0) GO TO 800
1659  GO TO 700
1660  602 msg = 'DLSODE- ITASK (=I1) illegal '
1661  CALL xerrwd (msg, 30, 2, 0, 1, itask, 0, 0, 0.0d0, 0.0d0)
1662  GO TO 700
1663  603 msg = .GT.'DLSODE- ISTATE 1 but DLSODE not initialized '
1664  CALL xerrwd (msg, 50, 3, 0, 0, 0, 0, 0, 0.0d0, 0.0d0)
1665  GO TO 700
1666  604 msg = .LT.'DLSODE- NEQ (=I1) 1 '
1667  CALL xerrwd (msg, 30, 4, 0, 1, neq(1), 0, 0, 0.0d0, 0.0d0)
1668  GO TO 700
1669  605 msg = 'DLSODE- ISTATE = 3 and NEQ increased (I1 to I2) '
1670  CALL xerrwd (msg, 50, 5, 0, 2, n, neq(1), 0, 0.0d0, 0.0d0)
1671  GO TO 700
1672  606 msg = 'DLSODE- ITOL (=I1) illegal '
1673  CALL xerrwd (msg, 30, 6, 0, 1, itol, 0, 0, 0.0d0, 0.0d0)
1674  GO TO 700
1675  607 msg = 'DLSODE- IOPT (=I1) illegal '
1676  CALL xerrwd (msg, 30, 7, 0, 1, iopt, 0, 0, 0.0d0, 0.0d0)
1677  GO TO 700
1678  608 msg = 'DLSODE- MF (=I1) illegal '
1679  CALL xerrwd (msg, 30, 8, 0, 1, mf, 0, 0, 0.0d0, 0.0d0)
1680  GO TO 700
1681  609 msg = .LT..GE.'DLSODE- ML (=I1) illegal.. 0 or NEQ (=I2)'
1682  CALL xerrwd (msg, 50, 9, 0, 2, ml, neq(1), 0, 0.0d0, 0.0d0)
1683  GO TO 700
1684  610 msg = .LT..GE.'DLSODE- MU (=I1) illegal.. 0 or NEQ (=I2)'
1685  CALL xerrwd (msg, 50, 10, 0, 2, mu, neq(1), 0, 0.0d0, 0.0d0)
1686  GO TO 700
1687  611 msg = .LT.'DLSODE- MAXORD (=I1) 0 '
1688  CALL xerrwd (msg, 30, 11, 0, 1, maxord, 0, 0, 0.0d0, 0.0d0)
1689  GO TO 700
1690  612 msg = .LT.'DLSODE- MXSTEP (=I1) 0 '
1691  CALL xerrwd (msg, 30, 12, 0, 1, mxstep, 0, 0, 0.0d0, 0.0d0)
1692  GO TO 700
1693  613 msg = .LT.'DLSODE- MXHNIL (=I1) 0 '
1694  CALL xerrwd (msg, 30, 13, 0, 1, mxhnil, 0, 0, 0.0d0, 0.0d0)
1695  GO TO 700
1696  614 msg = 'DLSODE- TOUT (=R1) behind T (=R2) '
1697  CALL xerrwd (msg, 40, 14, 0, 0, 0, 0, 2, tout, t)
1698  msg = ' Integration direction is given by H0 (=R1) '
1699  CALL xerrwd (msg, 50, 14, 0, 0, 0, 0, 1, h0, 0.0d0)
1700  GO TO 700
1701  615 msg = .LT.'DLSODE- HMAX (=R1) 0.0 '
1702  CALL xerrwd (msg, 30, 15, 0, 0, 0, 0, 1, hmax, 0.0d0)
1703  GO TO 700
1704  616 msg = .LT.'DLSODE- HMIN (=R1) 0.0 '
1705  CALL xerrwd (msg, 30, 16, 0, 0, 0, 0, 1, hmin, 0.0d0)
1706  GO TO 700
1707  617 CONTINUE
1708  msg='DLSODE- RWORK length needed, LENRW (=I1), exceeds LRW (=I2)'
1709  CALL xerrwd (msg, 60, 17, 0, 2, lenrw, lrw, 0, 0.0d0, 0.0d0)
1710  GO TO 700
1711  618 CONTINUE
1712  msg='DLSODE- IWORK length needed, LENIW (=I1), exceeds LIW (=I2)'
1713  CALL xerrwd (msg, 60, 18, 0, 2, leniw, liw, 0, 0.0d0, 0.0d0)
1714  GO TO 700
1715  619 msg = .LT.'DLSODE- RTOL(I1) is R1 0.0 '
1716  CALL xerrwd (msg, 40, 19, 0, 1, i, 0, 1, rtoli, 0.0d0)
1717  GO TO 700
1718  620 msg = .LT.'DLSODE- ATOL(I1) is R1 0.0 '
1719  CALL xerrwd (msg, 40, 20, 0, 1, i, 0, 1, atoli, 0.0d0)
1720  GO TO 700
1721  621 ewti = rwork(lewt+i-1)
1722  msg = .LE.'DLSODE- EWT(I1) is R1 0.0 '
1723  CALL xerrwd (msg, 40, 21, 0, 1, i, 0, 1, ewti, 0.0d0)
1724  GO TO 700
1725  622 CONTINUE
1726  msg='DLSODE- TOUT (=R1) too close to T(=R2) to start integration'
1727  CALL xerrwd (msg, 60, 22, 0, 0, 0, 0, 2, tout, t)
1728  GO TO 700
1729  623 CONTINUE
1730  msg='DLSODE- ITASK = I1 and TOUT (=R1) behind TCUR - HU (= R2) '
1731  CALL xerrwd (msg, 60, 23, 0, 1, itask, 0, 2, tout, tp)
1732  GO TO 700
1733  624 CONTINUE
1734  msg='DLSODE- ITASK = 4 OR 5 and TCRIT (=R1) behind TCUR (=R2) '
1735  CALL xerrwd (msg, 60, 24, 0, 0, 0, 0, 2, tcrit, tn)
1736  GO TO 700
1737  625 CONTINUE
1738  msg='DLSODE- ITASK = 4 or 5 and TCRIT (=R1) behind TOUT (=R2) '
1739  CALL xerrwd (msg, 60, 25, 0, 0, 0, 0, 2, tcrit, tout)
1740  GO TO 700
1741  626 msg = 'DLSODE- At start of problem, too much accuracy '
1742  CALL xerrwd (msg, 50, 26, 0, 0, 0, 0, 0, 0.0d0, 0.0d0)
1743  msg=' requested for precision of machine.. See TOLSF (=R1) '
1744  CALL xerrwd (msg, 60, 26, 0, 0, 0, 0, 1, tolsf, 0.0d0)
1745  rwork(14) = tolsf
1746  GO TO 700
1747  627 msg = 'DLSODE- Trouble in DINTDY. ITASK = I1, TOUT = R1'
1748  CALL xerrwd (msg, 50, 27, 0, 1, itask, 0, 1, tout, 0.0d0)
1749 C
1750  700 istate = -3
1751  RETURN
1752 C
1753  800 msg = 'DLSODE- Run aborted.. apparent infinite loop '
1754  CALL xerrwd (msg, 50, 303, 2, 0, 0, 0, 0, 0.0d0, 0.0d0)
1755  RETURN
1756 C----------------------- END OF SUBROUTINE DLSODE ----------------------
1757  END
charNDArray max(char d, const charNDArray &m)
Definition: chNDArray.cc:230
charNDArray min(char d, const charNDArray &m)
Definition: chNDArray.cc:207
subroutine dewset(N, ITOL, RTOL, ATOL, YCUR, EWT)
Definition: dewset.f:2
subroutine dintdy(T, K, YH, NYH, DKY, IFLAG)
Definition: dintdy.f:2
subroutine dlsode(F, NEQ, Y, T, TOUT, ITOL, RTOL, ATOL, ITASK, ISTATE, IOPT, RWORK, LRW, IWORK, LIW, JAC, MF)
Definition: dlsode.f:3
subroutine dprepj(NEQ, Y, YH, NYH, EWT, FTEM, SAVF, WM, IWM, F, JAC)
Definition: dprepj.f:3
subroutine dsolsy(WM, IWM, X, TEM)
Definition: dsolsy.f:2
subroutine dstode(NEQ, Y, YH, NYH, YH1, EWT, SAVF, ACOR, WM, IWM, F, JAC, PJAC, SLVS)
Definition: dstode.f:3
subroutine xerrwd(MSG, NMES, NERR, LEVEL, NI, I1, I2, NR, R1, R2)
Definition: xerrwd.f:4