GNU Octave 10.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 
Loading...
Searching...
No Matches
ddaspk.f
Go to the documentation of this file.
1C Work performed under the auspices of the U.S. Department of Energy
2C by Lawrence Livermore National Laboratory under contract number
3C W-7405-Eng-48.
4C
5 SUBROUTINE ddaspk (RES, NEQ, T, Y, YPRIME, TOUT, INFO, RTOL, ATOL,
6 * IDID, RWORK, LRW, IWORK, LIW, RPAR, IPAR, JAC, PSOL)
7C
8C***BEGIN PROLOGUE DDASPK
9C***DATE WRITTEN 890101 (YYMMDD)
10C***REVISION DATE 910624
11C***REVISION DATE 920929 (CJ in RES call, RES counter fix.)
12C***REVISION DATE 921215 (Warnings on poor iteration performance)
13C***REVISION DATE 921216 (NRMAX as optional input)
14C***REVISION DATE 930315 (Name change: DDINI to DDINIT)
15C***REVISION DATE 940822 (Replaced initial condition calculation)
16C***REVISION DATE 941101 (Added linesearch in I.C. calculations)
17C***REVISION DATE 941220 (Misc. corrections throughout)
18C***REVISION DATE 950125 (Added DINVWT routine)
19C***REVISION DATE 950714 (Misc. corrections throughout)
20C***REVISION DATE 950802 (Default NRMAX = 5, based on tests.)
21C***REVISION DATE 950808 (Optional error test added.)
22C***REVISION DATE 950814 (Added I.C. constraints and INFO(14))
23C***REVISION DATE 950828 (Various minor corrections.)
24C***REVISION DATE 951006 (Corrected WT scaling in DFNRMK.)
25C***REVISION DATE 960129 (Corrected RL bug in DLINSD, DLINSK.)
26C***REVISION DATE 960301 (Added NONNEG to SAVE statement.)
27C***CATEGORY NO. I1A2
28C***KEYWORDS DIFFERENTIAL/ALGEBRAIC, BACKWARD DIFFERENTIATION FORMULAS,
29C IMPLICIT DIFFERENTIAL SYSTEMS, KRYLOV ITERATION
30C***AUTHORS Linda R. Petzold, Peter N. Brown, Alan C. Hindmarsh, and
31C Clement W. Ulrich
32C Center for Computational Sciences & Engineering, L-316
33C Lawrence Livermore National Laboratory
34C P.O. Box 808,
35C Livermore, CA 94551
36C***PURPOSE This code solves a system of differential/algebraic
37C equations of the form
38C G(t,y,y') = 0 ,
39C using a combination of Backward Differentiation Formula
40C (BDF) methods and a choice of two linear system solution
41C methods: direct (dense or band) or Krylov (iterative).
42C This version is in double precision.
43C-----------------------------------------------------------------------
44C***DESCRIPTION
45C
46C *Usage:
47C
48C IMPLICIT DOUBLE PRECISION(A-H,O-Z)
49C INTEGER NEQ, INFO(N), IDID, LRW, LIW, IWORK(LIW), IPAR(*)
50C DOUBLE PRECISION T, Y(*), YPRIME(*), TOUT, RTOL(*), ATOL(*),
51C RWORK(LRW), RPAR(*)
52C EXTERNAL RES, JAC, PSOL
53C
54C CALL DDASPK (RES, NEQ, T, Y, YPRIME, TOUT, INFO, RTOL, ATOL,
55C * IDID, RWORK, LRW, IWORK, LIW, RPAR, IPAR, JAC, PSOL)
56C
57C Quantities which may be altered by the code are:
58C T, Y(*), YPRIME(*), INFO(*), RTOL, ATOL, IDID, RWORK(*), IWORK(*)
59C
60C
61C *Arguments:
62C
63C RES:EXT This is the name of a subroutine which you
64C provide to define the residual function G(t,y,y')
65C of the differential/algebraic system.
66C
67C NEQ:IN This is the number of equations in the system.
68C
69C T:INOUT This is the current value of the independent
70C variable.
71C
72C Y(*):INOUT This array contains the solution components at T.
73C
74C YPRIME(*):INOUT This array contains the derivatives of the solution
75C components at T.
76C
77C TOUT:IN This is a point at which a solution is desired.
78C
79C INFO(N):IN This is an integer array used to communicate details
80C of how the solution is to be carried out, such as
81C tolerance type, matrix structure, step size and
82C order limits, and choice of nonlinear system method.
83C N must be at least 20.
84C
85C RTOL,ATOL:INOUT These quantities represent absolute and relative
86C error tolerances (on local error) which you provide
87C to indicate how accurately you wish the solution to
88C be computed. You may choose them to be both scalars
89C or else both arrays of length NEQ.
90C
91C IDID:OUT This integer scalar is an indicator reporting what
92C the code did. You must monitor this variable to
93C decide what action to take next.
94C
95C RWORK:WORK A real work array of length LRW which provides the
96C code with needed storage space.
97C
98C LRW:IN The length of RWORK.
99C
100C IWORK:WORK An integer work array of length LIW which provides
101C the code with needed storage space.
102C
103C LIW:IN The length of IWORK.
104C
105C RPAR,IPAR:IN These are real and integer parameter arrays which
106C you can use for communication between your calling
107C program and the RES, JAC, and PSOL subroutines.
108C
109C JAC:EXT This is the name of a subroutine which you may
110C provide (optionally) for calculating Jacobian
111C (partial derivative) data involved in solving linear
112C systems within DDASPK.
113C
114C PSOL:EXT This is the name of a subroutine which you must
115C provide for solving linear systems if you selected
116C a Krylov method. The purpose of PSOL is to solve
117C linear systems involving a left preconditioner P.
118C
119C *Overview
120C
121C The DDASPK solver uses the backward differentiation formulas of
122C orders one through five to solve a system of the form G(t,y,y') = 0
123C for y = Y and y' = YPRIME. Values for Y and YPRIME at the initial
124C time must be given as input. These values should be consistent,
125C that is, if T, Y, YPRIME are the given initial values, they should
126C satisfy G(T,Y,YPRIME) = 0. However, if consistent values are not
127C known, in many cases you can have DDASPK solve for them -- see INFO(11).
128C (This and other options are described in more detail below.)
129C
130C Normally, DDASPK solves the system from T to TOUT. It is easy to
131C continue the solution to get results at additional TOUT. This is
132C the interval mode of operation. Intermediate results can also be
133C obtained easily by specifying INFO(3).
134C
135C On each step taken by DDASPK, a sequence of nonlinear algebraic
136C systems arises. These are solved by one of two types of
137C methods:
138C * a Newton iteration with a direct method for the linear
139C systems involved (INFO(12) = 0), or
140C * a Newton iteration with a preconditioned Krylov iterative
141C method for the linear systems involved (INFO(12) = 1).
142C
143C The direct method choices are dense and band matrix solvers,
144C with either a user-supplied or an internal difference quotient
145C Jacobian matrix, as specified by INFO(5) and INFO(6).
146C In the band case, INFO(6) = 1, you must supply half-bandwidths
147C in IWORK(1) and IWORK(2).
148C
149C The Krylov method is the Generalized Minimum Residual (GMRES)
150C method, in either complete or incomplete form, and with
151C scaling and preconditioning. The method is implemented
152C in an algorithm called SPIGMR. Certain options in the Krylov
153C method case are specified by INFO(13) and INFO(15).
154C
155C If the Krylov method is chosen, you may supply a pair of routines,
156C JAC and PSOL, to apply preconditioning to the linear system.
157C If the system is A*x = b, the matrix is A = dG/dY + CJ*dG/dYPRIME
158C (of order NEQ). This system can then be preconditioned in the form
159C (P-inverse)*A*x = (P-inverse)*b, with left preconditioner P.
160C (DDASPK does not allow right preconditioning.)
161C Then the Krylov method is applied to this altered, but equivalent,
162C linear system, hopefully with much better performance than without
163C preconditioning. (In addition, a diagonal scaling matrix based on
164C the tolerances is also introduced into the altered system.)
165C
166C The JAC routine evaluates any data needed for solving systems
167C with coefficient matrix P, and PSOL carries out that solution.
168C In any case, in order to improve convergence, you should try to
169C make P approximate the matrix A as much as possible, while keeping
170C the system P*x = b reasonably easy and inexpensive to solve for x,
171C given a vector b.
172C
173C
174C *Description
175C
176C------INPUT - WHAT TO DO ON THE FIRST CALL TO DDASPK-------------------
177C
178C
179C The first call of the code is defined to be the start of each new
180C problem. Read through the descriptions of all the following items,
181C provide sufficient storage space for designated arrays, set
182C appropriate variables for the initialization of the problem, and
183C give information about how you want the problem to be solved.
184C
185C
186C RES -- Provide a subroutine of the form
187C
188C SUBROUTINE RES (T, Y, YPRIME, CJ, DELTA, IRES, RPAR, IPAR)
189C
190C to define the system of differential/algebraic
191C equations which is to be solved. For the given values
192C of T, Y and YPRIME, the subroutine should return
193C the residual of the differential/algebraic system
194C DELTA = G(T,Y,YPRIME)
195C DELTA is a vector of length NEQ which is output from RES.
196C
197C Subroutine RES must not alter T, Y, YPRIME, or CJ.
198C You must declare the name RES in an EXTERNAL
199C statement in your program that calls DDASPK.
200C You must dimension Y, YPRIME, and DELTA in RES.
201C
202C The input argument CJ can be ignored, or used to rescale
203C constraint equations in the system (see Ref. 2, p. 145).
204C Note: In this respect, DDASPK is not downward-compatible
205C with DDASSL, which does not have the RES argument CJ.
206C
207C IRES is an integer flag which is always equal to zero
208C on input. Subroutine RES should alter IRES only if it
209C encounters an illegal value of Y or a stop condition.
210C Set IRES = -1 if an input value is illegal, and DDASPK
211C will try to solve the problem without getting IRES = -1.
212C If IRES = -2, DDASPK will return control to the calling
213C program with IDID = -11.
214C
215C RPAR and IPAR are real and integer parameter arrays which
216C you can use for communication between your calling program
217C and subroutine RES. They are not altered by DDASPK. If you
218C do not need RPAR or IPAR, ignore these parameters by treat-
219C ing them as dummy arguments. If you do choose to use them,
220C dimension them in your calling program and in RES as arrays
221C of appropriate length.
222C
223C NEQ -- Set it to the number of equations in the system (NEQ .GE. 1).
224C
225C T -- Set it to the initial point of the integration. (T must be
226C a variable.)
227C
228C Y(*) -- Set this array to the initial values of the NEQ solution
229C components at the initial point. You must dimension Y of
230C length at least NEQ in your calling program.
231C
232C YPRIME(*) -- Set this array to the initial values of the NEQ first
233C derivatives of the solution components at the initial
234C point. You must dimension YPRIME at least NEQ in your
235C calling program.
236C
237C TOUT - Set it to the first point at which a solution is desired.
238C You cannot take TOUT = T. Integration either forward in T
239C (TOUT .GT. T) or backward in T (TOUT .LT. T) is permitted.
240C
241C The code advances the solution from T to TOUT using step
242C sizes which are automatically selected so as to achieve the
243C desired accuracy. If you wish, the code will return with the
244C solution and its derivative at intermediate steps (the
245C intermediate-output mode) so that you can monitor them,
246C but you still must provide TOUT in accord with the basic
247C aim of the code.
248C
249C The first step taken by the code is a critical one because
250C it must reflect how fast the solution changes near the
251C initial point. The code automatically selects an initial
252C step size which is practically always suitable for the
253C problem. By using the fact that the code will not step past
254C TOUT in the first step, you could, if necessary, restrict the
255C length of the initial step.
256C
257C For some problems it may not be permissible to integrate
258C past a point TSTOP, because a discontinuity occurs there
259C or the solution or its derivative is not defined beyond
260C TSTOP. When you have declared a TSTOP point (see INFO(4)
261C and RWORK(1)), you have told the code not to integrate past
262C TSTOP. In this case any tout beyond TSTOP is invalid input.
263C
264C INFO(*) - Use the INFO array to give the code more details about
265C how you want your problem solved. This array should be
266C dimensioned of length 20, though DDASPK uses only the
267C first 15 entries. You must respond to all of the following
268C items, which are arranged as questions. The simplest use
269C of DDASPK corresponds to setting all entries of INFO to 0.
270C
271C INFO(1) - This parameter enables the code to initialize itself.
272C You must set it to indicate the start of every new
273C problem.
274C
275C **** Is this the first call for this problem ...
276C yes - set INFO(1) = 0
277C no - not applicable here.
278C See below for continuation calls. ****
279C
280C INFO(2) - How much accuracy you want of your solution
281C is specified by the error tolerances RTOL and ATOL.
282C The simplest use is to take them both to be scalars.
283C To obtain more flexibility, they can both be arrays.
284C The code must be told your choice.
285C
286C **** Are both error tolerances RTOL, ATOL scalars ...
287C yes - set INFO(2) = 0
288C and input scalars for both RTOL and ATOL
289C no - set INFO(2) = 1
290C and input arrays for both RTOL and ATOL ****
291C
292C INFO(3) - The code integrates from T in the direction of TOUT
293C by steps. If you wish, it will return the computed
294C solution and derivative at the next intermediate step
295C (the intermediate-output mode) or TOUT, whichever comes
296C first. This is a good way to proceed if you want to
297C see the behavior of the solution. If you must have
298C solutions at a great many specific TOUT points, this
299C code will compute them efficiently.
300C
301C **** Do you want the solution only at
302C TOUT (and not at the next intermediate step) ...
303C yes - set INFO(3) = 0
304C no - set INFO(3) = 1 ****
305C
306C INFO(4) - To handle solutions at a great many specific
307C values TOUT efficiently, this code may integrate past
308C TOUT and interpolate to obtain the result at TOUT.
309C Sometimes it is not possible to integrate beyond some
310C point TSTOP because the equation changes there or it is
311C not defined past TSTOP. Then you must tell the code
312C this stop condition.
313C
314C **** Can the integration be carried out without any
315C restrictions on the independent variable T ...
316C yes - set INFO(4) = 0
317C no - set INFO(4) = 1
318C and define the stopping point TSTOP by
319C setting RWORK(1) = TSTOP ****
320C
321C INFO(5) - used only when INFO(12) = 0 (direct methods).
322C To solve differential/algebraic systems you may wish
323C to use a matrix of partial derivatives of the
324C system of differential equations. If you do not
325C provide a subroutine to evaluate it analytically (see
326C description of the item JAC in the call list), it will
327C be approximated by numerical differencing in this code.
328C Although it is less trouble for you to have the code
329C compute partial derivatives by numerical differencing,
330C the solution will be more reliable if you provide the
331C derivatives via JAC. Usually numerical differencing is
332C more costly than evaluating derivatives in JAC, but
333C sometimes it is not - this depends on your problem.
334C
335C **** Do you want the code to evaluate the partial deriv-
336C atives automatically by numerical differences ...
337C yes - set INFO(5) = 0
338C no - set INFO(5) = 1
339C and provide subroutine JAC for evaluating the
340C matrix of partial derivatives ****
341C
342C INFO(6) - used only when INFO(12) = 0 (direct methods).
343C DDASPK will perform much better if the matrix of
344C partial derivatives, dG/dY + CJ*dG/dYPRIME (here CJ is
345C a scalar determined by DDASPK), is banded and the code
346C is told this. In this case, the storage needed will be
347C greatly reduced, numerical differencing will be performed
348C much cheaper, and a number of important algorithms will
349C execute much faster. The differential equation is said
350C to have half-bandwidths ML (lower) and MU (upper) if
351C equation i involves only unknowns Y(j) with
352C i-ML .le. j .le. i+MU .
353C For all i=1,2,...,NEQ. Thus, ML and MU are the widths
354C of the lower and upper parts of the band, respectively,
355C with the main diagonal being excluded. If you do not
356C indicate that the equation has a banded matrix of partial
357C derivatives the code works with a full matrix of NEQ**2
358C elements (stored in the conventional way). Computations
359C with banded matrices cost less time and storage than with
360C full matrices if 2*ML+MU .lt. NEQ. If you tell the
361C code that the matrix of partial derivatives has a banded
362C structure and you want to provide subroutine JAC to
363C compute the partial derivatives, then you must be careful
364C to store the elements of the matrix in the special form
365C indicated in the description of JAC.
366C
367C **** Do you want to solve the problem using a full (dense)
368C matrix (and not a special banded structure) ...
369C yes - set INFO(6) = 0
370C no - set INFO(6) = 1
371C and provide the lower (ML) and upper (MU)
372C bandwidths by setting
373C IWORK(1)=ML
374C IWORK(2)=MU ****
375C
376C INFO(7) - You can specify a maximum (absolute value of)
377C stepsize, so that the code will avoid passing over very
378C large regions.
379C
380C **** Do you want the code to decide on its own the maximum
381C stepsize ...
382C yes - set INFO(7) = 0
383C no - set INFO(7) = 1
384C and define HMAX by setting
385C RWORK(2) = HMAX ****
386C
387C INFO(8) - Differential/algebraic problems may occasionally
388C suffer from severe scaling difficulties on the first
389C step. If you know a great deal about the scaling of
390C your problem, you can help to alleviate this problem
391C by specifying an initial stepsize H0.
392C
393C **** Do you want the code to define its own initial
394C stepsize ...
395C yes - set INFO(8) = 0
396C no - set INFO(8) = 1
397C and define H0 by setting
398C RWORK(3) = H0 ****
399C
400C INFO(9) - If storage is a severe problem, you can save some
401C storage by restricting the maximum method order MAXORD.
402C The default value is 5. For each order decrease below 5,
403C the code requires NEQ fewer locations, but it is likely
404C to be slower. In any case, you must have
405C 1 .le. MAXORD .le. 5.
406C **** Do you want the maximum order to default to 5 ...
407C yes - set INFO(9) = 0
408C no - set INFO(9) = 1
409C and define MAXORD by setting
410C IWORK(3) = MAXORD ****
411C
412C INFO(10) - If you know that certain components of the
413C solutions to your equations are always nonnegative
414C (or nonpositive), it may help to set this
415C parameter. There are three options that are
416C available:
417C 1. To have constraint checking only in the initial
418C condition calculation.
419C 2. To enforce nonnegativity in Y during the integration.
420C 3. To enforce both options 1 and 2.
421C
422C When selecting option 2 or 3, it is probably best to try the
423C code without using this option first, and only use
424C this option if that does not work very well.
425C
426C **** Do you want the code to solve the problem without
427C invoking any special inequality constraints ...
428C yes - set INFO(10) = 0
429C no - set INFO(10) = 1 to have option 1 enforced
430C no - set INFO(10) = 2 to have option 2 enforced
431C no - set INFO(10) = 3 to have option 3 enforced ****
432C
433C If you have specified INFO(10) = 1 or 3, then you
434C will also need to identify how each component of Y
435C in the initial condition calculation is constrained.
436C You must set:
437C IWORK(40+I) = +1 if Y(I) must be .GE. 0,
438C IWORK(40+I) = +2 if Y(I) must be .GT. 0,
439C IWORK(40+I) = -1 if Y(I) must be .LE. 0, while
440C IWORK(40+I) = -2 if Y(I) must be .LT. 0, while
441C IWORK(40+I) = 0 if Y(I) is not constrained.
442C
443C INFO(11) - DDASPK normally requires the initial T, Y, and
444C YPRIME to be consistent. That is, you must have
445C G(T,Y,YPRIME) = 0 at the initial T. If you do not know
446C the initial conditions precisely, in some cases
447C DDASPK may be able to compute it.
448C
449C Denoting the differential variables in Y by Y_d
450C and the algebraic variables by Y_a, DDASPK can solve
451C one of two initialization problems:
452C 1. Given Y_d, calculate Y_a and Y'_d, or
453C 2. Given Y', calculate Y.
454C In either case, initial values for the given
455C components are input, and initial guesses for
456C the unknown components must also be provided as input.
457C
458C **** Are the initial T, Y, YPRIME consistent ...
459C
460C yes - set INFO(11) = 0
461C no - set INFO(11) = 1 to calculate option 1 above,
462C or set INFO(11) = 2 to calculate option 2 ****
463C
464C If you have specified INFO(11) = 1, then you
465C will also need to identify which are the
466C differential and which are the algebraic
467C components (algebraic components are components
468C whose derivatives do not appear explicitly
469C in the function G(T,Y,YPRIME)). You must set:
470C IWORK(LID+I) = +1 if Y(I) is a differential variable
471C IWORK(LID+I) = -1 if Y(I) is an algebraic variable,
472C where LID = 40 if INFO(10) = 0 or 2 and LID = 40+NEQ
473C if INFO(10) = 1 or 3.
474C
475C INFO(12) - Except for the addition of the RES argument CJ,
476C DDASPK by default is downward-compatible with DDASSL,
477C which uses only direct (dense or band) methods to solve
478C the linear systems involved. You must set INFO(12) to
479C indicate whether you want the direct methods or the
480C Krylov iterative method.
481C **** Do you want DDASPK to use standard direct methods
482C (dense or band) or the Krylov (iterative) method ...
483C direct methods - set INFO(12) = 0.
484C Krylov method - set INFO(12) = 1,
485C and check the settings of INFO(13) and INFO(15).
486C
487C INFO(13) - used when INFO(12) = 1 (Krylov methods).
488C DDASPK uses scalars MAXL, KMP, NRMAX, and EPLI for the
489C iterative solution of linear systems. INFO(13) allows
490C you to override the default values of these parameters.
491C These parameters and their defaults are as follows:
492C MAXL = maximum number of iterations in the SPIGMR
493C algorithm (MAXL .le. NEQ). The default is
494C MAXL = MIN(5,NEQ).
495C KMP = number of vectors on which orthogonalization is
496C done in the SPIGMR algorithm. The default is
497C KMP = MAXL, which corresponds to complete GMRES
498C iteration, as opposed to the incomplete form.
499C NRMAX = maximum number of restarts of the SPIGMR
500C algorithm per nonlinear iteration. The default is
501C NRMAX = 5.
502C EPLI = convergence test constant in SPIGMR algorithm.
503C The default is EPLI = 0.05.
504C Note that the length of RWORK depends on both MAXL
505C and KMP. See the definition of LRW below.
506C **** Are MAXL, KMP, and EPLI to be given their
507C default values ...
508C yes - set INFO(13) = 0
509C no - set INFO(13) = 1,
510C and set all of the following:
511C IWORK(24) = MAXL (1 .le. MAXL .le. NEQ)
512C IWORK(25) = KMP (1 .le. KMP .le. MAXL)
513C IWORK(26) = NRMAX (NRMAX .ge. 0)
514C RWORK(10) = EPLI (0 .lt. EPLI .lt. 1.0) ****
515C
516C INFO(14) - used with INFO(11) > 0 (initial condition
517C calculation is requested). In this case, you may
518C request control to be returned to the calling program
519C immediately after the initial condition calculation,
520C before proceeding to the integration of the system
521C (e.g. to examine the computed Y and YPRIME).
522C If this is done, and if the initialization succeeded
523C (IDID = 4), you should reset INFO(11) to 0 for the
524C next call, to prevent the solver from repeating the
525C initialization (and to avoid an infinite loop).
526C **** Do you want to proceed to the integration after
527C the initial condition calculation is done ...
528C yes - set INFO(14) = 0
529C no - set INFO(14) = 1 ****
530C
531C INFO(15) - used when INFO(12) = 1 (Krylov methods).
532C When using preconditioning in the Krylov method,
533C you must supply a subroutine, PSOL, which solves the
534C associated linear systems using P.
535C The usage of DDASPK is simpler if PSOL can carry out
536C the solution without any prior calculation of data.
537C However, if some partial derivative data is to be
538C calculated in advance and used repeatedly in PSOL,
539C then you must supply a JAC routine to do this,
540C and set INFO(15) to indicate that JAC is to be called
541C for this purpose. For example, P might be an
542C approximation to a part of the matrix A which can be
543C calculated and LU-factored for repeated solutions of
544C the preconditioner system. The arrays WP and IWP
545C (described under JAC and PSOL) can be used to
546C communicate data between JAC and PSOL.
547C **** Does PSOL operate with no prior preparation ...
548C yes - set INFO(15) = 0 (no JAC routine)
549C no - set INFO(15) = 1
550C and supply a JAC routine to evaluate and
551C preprocess any required Jacobian data. ****
552C
553C INFO(16) - option to exclude algebraic variables from
554C the error test.
555C **** Do you wish to control errors locally on
556C all the variables...
557C yes - set INFO(16) = 0
558C no - set INFO(16) = 1
559C If you have specified INFO(16) = 1, then you
560C will also need to identify which are the
561C differential and which are the algebraic
562C components (algebraic components are components
563C whose derivatives do not appear explicitly
564C in the function G(T,Y,YPRIME)). You must set:
565C IWORK(LID+I) = +1 if Y(I) is a differential
566C variable, and
567C IWORK(LID+I) = -1 if Y(I) is an algebraic
568C variable,
569C where LID = 40 if INFO(10) = 0 or 2 and
570C LID = 40 + NEQ if INFO(10) = 1 or 3.
571C
572C INFO(17) - used when INFO(11) > 0 (DDASPK is to do an
573C initial condition calculation).
574C DDASPK uses several heuristic control quantities in the
575C initial condition calculation. They have default values,
576C but can also be set by the user using INFO(17).
577C These parameters and their defaults are as follows:
578C MXNIT = maximum number of Newton iterations
579C per Jacobian or preconditioner evaluation.
580C The default is:
581C MXNIT = 5 in the direct case (INFO(12) = 0), and
582C MXNIT = 15 in the Krylov case (INFO(12) = 1).
583C MXNJ = maximum number of Jacobian or preconditioner
584C evaluations. The default is:
585C MXNJ = 6 in the direct case (INFO(12) = 0), and
586C MXNJ = 2 in the Krylov case (INFO(12) = 1).
587C MXNH = maximum number of values of the artificial
588C stepsize parameter H to be tried if INFO(11) = 1.
589C The default is MXNH = 5.
590C NOTE: the maximum number of Newton iterations
591C allowed in all is MXNIT*MXNJ*MXNH if INFO(11) = 1,
592C and MXNIT*MXNJ if INFO(11) = 2.
593C LSOFF = flag to turn off the linesearch algorithm
594C (LSOFF = 0 means linesearch is on, LSOFF = 1 means
595C it is turned off). The default is LSOFF = 0.
596C STPTOL = minimum scaled step in linesearch algorithm.
597C The default is STPTOL = (unit roundoff)**(2/3).
598C EPINIT = swing factor in the Newton iteration convergence
599C test. The test is applied to the residual vector,
600C premultiplied by the approximate Jacobian (in the
601C direct case) or the preconditioner (in the Krylov
602C case). For convergence, the weighted RMS norm of
603C this vector (scaled by the error weights) must be
604C less than EPINIT*EPCON, where EPCON = .33 is the
605C analogous test constant used in the time steps.
606C The default is EPINIT = .01.
607C **** Are the initial condition heuristic controls to be
608C given their default values...
609C yes - set INFO(17) = 0
610C no - set INFO(17) = 1,
611C and set all of the following:
612C IWORK(32) = MXNIT (.GT. 0)
613C IWORK(33) = MXNJ (.GT. 0)
614C IWORK(34) = MXNH (.GT. 0)
615C IWORK(35) = LSOFF ( = 0 or 1)
616C RWORK(14) = STPTOL (.GT. 0.0)
617C RWORK(15) = EPINIT (.GT. 0.0) ****
618C
619C INFO(18) - option to get extra printing in initial condition
620C calculation.
621C **** Do you wish to have extra printing...
622C no - set INFO(18) = 0
623C yes - set INFO(18) = 1 for minimal printing, or
624C set INFO(18) = 2 for full printing.
625C If you have specified INFO(18) .ge. 1, data
626C will be printed with the error handler routines.
627C To print to a non-default unit number L, include
628C the line CALL XSETUN(L) in your program. ****
629C
630C RTOL, ATOL -- You must assign relative (RTOL) and absolute (ATOL)
631C error tolerances to tell the code how accurately you
632C want the solution to be computed. They must be defined
633C as variables because the code may change them.
634C you have two choices --
635C Both RTOL and ATOL are scalars (INFO(2) = 0), or
636C both RTOL and ATOL are vectors (INFO(2) = 1).
637C In either case all components must be non-negative.
638C
639C The tolerances are used by the code in a local error
640C test at each step which requires roughly that
641C abs(local error in Y(i)) .le. EWT(i) ,
642C where EWT(i) = RTOL*abs(Y(i)) + ATOL is an error weight
643C quantity, for each vector component.
644C (More specifically, a root-mean-square norm is used to
645C measure the size of vectors, and the error test uses the
646C magnitude of the solution at the beginning of the step.)
647C
648C The true (global) error is the difference between the
649C true solution of the initial value problem and the
650C computed approximation. Practically all present day
651C codes, including this one, control the local error at
652C each step and do not even attempt to control the global
653C error directly.
654C
655C Usually, but not always, the true accuracy of
656C the computed Y is comparable to the error tolerances.
657C This code will usually, but not always, deliver a more
658C accurate solution if you reduce the tolerances and
659C integrate again. By comparing two such solutions you
660C can get a fairly reliable idea of the true error in the
661C solution at the larger tolerances.
662C
663C Setting ATOL = 0. results in a pure relative error test
664C on that component. Setting RTOL = 0. results in a pure
665C absolute error test on that component. A mixed test
666C with non-zero RTOL and ATOL corresponds roughly to a
667C relative error test when the solution component is
668C much bigger than ATOL and to an absolute error test
669C when the solution component is smaller than the
670C threshold ATOL.
671C
672C The code will not attempt to compute a solution at an
673C accuracy unreasonable for the machine being used. It
674C will advise you if you ask for too much accuracy and
675C inform you as to the maximum accuracy it believes
676C possible.
677C
678C RWORK(*) -- a real work array, which should be dimensioned in your
679C calling program with a length equal to the value of
680C LRW (or greater).
681C
682C LRW -- Set it to the declared length of the RWORK array. The
683C minimum length depends on the options you have selected,
684C given by a base value plus additional storage as described
685C below.
686C
687C If INFO(12) = 0 (standard direct method), the base value is
688C base = 50 + max(MAXORD+4,7)*NEQ.
689C The default value is MAXORD = 5 (see INFO(9)). With the
690C default MAXORD, base = 50 + 9*NEQ.
691C Additional storage must be added to the base value for
692C any or all of the following options:
693C if INFO(6) = 0 (dense matrix), add NEQ**2
694C if INFO(6) = 1 (banded matrix), then
695C if INFO(5) = 0, add (2*ML+MU+1)*NEQ + 2*(NEQ/(ML+MU+1)+1),
696C if INFO(5) = 1, add (2*ML+MU+1)*NEQ,
697C if INFO(16) = 1, add NEQ.
698C
699C If INFO(12) = 1 (Krylov method), the base value is
700C base = 50 + (MAXORD+5)*NEQ + (MAXL+3+MIN0(1,MAXL-KMP))*NEQ +
701C + (MAXL+3)*MAXL + 1 + LENWP.
702C See PSOL for description of LENWP. The default values are:
703C MAXORD = 5 (see INFO(9)), MAXL = min(5,NEQ) and KMP = MAXL
704C (see INFO(13)).
705C With the default values for MAXORD, MAXL and KMP,
706C base = 91 + 18*NEQ + LENWP.
707C Additional storage must be added to the base value for
708C any or all of the following options:
709C if INFO(16) = 1, add NEQ.
710C
711C
712C IWORK(*) -- an integer work array, which should be dimensioned in
713C your calling program with a length equal to the value
714C of LIW (or greater).
715C
716C LIW -- Set it to the declared length of the IWORK array. The
717C minimum length depends on the options you have selected,
718C given by a base value plus additional storage as described
719C below.
720C
721C If INFO(12) = 0 (standard direct method), the base value is
722C base = 40 + NEQ.
723C IF INFO(10) = 1 or 3, add NEQ to the base value.
724C If INFO(11) = 1 or INFO(16) =1, add NEQ to the base value.
725C
726C If INFO(12) = 1 (Krylov method), the base value is
727C base = 40 + LENIWP.
728C See PSOL for description of LENIWP.
729C IF INFO(10) = 1 or 3, add NEQ to the base value.
730C If INFO(11) = 1 or INFO(16) = 1, add NEQ to the base value.
731C
732C
733C RPAR, IPAR -- These are arrays of double precision and integer type,
734C respectively, which are available for you to use
735C for communication between your program that calls
736C DDASPK and the RES subroutine (and the JAC and PSOL
737C subroutines). They are not altered by DDASPK.
738C If you do not need RPAR or IPAR, ignore these
739C parameters by treating them as dummy arguments.
740C If you do choose to use them, dimension them in
741C your calling program and in RES (and in JAC and PSOL)
742C as arrays of appropriate length.
743C
744C JAC -- This is the name of a routine that you may supply
745C (optionally) that relates to the Jacobian matrix of the
746C nonlinear system that the code must solve at each T step.
747C The role of JAC (and its call sequence) depends on whether
748C a direct (INFO(12) = 0) or Krylov (INFO(12) = 1) method
749C is selected.
750C
751C **** INFO(12) = 0 (direct methods):
752C If you are letting the code generate partial derivatives
753C numerically (INFO(5) = 0), then JAC can be absent
754C (or perhaps a dummy routine to satisfy the loader).
755C Otherwise you must supply a JAC routine to compute
756C the matrix A = dG/dY + CJ*dG/dYPRIME. It must have
757C the form
758C
759C SUBROUTINE JAC (T, Y, YPRIME, PD, CJ, RPAR, IPAR)
760C
761C The JAC routine must dimension Y, YPRIME, and PD (and RPAR
762C and IPAR if used). CJ is a scalar which is input to JAC.
763C For the given values of T, Y, and YPRIME, the JAC routine
764C must evaluate the nonzero elements of the matrix A, and
765C store these values in the array PD. The elements of PD are
766C set to zero before each call to JAC, so that only nonzero
767C elements need to be defined.
768C The way you store the elements into the PD array depends
769C on the structure of the matrix indicated by INFO(6).
770C *** INFO(6) = 0 (full or dense matrix) ***
771C Give PD a first dimension of NEQ. When you evaluate the
772C nonzero partial derivatives of equation i (i.e. of G(i))
773C with respect to component j (of Y and YPRIME), you must
774C store the element in PD according to
775C PD(i,j) = dG(i)/dY(j) + CJ*dG(i)/dYPRIME(j).
776C *** INFO(6) = 1 (banded matrix with half-bandwidths ML, MU
777C as described under INFO(6)) ***
778C Give PD a first dimension of 2*ML+MU+1. When you
779C evaluate the nonzero partial derivatives of equation i
780C (i.e. of G(i)) with respect to component j (of Y and
781C YPRIME), you must store the element in PD according to
782C IROW = i - j + ML + MU + 1
783C PD(IROW,j) = dG(i)/dY(j) + CJ*dG(i)/dYPRIME(j).
784C
785C **** INFO(12) = 1 (Krylov method):
786C If you are not calculating Jacobian data in advance for use
787C in PSOL (INFO(15) = 0), JAC can be absent (or perhaps a
788C dummy routine to satisfy the loader). Otherwise, you may
789C supply a JAC routine to compute and preprocess any parts of
790C of the Jacobian matrix A = dG/dY + CJ*dG/dYPRIME that are
791C involved in the preconditioner matrix P.
792C It is to have the form
793C
794C SUBROUTINE JAC (RES, IRES, NEQ, T, Y, YPRIME, REWT, SAVR,
795C WK, H, CJ, WP, IWP, IER, RPAR, IPAR)
796C
797C The JAC routine must dimension Y, YPRIME, REWT, SAVR, WK,
798C and (if used) WP, IWP, RPAR, and IPAR.
799C The Y, YPRIME, and SAVR arrays contain the current values
800C of Y, YPRIME, and the residual G, respectively.
801C The array WK is work space of length NEQ.
802C H is the step size. CJ is a scalar, input to JAC, that is
803C normally proportional to 1/H. REWT is an array of
804C reciprocal error weights, 1/EWT(i), where EWT(i) is
805C RTOL*abs(Y(i)) + ATOL (unless you supplied routine DDAWTS
806C instead), for use in JAC if needed. For example, if JAC
807C computes difference quotient approximations to partial
808C derivatives, the REWT array may be useful in setting the
809C increments used. The JAC routine should do any
810C factorization operations called for, in preparation for
811C solving linear systems in PSOL. The matrix P should
812C be an approximation to the Jacobian,
813C A = dG/dY + CJ*dG/dYPRIME.
814C
815C WP and IWP are real and integer work arrays which you may
816C use for communication between your JAC routine and your
817C PSOL routine. These may be used to store elements of the
818C preconditioner P, or related matrix data (such as factored
819C forms). They are not altered by DDASPK.
820C If you do not need WP or IWP, ignore these parameters by
821C treating them as dummy arguments. If you do use them,
822C dimension them appropriately in your JAC and PSOL routines.
823C See the PSOL description for instructions on setting
824C the lengths of WP and IWP.
825C
826C On return, JAC should set the error flag IER as follows..
827C IER = 0 if JAC was successful,
828C IER .ne. 0 if JAC was unsuccessful (e.g. if Y or YPRIME
829C was illegal, or a singular matrix is found).
830C (If IER .ne. 0, a smaller stepsize will be tried.)
831C IER = 0 on entry to JAC, so need be reset only on a failure.
832C If RES is used within JAC, then a nonzero value of IRES will
833C override any nonzero value of IER (see the RES description).
834C
835C Regardless of the method type, subroutine JAC must not
836C alter T, Y(*), YPRIME(*), H, CJ, or REWT(*).
837C You must declare the name JAC in an EXTERNAL statement in
838C your program that calls DDASPK.
839C
840C PSOL -- This is the name of a routine you must supply if you have
841C selected a Krylov method (INFO(12) = 1) with preconditioning.
842C In the direct case (INFO(12) = 0), PSOL can be absent
843C (a dummy routine may have to be supplied to satisfy the
844C loader). Otherwise, you must provide a PSOL routine to
845C solve linear systems arising from preconditioning.
846C When supplied with INFO(12) = 1, the PSOL routine is to
847C have the form
848C
849C SUBROUTINE PSOL (NEQ, T, Y, YPRIME, SAVR, WK, CJ, WGHT,
850C WP, IWP, B, EPLIN, IER, RPAR, IPAR)
851C
852C The PSOL routine must solve linear systems of the form
853C P*x = b where P is the left preconditioner matrix.
854C
855C The right-hand side vector b is in the B array on input, and
856C PSOL must return the solution vector x in B.
857C The Y, YPRIME, and SAVR arrays contain the current values
858C of Y, YPRIME, and the residual G, respectively.
859C
860C Work space required by JAC and/or PSOL, and space for data to
861C be communicated from JAC to PSOL is made available in the form
862C of arrays WP and IWP, which are parts of the RWORK and IWORK
863C arrays, respectively. The lengths of these real and integer
864C work spaces WP and IWP must be supplied in LENWP and LENIWP,
865C respectively, as follows..
866C IWORK(27) = LENWP = length of real work space WP
867C IWORK(28) = LENIWP = length of integer work space IWP.
868C
869C WK is a work array of length NEQ for use by PSOL.
870C CJ is a scalar, input to PSOL, that is normally proportional
871C to 1/H (H = stepsize). If the old value of CJ
872C (at the time of the last JAC call) is needed, it must have
873C been saved by JAC in WP.
874C
875C WGHT is an array of weights, to be used if PSOL uses an
876C iterative method and performs a convergence test. (In terms
877C of the argument REWT to JAC, WGHT is REWT/sqrt(NEQ).)
878C If PSOL uses an iterative method, it should use EPLIN
879C (a heuristic parameter) as the bound on the weighted norm of
880C the residual for the computed solution. Specifically, the
881C residual vector R should satisfy
882C SQRT (SUM ( (R(i)*WGHT(i))**2 ) ) .le. EPLIN
883C
884C PSOL must not alter NEQ, T, Y, YPRIME, SAVR, CJ, WGHT, EPLIN.
885C
886C On return, PSOL should set the error flag IER as follows..
887C IER = 0 if PSOL was successful,
888C IER .lt. 0 if an unrecoverable error occurred, meaning
889C control will be passed to the calling routine,
890C IER .gt. 0 if a recoverable error occurred, meaning that
891C the step will be retried with the same step size
892C but with a call to JAC to update necessary data,
893C unless the Jacobian data is current, in which case
894C the step will be retried with a smaller step size.
895C IER = 0 on entry to PSOL so need be reset only on a failure.
896C
897C You must declare the name PSOL in an EXTERNAL statement in
898C your program that calls DDASPK.
899C
900C
901C OPTIONALLY REPLACEABLE SUBROUTINE:
902C
903C DDASPK uses a weighted root-mean-square norm to measure the
904C size of various error vectors. The weights used in this norm
905C are set in the following subroutine:
906C
907C SUBROUTINE DDAWTS (NEQ, IWT, RTOL, ATOL, Y, EWT, RPAR, IPAR)
908C DIMENSION RTOL(*), ATOL(*), Y(*), EWT(*), RPAR(*), IPAR(*)
909C
910C A DDAWTS routine has been included with DDASPK which sets the
911C weights according to
912C EWT(I) = RTOL*ABS(Y(I)) + ATOL
913C in the case of scalar tolerances (IWT = 0) or
914C EWT(I) = RTOL(I)*ABS(Y(I)) + ATOL(I)
915C in the case of array tolerances (IWT = 1). (IWT is INFO(2).)
916C In some special cases, it may be appropriate for you to define
917C your own error weights by writing a subroutine DDAWTS to be
918C called instead of the version supplied. However, this should
919C be attempted only after careful thought and consideration.
920C If you supply this routine, you may use the tolerances and Y
921C as appropriate, but do not overwrite these variables. You
922C may also use RPAR and IPAR to communicate data as appropriate.
923C ***Note: Aside from the values of the weights, the choice of
924C norm used in DDASPK (weighted root-mean-square) is not subject
925C to replacement by the user. In this respect, DDASPK is not
926C downward-compatible with the original DDASSL solver (in which
927C the norm routine was optionally user-replaceable).
928C
929C
930C------OUTPUT - AFTER ANY RETURN FROM DDASPK----------------------------
931C
932C The principal aim of the code is to return a computed solution at
933C T = TOUT, although it is also possible to obtain intermediate
934C results along the way. To find out whether the code achieved its
935C goal or if the integration process was interrupted before the task
936C was completed, you must check the IDID parameter.
937C
938C
939C T -- The output value of T is the point to which the solution
940C was successfully advanced.
941C
942C Y(*) -- contains the computed solution approximation at T.
943C
944C YPRIME(*) -- contains the computed derivative approximation at T.
945C
946C IDID -- reports what the code did, described as follows:
947C
948C *** TASK COMPLETED ***
949C Reported by positive values of IDID
950C
951C IDID = 1 -- a step was successfully taken in the
952C intermediate-output mode. The code has not
953C yet reached TOUT.
954C
955C IDID = 2 -- the integration to TSTOP was successfully
956C completed (T = TSTOP) by stepping exactly to TSTOP.
957C
958C IDID = 3 -- the integration to TOUT was successfully
959C completed (T = TOUT) by stepping past TOUT.
960C Y(*) and YPRIME(*) are obtained by interpolation.
961C
962C IDID = 4 -- the initial condition calculation, with
963C INFO(11) > 0, was successful, and INFO(14) = 1.
964C No integration steps were taken, and the solution
965C is not considered to have been started.
966C
967C *** TASK INTERRUPTED ***
968C Reported by negative values of IDID
969C
970C IDID = -1 -- a large amount of work has been expended
971C (about 500 steps).
972C
973C IDID = -2 -- the error tolerances are too stringent.
974C
975C IDID = -3 -- the local error test cannot be satisfied
976C because you specified a zero component in ATOL
977C and the corresponding computed solution component
978C is zero. Thus, a pure relative error test is
979C impossible for this component.
980C
981C IDID = -5 -- there were repeated failures in the evaluation
982C or processing of the preconditioner (in JAC).
983C
984C IDID = -6 -- DDASPK had repeated error test failures on the
985C last attempted step.
986C
987C IDID = -7 -- the nonlinear system solver in the time integration
988C could not converge.
989C
990C IDID = -8 -- the matrix of partial derivatives appears
991C to be singular (direct method).
992C
993C IDID = -9 -- the nonlinear system solver in the time integration
994C failed to achieve convergence, and there were repeated
995C error test failures in this step.
996C
997C IDID =-10 -- the nonlinear system solver in the time integration
998C failed to achieve convergence because IRES was equal
999C to -1.
1000C
1001C IDID =-11 -- IRES = -2 was encountered and control is
1002C being returned to the calling program.
1003C
1004C IDID =-12 -- DDASPK failed to compute the initial Y, YPRIME.
1005C
1006C IDID =-13 -- unrecoverable error encountered inside user's
1007C PSOL routine, and control is being returned to
1008C the calling program.
1009C
1010C IDID =-14 -- the Krylov linear system solver could not
1011C achieve convergence.
1012C
1013C IDID =-15,..,-32 -- Not applicable for this code.
1014C
1015C *** TASK TERMINATED ***
1016C reported by the value of IDID=-33
1017C
1018C IDID = -33 -- the code has encountered trouble from which
1019C it cannot recover. A message is printed
1020C explaining the trouble and control is returned
1021C to the calling program. For example, this occurs
1022C when invalid input is detected.
1023C
1024C RTOL, ATOL -- these quantities remain unchanged except when
1025C IDID = -2. In this case, the error tolerances have been
1026C increased by the code to values which are estimated to
1027C be appropriate for continuing the integration. However,
1028C the reported solution at T was obtained using the input
1029C values of RTOL and ATOL.
1030C
1031C RWORK, IWORK -- contain information which is usually of no interest
1032C to the user but necessary for subsequent calls.
1033C However, you may be interested in the performance data
1034C listed below. These quantities are accessed in RWORK
1035C and IWORK but have internal mnemonic names, as follows..
1036C
1037C RWORK(3)--contains H, the step size h to be attempted
1038C on the next step.
1039C
1040C RWORK(4)--contains TN, the current value of the
1041C independent variable, i.e. the farthest point
1042C integration has reached. This will differ
1043C from T if interpolation has been performed
1044C (IDID = 3).
1045C
1046C RWORK(7)--contains HOLD, the stepsize used on the last
1047C successful step. If INFO(11) = INFO(14) = 1,
1048C this contains the value of H used in the
1049C initial condition calculation.
1050C
1051C IWORK(7)--contains K, the order of the method to be
1052C attempted on the next step.
1053C
1054C IWORK(8)--contains KOLD, the order of the method used
1055C on the last step.
1056C
1057C IWORK(11)--contains NST, the number of steps (in T)
1058C taken so far.
1059C
1060C IWORK(12)--contains NRE, the number of calls to RES
1061C so far.
1062C
1063C IWORK(13)--contains NJE, the number of calls to JAC so
1064C far (Jacobian or preconditioner evaluations).
1065C
1066C IWORK(14)--contains NETF, the total number of error test
1067C failures so far.
1068C
1069C IWORK(15)--contains NCFN, the total number of nonlinear
1070C convergence failures so far (includes counts
1071C of singular iteration matrix or singular
1072C preconditioners).
1073C
1074C IWORK(16)--contains NCFL, the number of convergence
1075C failures of the linear iteration so far.
1076C
1077C IWORK(17)--contains LENIW, the length of IWORK actually
1078C required. This is defined on normal returns
1079C and on an illegal input return for
1080C insufficient storage.
1081C
1082C IWORK(18)--contains LENRW, the length of RWORK actually
1083C required. This is defined on normal returns
1084C and on an illegal input return for
1085C insufficient storage.
1086C
1087C IWORK(19)--contains NNI, the total number of nonlinear
1088C iterations so far (each of which calls a
1089C linear solver).
1090C
1091C IWORK(20)--contains NLI, the total number of linear
1092C (Krylov) iterations so far.
1093C
1094C IWORK(21)--contains NPS, the number of PSOL calls so
1095C far, for preconditioning solve operations or
1096C for solutions with the user-supplied method.
1097C
1098C Note: The various counters in IWORK do not include
1099C counts during a call made with INFO(11) > 0 and
1100C INFO(14) = 1.
1101C
1102C
1103C------INPUT - WHAT TO DO TO CONTINUE THE INTEGRATION -----------------
1104C (CALLS AFTER THE FIRST)
1105C
1106C This code is organized so that subsequent calls to continue the
1107C integration involve little (if any) additional effort on your
1108C part. You must monitor the IDID parameter in order to determine
1109C what to do next.
1110C
1111C Recalling that the principal task of the code is to integrate
1112C from T to TOUT (the interval mode), usually all you will need
1113C to do is specify a new TOUT upon reaching the current TOUT.
1114C
1115C Do not alter any quantity not specifically permitted below. In
1116C particular do not alter NEQ, T, Y(*), YPRIME(*), RWORK(*),
1117C IWORK(*), or the differential equation in subroutine RES. Any
1118C such alteration constitutes a new problem and must be treated
1119C as such, i.e. you must start afresh.
1120C
1121C You cannot change from array to scalar error control or vice
1122C versa (INFO(2)), but you can change the size of the entries of
1123C RTOL or ATOL. Increasing a tolerance makes the equation easier
1124C to integrate. Decreasing a tolerance will make the equation
1125C harder to integrate and should generally be avoided.
1126C
1127C You can switch from the intermediate-output mode to the
1128C interval mode (INFO(3)) or vice versa at any time.
1129C
1130C If it has been necessary to prevent the integration from going
1131C past a point TSTOP (INFO(4), RWORK(1)), keep in mind that the
1132C code will not integrate to any TOUT beyond the currently
1133C specified TSTOP. Once TSTOP has been reached, you must change
1134C the value of TSTOP or set INFO(4) = 0. You may change INFO(4)
1135C or TSTOP at any time but you must supply the value of TSTOP in
1136C RWORK(1) whenever you set INFO(4) = 1.
1137C
1138C Do not change INFO(5), INFO(6), INFO(12-17) or their associated
1139C IWORK/RWORK locations unless you are going to restart the code.
1140C
1141C *** FOLLOWING A COMPLETED TASK ***
1142C
1143C If..
1144C IDID = 1, call the code again to continue the integration
1145C another step in the direction of TOUT.
1146C
1147C IDID = 2 or 3, define a new TOUT and call the code again.
1148C TOUT must be different from T. You cannot change
1149C the direction of integration without restarting.
1150C
1151C IDID = 4, reset INFO(11) = 0 and call the code again to begin
1152C the integration. (If you leave INFO(11) > 0 and
1153C INFO(14) = 1, you may generate an infinite loop.)
1154C In this situation, the next call to DASPK is
1155C considered to be the first call for the problem,
1156C in that all initializations are done.
1157C
1158C *** FOLLOWING AN INTERRUPTED TASK ***
1159C
1160C To show the code that you realize the task was interrupted and
1161C that you want to continue, you must take appropriate action and
1162C set INFO(1) = 1.
1163C
1164C If..
1165C IDID = -1, the code has taken about 500 steps. If you want to
1166C continue, set INFO(1) = 1 and call the code again.
1167C An additional 500 steps will be allowed.
1168C
1169C
1170C IDID = -2, the error tolerances RTOL, ATOL have been increased
1171C to values the code estimates appropriate for
1172C continuing. You may want to change them yourself.
1173C If you are sure you want to continue with relaxed
1174C error tolerances, set INFO(1) = 1 and call the code
1175C again.
1176C
1177C IDID = -3, a solution component is zero and you set the
1178C corresponding component of ATOL to zero. If you
1179C are sure you want to continue, you must first alter
1180C the error criterion to use positive values of ATOL
1181C for those components corresponding to zero solution
1182C components, then set INFO(1) = 1 and call the code
1183C again.
1184C
1185C IDID = -4 --- cannot occur with this code.
1186C
1187C IDID = -5, your JAC routine failed with the Krylov method. Check
1188C for errors in JAC and restart the integration.
1189C
1190C IDID = -6, repeated error test failures occurred on the last
1191C attempted step in DDASPK. A singularity in the
1192C solution may be present. If you are absolutely
1193C certain you want to continue, you should restart
1194C the integration. (Provide initial values of Y and
1195C YPRIME which are consistent.)
1196C
1197C IDID = -7, repeated convergence test failures occurred on the last
1198C attempted step in DDASPK. An inaccurate or ill-
1199C conditioned Jacobian or preconditioner may be the
1200C problem. If you are absolutely certain you want
1201C to continue, you should restart the integration.
1202C
1203C
1204C IDID = -8, the matrix of partial derivatives is singular, with
1205C the use of direct methods. Some of your equations
1206C may be redundant. DDASPK cannot solve the problem
1207C as stated. It is possible that the redundant
1208C equations could be removed, and then DDASPK could
1209C solve the problem. It is also possible that a
1210C solution to your problem either does not exist
1211C or is not unique.
1212C
1213C IDID = -9, DDASPK had multiple convergence test failures, preceded
1214C by multiple error test failures, on the last
1215C attempted step. It is possible that your problem is
1216C ill-posed and cannot be solved using this code. Or,
1217C there may be a discontinuity or a singularity in the
1218C solution. If you are absolutely certain you want to
1219C continue, you should restart the integration.
1220C
1221C IDID = -10, DDASPK had multiple convergence test failures
1222C because IRES was equal to -1. If you are
1223C absolutely certain you want to continue, you
1224C should restart the integration.
1225C
1226C IDID = -11, there was an unrecoverable error (IRES = -2) from RES
1227C inside the nonlinear system solver. Determine the
1228C cause before trying again.
1229C
1230C IDID = -12, DDASPK failed to compute the initial Y and YPRIME
1231C vectors. This could happen because the initial
1232C approximation to Y or YPRIME was not very good, or
1233C because no consistent values of these vectors exist.
1234C The problem could also be caused by an inaccurate or
1235C singular iteration matrix, or a poor preconditioner.
1236C
1237C IDID = -13, there was an unrecoverable error encountered inside
1238C your PSOL routine. Determine the cause before
1239C trying again.
1240C
1241C IDID = -14, the Krylov linear system solver failed to achieve
1242C convergence. This may be due to ill-conditioning
1243C in the iteration matrix, or a singularity in the
1244C preconditioner (if one is being used).
1245C Another possibility is that there is a better
1246C choice of Krylov parameters (see INFO(13)).
1247C Possibly the failure is caused by redundant equations
1248C in the system, or by inconsistent equations.
1249C In that case, reformulate the system to make it
1250C consistent and non-redundant.
1251C
1252C IDID = -15,..,-32 --- Cannot occur with this code.
1253C
1254C *** FOLLOWING A TERMINATED TASK ***
1255C
1256C If IDID = -33, you cannot continue the solution of this problem.
1257C An attempt to do so will result in your run being
1258C terminated.
1259C
1260C ---------------------------------------------------------------------
1261C
1262C***REFERENCES
1263C 1. L. R. Petzold, A Description of DASSL: A Differential/Algebraic
1264C System Solver, in Scientific Computing, R. S. Stepleman et al.
1265C (Eds.), North-Holland, Amsterdam, 1983, pp. 65-68.
1266C 2. K. E. Brenan, S. L. Campbell, and L. R. Petzold, Numerical
1267C Solution of Initial-Value Problems in Differential-Algebraic
1268C Equations, Elsevier, New York, 1989.
1269C 3. P. N. Brown and A. C. Hindmarsh, Reduced Storage Matrix Methods
1270C in Stiff ODE Systems, J. Applied Mathematics and Computation,
1271C 31 (1989), pp. 40-91.
1272C 4. P. N. Brown, A. C. Hindmarsh, and L. R. Petzold, Using Krylov
1273C Methods in the Solution of Large-Scale Differential-Algebraic
1274C Systems, SIAM J. Sci. Comp., 15 (1994), pp. 1467-1488.
1275C 5. P. N. Brown, A. C. Hindmarsh, and L. R. Petzold, Consistent
1276C Initial Condition Calculation for Differential-Algebraic
1277C Systems, LLNL Report UCRL-JC-122175, August 1995; submitted to
1278C SIAM J. Sci. Comp.
1279C
1280C***ROUTINES CALLED
1281C
1282C The following are all the subordinate routines used by DDASPK.
1283C
1284C DDASIC computes consistent initial conditions.
1285C DYYPNW updates Y and YPRIME in linesearch for initial condition
1286C calculation.
1287C DDSTP carries out one step of the integration.
1288C DCNSTR/DCNST0 check the current solution for constraint violations.
1289C DDAWTS sets error weight quantities.
1290C DINVWT tests and inverts the error weights.
1291C DDATRP performs interpolation to get an output solution.
1292C DDWNRM computes the weighted root-mean-square norm of a vector.
1293C D1MACH provides the unit roundoff of the computer.
1294C XERRWD/XSETF/XSETUN/IXSAV is a package to handle error messages.
1295C DDASID nonlinear equation driver to initialize Y and YPRIME using
1296C direct linear system solver methods. Interfaces to Newton
1297C solver (direct case).
1298C DNSID solves the nonlinear system for unknown initial values by
1299C modified Newton iteration and direct linear system methods.
1300C DLINSD carries out linesearch algorithm for initial condition
1301C calculation (direct case).
1302C DFNRMD calculates weighted norm of preconditioned residual in
1303C initial condition calculation (direct case).
1304C DNEDD nonlinear equation driver for direct linear system solver
1305C methods. Interfaces to Newton solver (direct case).
1306C DMATD assembles the iteration matrix (direct case).
1307C DNSD solves the associated nonlinear system by modified
1308C Newton iteration and direct linear system methods.
1309C DSLVD interfaces to linear system solver (direct case).
1310C DDASIK nonlinear equation driver to initialize Y and YPRIME using
1311C Krylov iterative linear system methods. Interfaces to
1312C Newton solver (Krylov case).
1313C DNSIK solves the nonlinear system for unknown initial values by
1314C Newton iteration and Krylov iterative linear system methods.
1315C DLINSK carries out linesearch algorithm for initial condition
1316C calculation (Krylov case).
1317C DFNRMK calculates weighted norm of preconditioned residual in
1318C initial condition calculation (Krylov case).
1319C DNEDK nonlinear equation driver for iterative linear system solver
1320C methods. Interfaces to Newton solver (Krylov case).
1321C DNSK solves the associated nonlinear system by Inexact Newton
1322C iteration and (linear) Krylov iteration.
1323C DSLVK interfaces to linear system solver (Krylov case).
1324C DSPIGM solves a linear system by SPIGMR algorithm.
1325C DATV computes matrix-vector product in Krylov algorithm.
1326C DORTH performs orthogonalization of Krylov basis vectors.
1327C DHEQR performs QR factorization of Hessenberg matrix.
1328C DHELS finds least-squares solution of Hessenberg linear system.
1329C DGETRF, DGETRS, DGBTRF, DGBTRS are LAPACK routines for solving
1330C linear systems (dense or band direct methods).
1331C DAXPY, DCOPY, DDOT, DNRM2, DSCAL are Basic Linear Algebra (BLAS)
1332C routines.
1333C
1334C The routines called directly by DDASPK are:
1335C DCNST0, DDAWTS, DINVWT, D1MACH, DDWNRM, DDASIC, DDATRP, DDSTP,
1336C XERRWD
1337C
1338C***END PROLOGUE DDASPK
1339C
1340C
1341 IMPLICIT DOUBLE PRECISION(a-h,o-z)
1342 LOGICAL DONE, LAVL, LCFN, LCFL, LWARN
1343 dimension y(*),yprime(*)
1344 dimension info(20)
1345 dimension rwork(lrw),iwork(liw)
1346 dimension rtol(*),atol(*)
1347 dimension rpar(*),ipar(*)
1348 CHARACTER MSG*80
1349 EXTERNAL res, jac, psol, ddasid, ddasik, dnedd, dnedk
1350C
1351C Set pointers into IWORK.
1352C
1353 parameter(lml=1, lmu=2, lmtype=4,
1354 * liwm=1, lmxord=3, ljcalc=5, lphase=6, lk=7, lkold=8,
1355 * lns=9, lnstl=10, lnst=11, lnre=12, lnje=13, letf=14, lncfn=15,
1356 * lncfl=16, lniw=17, lnrw=18, lnni=19, lnli=20, lnps=21,
1357 * lnpd=22, lmiter=23, lmaxl=24, lkmp=25, lnrmax=26, llnwp=27,
1358 * llniwp=28, llocwp=29, llciwp=30, lkprin=31,
1359 * lmxnit=32, lmxnj=33, lmxnh=34, llsoff=35, licns=41)
1360C
1361C Set pointers into RWORK.
1362C
1363 parameter(ltstop=1, lhmax=2, lh=3, ltn=4, lcj=5, lcjold=6,
1364 * lhold=7, ls=8, lround=9, lepli=10, lsqrn=11, lrsqrn=12,
1365 * lepcon=13, lstol=14, lepin=15,
1366 * lalpha=21, lbeta=27, lgamma=33, lpsi=39, lsigma=45, ldelta=51)
1367C
1368 SAVE lid, lenid, nonneg
1369C
1370C
1371C***FIRST EXECUTABLE STATEMENT DDASPK
1372C
1373C
1374 IF(info(1).NE.0) GO TO 100
1375C
1376C-----------------------------------------------------------------------
1377C This block is executed for the initial call only.
1378C It contains checking of inputs and initializations.
1379C-----------------------------------------------------------------------
1380C
1381C First check INFO array to make sure all elements of INFO
1382C Are within the proper range. (INFO(1) is checked later, because
1383C it must be tested on every call.) ITEMP holds the location
1384C within INFO which may be out of range.
1385C
1386 DO 10 i=2,9
1387 itemp = i
1388 IF (info(i) .NE. 0 .AND. info(i) .NE. 1) GO TO 701
1389 10 CONTINUE
1390 itemp = 10
1391 IF(info(10).LT.0 .OR. info(10).GT.3) GO TO 701
1392 itemp = 11
1393 IF(info(11).LT.0 .OR. info(11).GT.2) GO TO 701
1394 DO 15 i=12,17
1395 itemp = i
1396 IF (info(i) .NE. 0 .AND. info(i) .NE. 1) GO TO 701
1397 15 CONTINUE
1398 itemp = 18
1399 IF(info(18).LT.0 .OR. info(18).GT.2) GO TO 701
1400
1401C
1402C Check NEQ to see if it is positive.
1403C
1404 IF (neq .LE. 0) GO TO 702
1405C
1406C Check and compute maximum order.
1407C
1408 mxord=5
1409 IF (info(9) .NE. 0) THEN
1410 mxord=iwork(lmxord)
1411 IF (mxord .LT. 1 .OR. mxord .GT. 5) GO TO 703
1412 ENDIF
1413 iwork(lmxord)=mxord
1414C
1415C Set and/or check inputs for constraint checking (INFO(10) .NE. 0).
1416C Set values for ICNFLG, NONNEG, and pointer LID.
1417C
1418 icnflg = 0
1419 nonneg = 0
1420 lid = licns
1421 IF (info(10) .EQ. 0) GO TO 20
1422 IF (info(10) .EQ. 1) THEN
1423 icnflg = 1
1424 nonneg = 0
1425 lid = licns + neq
1426 ELSEIF (info(10) .EQ. 2) THEN
1427 icnflg = 0
1428 nonneg = 1
1429 ELSE
1430 icnflg = 1
1431 nonneg = 1
1432 lid = licns + neq
1433 ENDIF
1434C
1435 20 CONTINUE
1436C
1437C Set and/or check inputs for Krylov solver (INFO(12) .NE. 0).
1438C If indicated, set default values for MAXL, KMP, NRMAX, and EPLI.
1439C Otherwise, verify inputs required for iterative solver.
1440C
1441 IF (info(12) .EQ. 0) GO TO 25
1442C
1443 iwork(lmiter) = info(12)
1444 IF (info(13) .EQ. 0) THEN
1445 iwork(lmaxl) = min(5,neq)
1446 iwork(lkmp) = iwork(lmaxl)
1447 iwork(lnrmax) = 5
1448 rwork(lepli) = 0.05d0
1449 ELSE
1450 IF(iwork(lmaxl) .LT. 1 .OR. iwork(lmaxl) .GT. neq) GO TO 720
1451 IF(iwork(lkmp) .LT. 1 .OR. iwork(lkmp) .GT. iwork(lmaxl))
1452 1 GO TO 721
1453 IF(iwork(lnrmax) .LT. 0) GO TO 722
1454 IF(rwork(lepli).LE.0.0d0 .OR. rwork(lepli).GE.1.0d0)GO TO 723
1455 ENDIF
1456C
1457 25 CONTINUE
1458C
1459C Set and/or check controls for the initial condition calculation
1460C (INFO(11) .GT. 0). If indicated, set default values.
1461C Otherwise, verify inputs required for iterative solver.
1462C
1463 IF (info(11) .EQ. 0) GO TO 30
1464 IF (info(17) .EQ. 0) THEN
1465 iwork(lmxnit) = 5
1466 IF (info(12) .GT. 0) iwork(lmxnit) = 15
1467 iwork(lmxnj) = 6
1468 IF (info(12) .GT. 0) iwork(lmxnj) = 2
1469 iwork(lmxnh) = 5
1470 iwork(llsoff) = 0
1471 rwork(lepin) = 0.01d0
1472 ELSE
1473 IF (iwork(lmxnit) .LE. 0) GO TO 725
1474 IF (iwork(lmxnj) .LE. 0) GO TO 725
1475 IF (iwork(lmxnh) .LE. 0) GO TO 725
1476 lsoff = iwork(llsoff)
1477 IF (lsoff .LT. 0 .OR. lsoff .GT. 1) GO TO 725
1478 IF (rwork(lepin) .LE. 0.0d0) GO TO 725
1479 ENDIF
1480C
1481 30 CONTINUE
1482C
1483C Below is the computation and checking of the work array lengths
1484C LENIW and LENRW, using direct methods (INFO(12) = 0) or
1485C the Krylov methods (INFO(12) = 1).
1486C
1487 lenic = 0
1488 IF (info(10) .EQ. 1 .OR. info(10) .EQ. 3) lenic = neq
1489 lenid = 0
1490 IF (info(11) .EQ. 1 .OR. info(16) .EQ. 1) lenid = neq
1491 IF (info(12) .EQ. 0) THEN
1492C
1493C Compute MTYPE, etc. Check ML and MU.
1494C
1495 ncphi = max(mxord + 1, 4)
1496 IF(info(6).EQ.0) THEN
1497 lenpd = neq**2
1498 lenrw = 50 + (ncphi+3)*neq + lenpd
1499 IF(info(5).EQ.0) THEN
1500 iwork(lmtype)=2
1501 ELSE
1502 iwork(lmtype)=1
1503 ENDIF
1504 ELSE
1505 IF(iwork(lml).LT.0.OR.iwork(lml).GE.neq)GO TO 717
1506 IF(iwork(lmu).LT.0.OR.iwork(lmu).GE.neq)GO TO 718
1507 lenpd=(2*iwork(lml)+iwork(lmu)+1)*neq
1508 IF(info(5).EQ.0) THEN
1509 iwork(lmtype)=5
1510 mband=iwork(lml)+iwork(lmu)+1
1511 msave=(neq/mband)+1
1512 lenrw = 50 + (ncphi+3)*neq + lenpd + 2*msave
1513 ELSE
1514 iwork(lmtype)=4
1515 lenrw = 50 + (ncphi+3)*neq + lenpd
1516 ENDIF
1517 ENDIF
1518C
1519C Compute LENIW, LENWP, LENIWP.
1520C
1521 leniw = 40 + lenic + lenid + neq
1522 lenwp = 0
1523 leniwp = 0
1524C
1525 ELSE IF (info(12) .EQ. 1) THEN
1526 maxl = iwork(lmaxl)
1527 lenwp = iwork(llnwp)
1528 leniwp = iwork(llniwp)
1529 lenpd = (maxl+3+min0(1,maxl-iwork(lkmp)))*neq
1530 1 + (maxl+3)*maxl + 1 + lenwp
1531 lenrw = 50 + (iwork(lmxord)+5)*neq + lenpd
1532 leniw = 40 + lenic + lenid + leniwp
1533C
1534 ENDIF
1535 IF(info(16) .NE. 0) lenrw = lenrw + neq
1536C
1537C Check lengths of RWORK and IWORK.
1538C
1539 iwork(lniw)=leniw
1540 iwork(lnrw)=lenrw
1541 iwork(lnpd)=lenpd
1542 iwork(llocwp) = lenpd-lenwp+1
1543 IF(lrw.LT.lenrw)GO TO 704
1544 IF(liw.LT.leniw)GO TO 705
1545C
1546C Check ICNSTR for legality.
1547C
1548 IF (lenic .GT. 0) THEN
1549 DO 40 i = 1,neq
1550 ici = iwork(licns-1+i)
1551 IF (ici .LT. -2 .OR. ici .GT. 2) GO TO 726
1552 40 CONTINUE
1553 ENDIF
1554C
1555C Check Y for consistency with constraints.
1556C
1557 IF (lenic .GT. 0) THEN
1558 CALL dcnst0(neq,y,iwork(licns),iret)
1559 IF (iret .NE. 0) GO TO 727
1560 ENDIF
1561C
1562C Check ID for legality.
1563C
1564 IF (lenid .GT. 0) THEN
1565 DO 50 i = 1,neq
1566 idi = iwork(lid-1+i)
1567 IF (idi .NE. 1 .AND. idi .NE. -1) GO TO 724
1568 50 CONTINUE
1569 ENDIF
1570C
1571C Check to see that TOUT is different from T.
1572C
1573 IF(tout .EQ. t)GO TO 719
1574C
1575C Check HMAX.
1576C
1577 IF(info(7) .NE. 0) THEN
1578 hmax = rwork(lhmax)
1579 IF (hmax .LE. 0.0d0) GO TO 710
1580 ENDIF
1581C
1582C Initialize counters and other flags.
1583C
1584 iwork(lnst)=0
1585 iwork(lnre)=0
1586 iwork(lnje)=0
1587 iwork(letf)=0
1588 iwork(lncfn)=0
1589 iwork(lnni)=0
1590 iwork(lnli)=0
1591 iwork(lnps)=0
1592 iwork(lncfl)=0
1593 iwork(lkprin)=info(18)
1594 idid=1
1595 GO TO 200
1596C
1597C-----------------------------------------------------------------------
1598C This block is for continuation calls only.
1599C Here we check INFO(1), and if the last step was interrupted,
1600C we check whether appropriate action was taken.
1601C-----------------------------------------------------------------------
1602C
1603100 CONTINUE
1604 IF(info(1).EQ.1)GO TO 110
1605 itemp = 1
1606 IF(info(1).NE.-1)GO TO 701
1607C
1608C If we are here, the last step was interrupted by an error
1609C condition from DDSTP, and appropriate action was not taken.
1610C This is a fatal error.
1611C
1612 msg = 'DASPK-- THE LAST STEP TERMINATED WITH A NEGATIVE'
1613 CALL xerrwd(msg,49,201,0,0,0,0,0,0.0d0,0.0d0)
1614 msg = 'DASPK-- VALUE (=I1) OF IDID AND NO APPROPRIATE'
1615 CALL xerrwd(msg,47,202,0,1,idid,0,0,0.0d0,0.0d0)
1616 msg = 'DASPK-- ACTION WAS TAKEN. RUN TERMINATED'
1617 CALL xerrwd(msg,41,203,1,0,0,0,0,0.0d0,0.0d0)
1618 RETURN
1619110 CONTINUE
1620C
1621C-----------------------------------------------------------------------
1622C This block is executed on all calls.
1623C
1624C Counters are saved for later checks of performance.
1625C Then the error tolerance parameters are checked, and the
1626C work array pointers are set.
1627C-----------------------------------------------------------------------
1628C
1629200 CONTINUE
1630C
1631C Save counters for use later.
1632C
1633 iwork(lnstl)=iwork(lnst)
1634 nli0 = iwork(lnli)
1635 nni0 = iwork(lnni)
1636 ncfn0 = iwork(lncfn)
1637 ncfl0 = iwork(lncfl)
1638 nwarn = 0
1639C
1640C Check RTOL and ATOL.
1641C
1642 nzflg = 0
1643 rtoli = rtol(1)
1644 atoli = atol(1)
1645 DO 210 i=1,neq
1646 IF (info(2) .EQ. 1) rtoli = rtol(i)
1647 IF (info(2) .EQ. 1) atoli = atol(i)
1648 IF (rtoli .GT. 0.0d0 .OR. atoli .GT. 0.0d0) nzflg = 1
1649 IF (rtoli .LT. 0.0d0) GO TO 706
1650 IF (atoli .LT. 0.0d0) GO TO 707
1651210 CONTINUE
1652 IF (nzflg .EQ. 0) GO TO 708
1653C
1654C Set pointers to RWORK and IWORK segments.
1655C For direct methods, SAVR is not used.
1656C
1657 iwork(llciwp) = lid + lenid
1658 lsavr = ldelta
1659 IF (info(12) .NE. 0) lsavr = ldelta + neq
1660 le = lsavr + neq
1661 lwt = le + neq
1662 lvt = lwt
1663 IF (info(16) .NE. 0) lvt = lwt + neq
1664 lphi = lvt + neq
1665 lwm = lphi + (iwork(lmxord)+1)*neq
1666 IF (info(1) .EQ. 1) GO TO 400
1667C
1668C-----------------------------------------------------------------------
1669C This block is executed on the initial call only.
1670C Set the initial step size, the error weight vector, and PHI.
1671C Compute unknown initial components of Y and YPRIME, if requested.
1672C-----------------------------------------------------------------------
1673C
1674300 CONTINUE
1675 tn=t
1676 idid=1
1677C
1678C Set error weight array WT and altered weight array VT.
1679C
1680 CALL ddawts(neq,info(2),rtol,atol,y,rwork(lwt),rpar,ipar)
1681 CALL dinvwt(neq,rwork(lwt),ier)
1682 IF (ier .NE. 0) GO TO 713
1683 IF (info(16) .NE. 0) THEN
1684 DO 305 i = 1, neq
1685 305 rwork(lvt+i-1) = max(iwork(lid+i-1),0)*rwork(lwt+i-1)
1686 ENDIF
1687C
1688C Compute unit roundoff and HMIN.
1689C
1690 uround = d1mach(4)
1691 rwork(lround) = uround
1692 hmin = 4.0d0*uround*max(abs(t),abs(tout))
1693C
1694C Set/check STPTOL control for initial condition calculation.
1695C
1696 IF (info(11) .NE. 0) THEN
1697 IF( info(17) .EQ. 0) THEN
1698 rwork(lstol) = uround**.6667d0
1699 ELSE
1700 IF (rwork(lstol) .LE. 0.0d0) GO TO 725
1701 ENDIF
1702 ENDIF
1703C
1704C Compute EPCON and square root of NEQ and its reciprocal, used
1705C inside iterative solver.
1706C
1707 rwork(lepcon) = 0.33d0
1708 floatn = neq
1709 rwork(lsqrn) = sqrt(floatn)
1710 rwork(lrsqrn) = 1.d0/rwork(lsqrn)
1711C
1712C Check initial interval to see that it is long enough.
1713C
1714 tdist = abs(tout - t)
1715 IF(tdist .LT. hmin) GO TO 714
1716C
1717C Check H0, if this was input.
1718C
1719 IF (info(8) .EQ. 0) GO TO 310
1720 h0 = rwork(lh)
1721 IF ((tout - t)*h0 .LT. 0.0d0) GO TO 711
1722 IF (h0 .EQ. 0.0d0) GO TO 712
1723 GO TO 320
1724310 CONTINUE
1725C
1726C Compute initial stepsize, to be used by either
1727C DDSTP or DDASIC, depending on INFO(11).
1728C
1729 h0 = 0.001d0*tdist
1730 ypnorm = ddwnrm(neq,yprime,rwork(lvt),rpar,ipar)
1731 IF (ypnorm .GT. 0.5d0/h0) h0 = 0.5d0/ypnorm
1732 h0 = sign(h0,tout-t)
1733C
1734C Adjust H0 if necessary to meet HMAX bound.
1735C
1736320 IF (info(7) .EQ. 0) GO TO 330
1737 rh = abs(h0)/rwork(lhmax)
1738 IF (rh .GT. 1.0d0) h0 = h0/rh
1739C
1740C Check against TSTOP, if applicable.
1741C
1742330 IF (info(4) .EQ. 0) GO TO 340
1743 tstop = rwork(ltstop)
1744 IF ((tstop - t)*h0 .LT. 0.0d0) GO TO 715
1745 IF ((t + h0 - tstop)*h0 .GT. 0.0d0) h0 = tstop - t
1746 IF ((tstop - tout)*h0 .LT. 0.0d0) GO TO 709
1747C
1748340 IF (info(11) .EQ. 0) GO TO 370
1749C
1750C Compute unknown components of initial Y and YPRIME, depending
1751C on INFO(11) and INFO(12). INFO(12) represents the nonlinear
1752C solver type (direct/Krylov). Pass the name of the specific
1753C nonlinear solver, depending on INFO(12). The location of the work
1754C arrays SAVR, YIC, YPIC, PWK also differ in the two cases.
1755C
1756 nwt = 1
1757 epconi = rwork(lepin)*rwork(lepcon)
1758350 IF (info(12) .EQ. 0) THEN
1759 lyic = lphi + 2*neq
1760 lypic = lyic + neq
1761 lpwk = lypic
1762 CALL ddasic(tn,y,yprime,neq,info(11),iwork(lid),
1763 * res,jac,psol,h0,rwork(lwt),nwt,idid,rpar,ipar,
1764 * rwork(lphi),rwork(lsavr),rwork(ldelta),rwork(le),
1765 * rwork(lyic),rwork(lypic),rwork(lpwk),rwork(lwm),iwork(liwm),
1766 * hmin,rwork(lround),rwork(lepli),rwork(lsqrn),rwork(lrsqrn),
1767 * epconi,rwork(lstol),info(15),icnflg,iwork(licns),ddasid)
1768 ELSE IF (info(12) .EQ. 1) THEN
1769 lyic = lwm
1770 lypic = lyic + neq
1771 lpwk = lypic + neq
1772 CALL ddasic(tn,y,yprime,neq,info(11),iwork(lid),
1773 * res,jac,psol,h0,rwork(lwt),nwt,idid,rpar,ipar,
1774 * rwork(lphi),rwork(lsavr),rwork(ldelta),rwork(le),
1775 * rwork(lyic),rwork(lypic),rwork(lpwk),rwork(lwm),iwork(liwm),
1776 * hmin,rwork(lround),rwork(lepli),rwork(lsqrn),rwork(lrsqrn),
1777 * epconi,rwork(lstol),info(15),icnflg,iwork(licns),ddasik)
1778 ENDIF
1779C
1780 IF (idid .LT. 0) GO TO 600
1781C
1782C DDASIC was successful. If this was the first call to DDASIC,
1783C update the WT array (with the current Y) and call it again.
1784C
1785 IF (nwt .EQ. 2) GO TO 355
1786 nwt = 2
1787 CALL ddawts(neq,info(2),rtol,atol,y,rwork(lwt),rpar,ipar)
1788 CALL dinvwt(neq,rwork(lwt),ier)
1789 IF (ier .NE. 0) GO TO 713
1790 GO TO 350
1791C
1792C If INFO(14) = 1, return now with IDID = 4.
1793C
1794355 IF (info(14) .EQ. 1) THEN
1795 idid = 4
1796 h = h0
1797 IF (info(11) .EQ. 1) rwork(lhold) = h0
1798 GO TO 590
1799 ENDIF
1800C
1801C Update the WT and VT arrays one more time, with the new Y.
1802C
1803 CALL ddawts(neq,info(2),rtol,atol,y,rwork(lwt),rpar,ipar)
1804 CALL dinvwt(neq,rwork(lwt),ier)
1805 IF (ier .NE. 0) GO TO 713
1806 IF (info(16) .NE. 0) THEN
1807 DO 357 i = 1, neq
1808 357 rwork(lvt+i-1) = max(iwork(lid+i-1),0)*rwork(lwt+i-1)
1809 ENDIF
1810C
1811C Reset the initial stepsize to be used by DDSTP.
1812C Use H0, if this was input. Otherwise, recompute H0,
1813C and adjust it if necessary to meet HMAX bound.
1814C
1815 IF (info(8) .NE. 0) THEN
1816 h0 = rwork(lh)
1817 GO TO 360
1818 ENDIF
1819C
1820 h0 = 0.001d0*tdist
1821 ypnorm = ddwnrm(neq,yprime,rwork(lvt),rpar,ipar)
1822 IF (ypnorm .GT. 0.5d0/h0) h0 = 0.5d0/ypnorm
1823 h0 = sign(h0,tout-t)
1824C
1825360 IF (info(7) .NE. 0) THEN
1826 rh = abs(h0)/rwork(lhmax)
1827 IF (rh .GT. 1.0d0) h0 = h0/rh
1828 ENDIF
1829C
1830C Check against TSTOP, if applicable.
1831C
1832 IF (info(4) .NE. 0) THEN
1833 tstop = rwork(ltstop)
1834 IF ((t + h0 - tstop)*h0 .GT. 0.0d0) h0 = tstop - t
1835 ENDIF
1836C
1837C Load H and RWORK(LH) with H0.
1838C
1839370 h = h0
1840 rwork(lh) = h
1841C
1842C Load Y and H*YPRIME into PHI(*,1) and PHI(*,2).
1843C
1844 itemp = lphi + neq
1845 DO 380 i = 1,neq
1846 rwork(lphi + i - 1) = y(i)
1847380 rwork(itemp + i - 1) = h*yprime(i)
1848C
1849 GO TO 500
1850C
1851C-----------------------------------------------------------------------
1852C This block is for continuation calls only.
1853C Its purpose is to check stop conditions before taking a step.
1854C Adjust H if necessary to meet HMAX bound.
1855C-----------------------------------------------------------------------
1856C
1857400 CONTINUE
1858 uround=rwork(lround)
1859 done = .false.
1860 tn=rwork(ltn)
1861 h=rwork(lh)
1862 IF(info(7) .EQ. 0) GO TO 410
1863 rh = abs(h)/rwork(lhmax)
1864 IF(rh .GT. 1.0d0) h = h/rh
1865410 CONTINUE
1866 IF(t .EQ. tout) GO TO 719
1867 IF((t - tout)*h .GT. 0.0d0) GO TO 711
1868 IF(info(4) .EQ. 1) GO TO 430
1869 IF(info(3) .EQ. 1) GO TO 420
1870 IF((tn-tout)*h.LT.0.0d0)GO TO 490
1871 CALL ddatrp(tn,tout,y,yprime,neq,iwork(lkold),
1872 * rwork(lphi),rwork(lpsi))
1873 t=tout
1874 idid = 3
1875 done = .true.
1876 GO TO 490
1877420 IF((tn-t)*h .LE. 0.0d0) GO TO 490
1878 IF((tn - tout)*h .GT. 0.0d0) GO TO 425
1879 CALL ddatrp(tn,tn,y,yprime,neq,iwork(lkold),
1880 * rwork(lphi),rwork(lpsi))
1881 t = tn
1882 idid = 1
1883 done = .true.
1884 GO TO 490
1885425 CONTINUE
1886 CALL ddatrp(tn,tout,y,yprime,neq,iwork(lkold),
1887 * rwork(lphi),rwork(lpsi))
1888 t = tout
1889 idid = 3
1890 done = .true.
1891 GO TO 490
1892430 IF(info(3) .EQ. 1) GO TO 440
1893 tstop=rwork(ltstop)
1894 IF((tn-tstop)*h.GT.0.0d0) GO TO 715
1895 IF((tstop-tout)*h.LT.0.0d0)GO TO 709
1896 IF((tn-tout)*h.LT.0.0d0)GO TO 450
1897 CALL ddatrp(tn,tout,y,yprime,neq,iwork(lkold),
1898 * rwork(lphi),rwork(lpsi))
1899 t=tout
1900 idid = 3
1901 done = .true.
1902 GO TO 490
1903440 tstop = rwork(ltstop)
1904 IF((tn-tstop)*h .GT. 0.0d0) GO TO 715
1905 IF((tstop-tout)*h .LT. 0.0d0) GO TO 709
1906 IF((tn-t)*h .LE. 0.0d0) GO TO 450
1907 IF((tn - tout)*h .GT. 0.0d0) GO TO 445
1908 CALL ddatrp(tn,tn,y,yprime,neq,iwork(lkold),
1909 * rwork(lphi),rwork(lpsi))
1910 t = tn
1911 idid = 1
1912 done = .true.
1913 GO TO 490
1914445 CONTINUE
1915 CALL ddatrp(tn,tout,y,yprime,neq,iwork(lkold),
1916 * rwork(lphi),rwork(lpsi))
1917 t = tout
1918 idid = 3
1919 done = .true.
1920 GO TO 490
1921450 CONTINUE
1922C
1923C Check whether we are within roundoff of TSTOP.
1924C
1925 IF(abs(tn-tstop).GT.100.0d0*uround*
1926 * (abs(tn)+abs(h)))GO TO 460
1927 CALL ddatrp(tn,tstop,y,yprime,neq,iwork(lkold),
1928 * rwork(lphi),rwork(lpsi))
1929 idid=2
1930 t=tstop
1931 done = .true.
1932 GO TO 490
1933460 tnext=tn+h
1934 IF((tnext-tstop)*h.LE.0.0d0)GO TO 490
1935 h=tstop-tn
1936 rwork(lh)=h
1937C
1938490 IF (done) GO TO 590
1939C
1940C-----------------------------------------------------------------------
1941C The next block contains the call to the one-step integrator DDSTP.
1942C This is a looping point for the integration steps.
1943C Check for too many steps.
1944C Check for poor Newton/Krylov performance.
1945C Update WT. Check for too much accuracy requested.
1946C Compute minimum stepsize.
1947C-----------------------------------------------------------------------
1948C
1949500 CONTINUE
1950C
1951C Check for too many steps.
1952C
1953 IF((iwork(lnst)-iwork(lnstl)).LT.500) GO TO 505
1954 idid=-1
1955 GO TO 527
1956C
1957C Check for poor Newton/Krylov performance.
1958C
1959505 IF (info(12) .EQ. 0) GO TO 510
1960 nstd = iwork(lnst) - iwork(lnstl)
1961 nnid = iwork(lnni) - nni0
1962 IF (nstd .LT. 10 .OR. nnid .EQ. 0) GO TO 510
1963 avlin = real(iwork(lnli) - nli0)/real(nnid)
1964 rcfn = real(iwork(lncfn) - ncfn0)/real(nstd)
1965 rcfl = real(iwork(lncfl) - ncfl0)/real(nnid)
1966 fmaxl = iwork(lmaxl)
1967 lavl = avlin .GT. fmaxl
1968 lcfn = rcfn .GT. 0.9d0
1969 lcfl = rcfl .GT. 0.9d0
1970 lwarn = lavl .OR. lcfn .OR. lcfl
1971 IF (.NOT.lwarn) GO TO 510
1972 nwarn = nwarn + 1
1973 IF (nwarn .GT. 10) GO TO 510
1974 IF (lavl) THEN
1975 msg = 'DASPK-- Warning. Poor iterative algorithm performance '
1976 CALL xerrwd (msg, 56, 501, 0, 0, 0, 0, 0, 0.0d0, 0.0d0)
1977 msg = ' at T = R1. Average no. of linear iterations = R2 '
1978 CALL xerrwd (msg, 56, 501, 0, 0, 0, 0, 2, tn, avlin)
1979 ENDIF
1980 IF (lcfn) THEN
1981 msg = 'DASPK-- Warning. Poor iterative algorithm performance '
1982 CALL xerrwd (msg, 56, 502, 0, 0, 0, 0, 0, 0.0d0, 0.0d0)
1983 msg = ' at T = R1. Nonlinear convergence failure rate = R2'
1984 CALL xerrwd (msg, 56, 502, 0, 0, 0, 0, 2, tn, rcfn)
1985 ENDIF
1986 IF (lcfl) THEN
1987 msg = 'DASPK-- Warning. Poor iterative algorithm performance '
1988 CALL xerrwd (msg, 56, 503, 0, 0, 0, 0, 0, 0.0d0, 0.0d0)
1989 msg = ' at T = R1. Linear convergence failure rate = R2 '
1990 CALL xerrwd (msg, 56, 503, 0, 0, 0, 0, 2, tn, rcfl)
1991 ENDIF
1992C
1993C Update WT and VT, if this is not the first call.
1994C
1995510 CALL ddawts(neq,info(2),rtol,atol,rwork(lphi),rwork(lwt),
1996 * rpar,ipar)
1997 CALL dinvwt(neq,rwork(lwt),ier)
1998 IF (ier .NE. 0) THEN
1999 idid = -3
2000 GO TO 527
2001 ENDIF
2002 IF (info(16) .NE. 0) THEN
2003 DO 515 i = 1, neq
2004 515 rwork(lvt+i-1) = max(iwork(lid+i-1),0)*rwork(lwt+i-1)
2005 ENDIF
2006C
2007C Test for too much accuracy requested.
2008C
2009 r = ddwnrm(neq,rwork(lphi),rwork(lwt),rpar,ipar)*100.0d0*uround
2010 IF (r .LE. 1.0d0) GO TO 525
2011C
2012C Multiply RTOL and ATOL by R and return.
2013C
2014 IF(info(2).EQ.1)GO TO 523
2015 rtol(1)=r*rtol(1)
2016 atol(1)=r*atol(1)
2017 idid=-2
2018 GO TO 527
2019523 DO 524 i=1,neq
2020 rtol(i)=r*rtol(i)
2021524 atol(i)=r*atol(i)
2022 idid=-2
2023 GO TO 527
2024525 CONTINUE
2025C
2026C Compute minimum stepsize.
2027C
2028 hmin=4.0d0*uround*max(abs(tn),abs(tout))
2029C
2030C Test H vs. HMAX
2031 IF (info(7) .NE. 0) THEN
2032 rh = abs(h)/rwork(lhmax)
2033 IF (rh .GT. 1.0d0) h = h/rh
2034 ENDIF
2035C
2036C Call the one-step integrator.
2037C Note that INFO(12) represents the nonlinear solver type.
2038C Pass the required nonlinear solver, depending upon INFO(12).
2039C
2040 IF (info(12) .EQ. 0) THEN
2041 CALL ddstp(tn,y,yprime,neq,
2042 * res,jac,psol,h,rwork(lwt),rwork(lvt),info(1),idid,rpar,ipar,
2043 * rwork(lphi),rwork(lsavr),rwork(ldelta),rwork(le),
2044 * rwork(lwm),iwork(liwm),
2045 * rwork(lalpha),rwork(lbeta),rwork(lgamma),
2046 * rwork(lpsi),rwork(lsigma),
2047 * rwork(lcj),rwork(lcjold),rwork(lhold),rwork(ls),hmin,
2048 * rwork(lround), rwork(lepli),rwork(lsqrn),rwork(lrsqrn),
2049 * rwork(lepcon), iwork(lphase),iwork(ljcalc),info(15),
2050 * iwork(lk), iwork(lkold),iwork(lns),nonneg,info(12),
2051 * dnedd)
2052 ELSE IF (info(12) .EQ. 1) THEN
2053 CALL ddstp(tn,y,yprime,neq,
2054 * res,jac,psol,h,rwork(lwt),rwork(lvt),info(1),idid,rpar,ipar,
2055 * rwork(lphi),rwork(lsavr),rwork(ldelta),rwork(le),
2056 * rwork(lwm),iwork(liwm),
2057 * rwork(lalpha),rwork(lbeta),rwork(lgamma),
2058 * rwork(lpsi),rwork(lsigma),
2059 * rwork(lcj),rwork(lcjold),rwork(lhold),rwork(ls),hmin,
2060 * rwork(lround), rwork(lepli),rwork(lsqrn),rwork(lrsqrn),
2061 * rwork(lepcon), iwork(lphase),iwork(ljcalc),info(15),
2062 * iwork(lk), iwork(lkold),iwork(lns),nonneg,info(12),
2063 * dnedk)
2064 ENDIF
2065C
2066527 IF(idid.LT.0)GO TO 600
2067C
2068C-----------------------------------------------------------------------
2069C This block handles the case of a successful return from DDSTP
2070C (IDID=1). Test for stop conditions.
2071C-----------------------------------------------------------------------
2072C
2073 IF(info(4).NE.0)GO TO 540
2074 IF(info(3).NE.0)GO TO 530
2075 IF((tn-tout)*h.LT.0.0d0)GO TO 500
2076 CALL ddatrp(tn,tout,y,yprime,neq,
2077 * iwork(lkold),rwork(lphi),rwork(lpsi))
2078 idid=3
2079 t=tout
2080 GO TO 580
2081530 IF((tn-tout)*h.GE.0.0d0)GO TO 535
2082 t=tn
2083 idid=1
2084 GO TO 580
2085535 CALL ddatrp(tn,tout,y,yprime,neq,
2086 * iwork(lkold),rwork(lphi),rwork(lpsi))
2087 idid=3
2088 t=tout
2089 GO TO 580
2090540 IF(info(3).NE.0)GO TO 550
2091 IF((tn-tout)*h.LT.0.0d0)GO TO 542
2092 CALL ddatrp(tn,tout,y,yprime,neq,
2093 * iwork(lkold),rwork(lphi),rwork(lpsi))
2094 t=tout
2095 idid=3
2096 GO TO 580
2097542 IF(abs(tn-tstop).LE.100.0d0*uround*
2098 * (abs(tn)+abs(h)))GO TO 545
2099 tnext=tn+h
2100 IF((tnext-tstop)*h.LE.0.0d0)GO TO 500
2101 h=tstop-tn
2102 GO TO 500
2103545 CALL ddatrp(tn,tstop,y,yprime,neq,
2104 * iwork(lkold),rwork(lphi),rwork(lpsi))
2105 idid=2
2106 t=tstop
2107 GO TO 580
2108550 IF((tn-tout)*h.GE.0.0d0)GO TO 555
2109 IF(abs(tn-tstop).LE.100.0d0*uround*(abs(tn)+abs(h)))GO TO 552
2110 t=tn
2111 idid=1
2112 GO TO 580
2113552 CALL ddatrp(tn,tstop,y,yprime,neq,
2114 * iwork(lkold),rwork(lphi),rwork(lpsi))
2115 idid=2
2116 t=tstop
2117 GO TO 580
2118555 CALL ddatrp(tn,tout,y,yprime,neq,
2119 * iwork(lkold),rwork(lphi),rwork(lpsi))
2120 t=tout
2121 idid=3
2122580 CONTINUE
2123C
2124C-----------------------------------------------------------------------
2125C All successful returns from DDASPK are made from this block.
2126C-----------------------------------------------------------------------
2127C
2128590 CONTINUE
2129 rwork(ltn)=tn
2130 rwork(lh)=h
2131 RETURN
2132C
2133C-----------------------------------------------------------------------
2134C This block handles all unsuccessful returns other than for
2135C illegal input.
2136C-----------------------------------------------------------------------
2137C
2138600 CONTINUE
2139 itemp = -idid
2140 GO TO (610,620,630,700,655,640,650,660,670,675,
2141 * 680,685,690,695), itemp
2142C
2143C The maximum number of steps was taken before
2144C reaching tout.
2145C
2146610 msg = 'DASPK-- AT CURRENT T (=R1) 500 STEPS'
2147 CALL xerrwd(msg,38,610,0,0,0,0,1,tn,0.0d0)
2148 msg = 'DASPK-- TAKEN ON THIS CALL BEFORE REACHING TOUT'
2149 CALL xerrwd(msg,48,611,0,0,0,0,0,0.0d0,0.0d0)
2150 GO TO 700
2151C
2152C Too much accuracy for machine precision.
2153C
2154620 msg = 'DASPK-- AT T (=R1) TOO MUCH ACCURACY REQUESTED'
2155 CALL xerrwd(msg,47,620,0,0,0,0,1,tn,0.0d0)
2156 msg = 'DASPK-- FOR PRECISION OF MACHINE. RTOL AND ATOL'
2157 CALL xerrwd(msg,48,621,0,0,0,0,0,0.0d0,0.0d0)
2158 msg = 'DASPK-- WERE INCREASED TO APPROPRIATE VALUES'
2159 CALL xerrwd(msg,45,622,0,0,0,0,0,0.0d0,0.0d0)
2160 GO TO 700
2161C
2162C WT(I) .LE. 0.0D0 for some I (not at start of problem).
2163C
2164630 msg = 'DASPK-- AT T (=R1) SOME ELEMENT OF WT'
2165 CALL xerrwd(msg,38,630,0,0,0,0,1,tn,0.0d0)
2166 msg = .LE.'DASPK-- HAS BECOME 0.0'
2167 CALL xerrwd(msg,28,631,0,0,0,0,0,0.0d0,0.0d0)
2168 GO TO 700
2169C
2170C Error test failed repeatedly or with H=HMIN.
2171C
2172640 msg = 'DASPK-- AT T (=R1) AND STEPSIZE H (=R2) THE'
2173 CALL xerrwd(msg,44,640,0,0,0,0,2,tn,h)
2174 msg='DASPK-- ERROR TEST FAILED REPEATEDLY OR WITH ABS(H)=HMIN'
2175 CALL xerrwd(msg,57,641,0,0,0,0,0,0.0d0,0.0d0)
2176 GO TO 700
2177C
2178C Nonlinear solver failed to converge repeatedly or with H=HMIN.
2179C
2180650 msg = 'DASPK-- AT T (=R1) AND STEPSIZE H (=R2) THE'
2181 CALL xerrwd(msg,44,650,0,0,0,0,2,tn,h)
2182 msg = 'DASPK-- NONLINEAR SOLVER FAILED TO CONVERGE'
2183 CALL xerrwd(msg,44,651,0,0,0,0,0,0.0d0,0.0d0)
2184 msg = 'DASPK-- REPEATEDLY OR WITH ABS(H)=HMIN'
2185 CALL xerrwd(msg,40,652,0,0,0,0,0,0.0d0,0.0d0)
2186 GO TO 700
2187C
2188C The preconditioner had repeated failures.
2189C
2190655 msg = 'DASPK-- AT T (=R1) AND STEPSIZE H (=R2) THE'
2191 CALL xerrwd(msg,44,655,0,0,0,0,2,tn,h)
2192 msg = 'DASPK-- PRECONDITIONER HAD REPEATED FAILURES.'
2193 CALL xerrwd(msg,46,656,0,0,0,0,0,0.0d0,0.0d0)
2194 GO TO 700
2195C
2196C The iteration matrix is singular.
2197C
2198660 msg = 'DASPK-- AT T (=R1) AND STEPSIZE H (=R2) THE'
2199 CALL xerrwd(msg,44,660,0,0,0,0,2,tn,h)
2200 msg = 'DASPK-- ITERATION MATRIX IS SINGULAR.'
2201 CALL xerrwd(msg,38,661,0,0,0,0,0,0.0d0,0.0d0)
2202 GO TO 700
2203C
2204C Nonlinear system failure preceded by error test failures.
2205C
2206670 msg = 'DASPK-- AT T (=R1) AND STEPSIZE H (=R2) THE'
2207 CALL xerrwd(msg,44,670,0,0,0,0,2,tn,h)
2208 msg = 'DASPK-- NONLINEAR SOLVER COULD NOT CONVERGE.'
2209 CALL xerrwd(msg,45,671,0,0,0,0,0,0.0d0,0.0d0)
2210 msg = 'DASPK-- ALSO, THE ERROR TEST FAILED REPEATEDLY.'
2211 CALL xerrwd(msg,49,672,0,0,0,0,0,0.0d0,0.0d0)
2212 GO TO 700
2213C
2214C Nonlinear system failure because IRES = -1.
2215C
2216675 msg = 'DASPK-- AT T (=R1) AND STEPSIZE H (=R2) THE'
2217 CALL xerrwd(msg,44,675,0,0,0,0,2,tn,h)
2218 msg = 'DASPK-- NONLINEAR SYSTEM SOLVER COULD NOT CONVERGE'
2219 CALL xerrwd(msg,51,676,0,0,0,0,0,0.0d0,0.0d0)
2220 msg = 'DASPK-- BECAUSE IRES WAS EQUAL TO MINUS ONE'
2221 CALL xerrwd(msg,44,677,0,0,0,0,0,0.0d0,0.0d0)
2222 GO TO 700
2223C
2224C Failure because IRES = -2.
2225C
2226680 msg = 'DASPK-- AT T (=R1) AND STEPSIZE H (=R2)'
2227 CALL xerrwd(msg,40,680,0,0,0,0,2,tn,h)
2228 msg = 'DASPK-- IRES WAS EQUAL TO MINUS TWO'
2229 CALL xerrwd(msg,36,681,0,0,0,0,0,0.0d0,0.0d0)
2230 GO TO 700
2231C
2232C Failed to compute initial YPRIME.
2233C
2234685 msg = 'DASPK-- AT T (=R1) AND STEPSIZE H (=R2) THE'
2235 CALL xerrwd(msg,44,685,0,0,0,0,0,0.0d0,0.0d0)
2236 msg = 'DASPK-- INITIAL (Y,YPRIME) COULD NOT BE COMPUTED'
2237 CALL xerrwd(msg,49,686,0,0,0,0,2,tn,h0)
2238 GO TO 700
2239C
2240C Failure because IER was negative from PSOL.
2241C
2242690 msg = 'DASPK-- AT T (=R1) AND STEPSIZE H (=R2)'
2243 CALL xerrwd(msg,40,690,0,0,0,0,2,tn,h)
2244 msg = 'DASPK-- IER WAS NEGATIVE FROM PSOL'
2245 CALL xerrwd(msg,35,691,0,0,0,0,0,0.0d0,0.0d0)
2246 GO TO 700
2247C
2248C Failure because the linear system solver could not converge.
2249C
2250695 msg = 'DASPK-- AT T (=R1) AND STEPSIZE H (=R2) THE'
2251 CALL xerrwd(msg,44,695,0,0,0,0,2,tn,h)
2252 msg = 'DASPK-- LINEAR SYSTEM SOLVER COULD NOT CONVERGE.'
2253 CALL xerrwd(msg,50,696,0,0,0,0,0,0.0d0,0.0d0)
2254 GO TO 700
2255C
2256C
2257700 CONTINUE
2258 info(1)=-1
2259 t=tn
2260 rwork(ltn)=tn
2261 rwork(lh)=h
2262 RETURN
2263C
2264C-----------------------------------------------------------------------
2265C This block handles all error returns due to illegal input,
2266C as detected before calling DDSTP.
2267C First the error message routine is called. If this happens
2268C twice in succession, execution is terminated.
2269C-----------------------------------------------------------------------
2270C
2271701 msg = 'DASPK-- ELEMENT (=I1) OF INFO VECTOR IS NOT VALID'
2272 CALL xerrwd(msg,50,1,0,1,itemp,0,0,0.0d0,0.0d0)
2273 GO TO 750
2274702 msg = .LE.'DASPK-- NEQ (=I1) 0'
2275 CALL xerrwd(msg,25,2,0,1,neq,0,0,0.0d0,0.0d0)
2276 GO TO 750
2277703 msg = 'DASPK-- MAXORD (=I1) NOT IN RANGE'
2278 CALL xerrwd(msg,34,3,0,1,mxord,0,0,0.0d0,0.0d0)
2279 GO TO 750
2280704 msg='DASPK-- RWORK LENGTH NEEDED, LENRW (=I1), EXCEEDS LRW (=I2)'
2281 CALL xerrwd(msg,60,4,0,2,lenrw,lrw,0,0.0d0,0.0d0)
2282 GO TO 750
2283705 msg='DASPK-- IWORK LENGTH NEEDED, LENIW (=I1), EXCEEDS LIW (=I2)'
2284 CALL xerrwd(msg,60,5,0,2,leniw,liw,0,0.0d0,0.0d0)
2285 GO TO 750
2286706 msg = .LT.'DASPK-- SOME ELEMENT OF RTOL IS 0'
2287 CALL xerrwd(msg,39,6,0,0,0,0,0,0.0d0,0.0d0)
2288 GO TO 750
2289707 msg = .LT.'DASPK-- SOME ELEMENT OF ATOL IS 0'
2290 CALL xerrwd(msg,39,7,0,0,0,0,0,0.0d0,0.0d0)
2291 GO TO 750
2292708 msg = 'DASPK-- ALL ELEMENTS OF RTOL AND ATOL ARE ZERO'
2293 CALL xerrwd(msg,47,8,0,0,0,0,0,0.0d0,0.0d0)
2294 GO TO 750
2295709 msg='DASPK-- INFO(4) = 1 AND TSTOP (=R1) BEHIND TOUT (=R2)'
2296 CALL xerrwd(msg,54,9,0,0,0,0,2,tstop,tout)
2297 GO TO 750
2298710 msg = .LT.'DASPK-- HMAX (=R1) 0.0'
2299 CALL xerrwd(msg,28,10,0,0,0,0,1,hmax,0.0d0)
2300 GO TO 750
2301711 msg = 'DASPK-- TOUT (=R1) BEHIND T (=R2)'
2302 CALL xerrwd(msg,34,11,0,0,0,0,2,tout,t)
2303 GO TO 750
2304712 msg = 'DASPK-- INFO(8)=1 AND H0=0.0'
2305 CALL xerrwd(msg,29,12,0,0,0,0,0,0.0d0,0.0d0)
2306 GO TO 750
2307713 msg = .LE.'DASPK-- SOME ELEMENT OF WT IS 0.0'
2308 CALL xerrwd(msg,39,13,0,0,0,0,0,0.0d0,0.0d0)
2309 GO TO 750
2310714 msg='DASPK-- TOUT (=R1) TOO CLOSE TO T (=R2) TO START INTEGRATION'
2311 CALL xerrwd(msg,60,14,0,0,0,0,2,tout,t)
2312 GO TO 750
2313715 msg = 'DASPK-- INFO(4)=1 AND TSTOP (=R1) BEHIND T (=R2)'
2314 CALL xerrwd(msg,49,15,0,0,0,0,2,tstop,t)
2315 GO TO 750
2316717 msg = .LT..GT.'DASPK-- ML (=I1) ILLEGAL. EITHER 0 OR NEQ'
2317 CALL xerrwd(msg,52,17,0,1,iwork(lml),0,0,0.0d0,0.0d0)
2318 GO TO 750
2319718 msg = .LT..GT.'DASPK-- MU (=I1) ILLEGAL. EITHER 0 OR NEQ'
2320 CALL xerrwd(msg,52,18,0,1,iwork(lmu),0,0,0.0d0,0.0d0)
2321 GO TO 750
2322719 msg = 'DASPK-- TOUT (=R1) IS EQUAL TO T (=R2)'
2323 CALL xerrwd(msg,39,19,0,0,0,0,2,tout,t)
2324 GO TO 750
2325720 msg = .LT..GT.'DASPK-- MAXL (=I1) ILLEGAL. EITHER 1 OR NEQ'
2326 CALL xerrwd(msg,54,20,0,1,iwork(lmaxl),0,0,0.0d0,0.0d0)
2327 GO TO 750
2328721 msg = .LT..GT.'DASPK-- KMP (=I1) ILLEGAL. EITHER 1 OR MAXL'
2329 CALL xerrwd(msg,54,21,0,1,iwork(lkmp),0,0,0.0d0,0.0d0)
2330 GO TO 750
2331722 msg = .LT.'DASPK-- NRMAX (=I1) ILLEGAL. 0'
2332 CALL xerrwd(msg,36,22,0,1,iwork(lnrmax),0,0,0.0d0,0.0d0)
2333 GO TO 750
2334723 msg = .LE..GE.'DASPK-- EPLI (=R1) ILLEGAL. EITHER 0.D0 OR 1.D0'
2335 CALL xerrwd(msg,58,23,0,0,0,0,1,rwork(lepli),0.0d0)
2336 GO TO 750
2337724 msg = .NE.'DASPK-- ILLEGAL IWORK VALUE FOR INFO(11) 0'
2338 CALL xerrwd(msg,48,24,0,0,0,0,0,0.0d0,0.0d0)
2339 GO TO 750
2340725 msg = 'DASPK-- ONE OF THE INPUTS FOR INFO(17) = 1 IS ILLEGAL'
2341 CALL xerrwd(msg,54,25,0,0,0,0,0,0.0d0,0.0d0)
2342 GO TO 750
2343726 msg = .NE.'DASPK-- ILLEGAL IWORK VALUE FOR INFO(10) 0'
2344 CALL xerrwd(msg,48,26,0,0,0,0,0,0.0d0,0.0d0)
2345 GO TO 750
2346727 msg = 'DASPK-- Y(I) AND IWORK(40+I) (I=I1) INCONSISTENT'
2347 CALL xerrwd(msg,49,27,0,1,iret,0,0,0.0d0,0.0d0)
2348 GO TO 750
2349750 IF(info(1).EQ.-1) GO TO 760
2350 info(1)=-1
2351 idid=-33
2352 RETURN
2353760 msg = 'DASPK-- REPEATED OCCURRENCES OF ILLEGAL INPUT'
2354 CALL xerrwd(msg,46,701,0,0,0,0,0,0.0d0,0.0d0)
2355770 msg = 'DASPK-- RUN TERMINATED. APPARENT INFINITE LOOP'
2356 CALL xerrwd(msg,47,702,1,0,0,0,0,0.0d0,0.0d0)
2357 RETURN
2358C
2359C------END OF SUBROUTINE DDASPK-----------------------------------------
2360 END
charNDArray max(char d, const charNDArray &m)
Definition chNDArray.cc:230
charNDArray min(char d, const charNDArray &m)
Definition chNDArray.cc:207
double precision function d1mach(i)
Definition d1mach.f:23
ColumnVector real(const ComplexColumnVector &a)
subroutine dcnst0(neq, y, icnstr, iret)
Definition dcnst0.f:6
subroutine ddasic(x, y, yprime, neq, icopt, id, res, jac, psol, h, wt, nic, idid, rpar, ipar, phi, savr, delta, e, yic, ypic, pwk, wm, iwm, hmin, uround, epli, sqrtn, rsqrtn, epconi, stptol, jflg, icnflg, icnstr, nlsic)
Definition ddasic.f:9
subroutine ddasid(x, y, yprime, neq, icopt, id, res, jacd, pdum, h, wt, jsdum, rpar, ipar, dumsvr, delta, r, yic, ypic, dumpwk, wm, iwm, cj, uround, dume, dums, dumr, epcon, ratemx, stptol, jfdum, icnflg, icnstr, iernls)
Definition ddasid.f:9
subroutine ddasik(x, y, yprime, neq, icopt, id, res, jack, psol, h, wt, jskip, rpar, ipar, savr, delta, r, yic, ypic, pwk, wm, iwm, cj, uround, epli, sqrtn, rsqrtn, epcon, ratemx, stptol, jflg, icnflg, icnstr, iernls)
Definition ddasik.f:9
subroutine ddaspk(res, neq, t, y, yprime, tout, info, rtol, atol, idid, rwork, lrw, iwork, liw, rpar, ipar, jac, psol)
Definition ddaspk.f:7
subroutine ddatrp(x, xout, yout, ypout, neq, kold, phi, psi)
Definition ddatrp.f:2
subroutine ddawts(neq, iwt, rtol, atol, y, wt, rpar, ipar)
Definition ddawts.f:2
subroutine ddstp(x, y, yprime, neq, res, jac, psol, h, wt, vt, jstart, idid, rpar, ipar, phi, savr, delta, e, wm, iwm, alpha, beta, gamma, psi, sigma, cj, cjold, hold, s, hmin, uround, epli, sqrtn, rsqrtn, epcon, iphase, jcalc, jflg, k, kold, ns, nonneg, ntype, nls)
Definition ddstp.f:10
double precision function ddwnrm(neq, v, rwt, rpar, ipar)
Definition ddwnrm.f:6
subroutine dinvwt(neq, wt, ier)
Definition dinvwt.f:6
subroutine dnedd(x, y, yprime, neq, res, jacd, pdum, h, wt, jstart, idid, rpar, ipar, phi, gamma, dumsvr, delta, e, wm, iwm, cj, cjold, cjlast, s, uround, dume, dums, dumr, epcon, jcalc, jfdum, kp1, nonneg, ntype, iernls)
Definition dnedd.f:9
subroutine dnedk(x, y, yprime, neq, res, jack, psol, h, wt, jstart, idid, rpar, ipar, phi, gamma, savr, delta, e, wm, iwm, cj, cjold, cjlast, s, uround, epli, sqrtn, rsqrtn, epcon, jcalc, jflg, kp1, nonneg, ntype, iernls)
Definition dnedk.f:9
double lgamma(double x)
Definition lo-specfun.h:336
subroutine xerrwd(msg, nmes, nerr, level, ni, i1, i2, nr, r1, r2)
Definition xerrwd.f:4