Original WRF subgrid support version from John Michalakes without fire
[wrffire.git] / wrfv2_fire / chem / module_cmu_svode_solver.F
blob30f19daf9f02d40472b08fe25e3c5797e8c8717a
1 !**********************************************************************************  
2 ! This computer software was prepared by Battelle Memorial Institute, hereinafter
3 ! the Contractor, under Contract No. DE-AC05-76RL0 1830 with the Department of 
4 ! Energy (DOE). NEITHER THE GOVERNMENT NOR THE CONTRACTOR MAKES ANY WARRANTY,
5 ! EXPRESS OR IMPLIED, OR ASSUMES ANY LIABILITY FOR THE USE OF THIS SOFTWARE.
7 ! MOSAIC module: see module_mosaic_driver.F for information and terms of use
8 !**********************************************************************************  
10       module module_cmu_svode_solver
13       use module_peg_util, only:  peg_error_fatal
16       contains
19 !-----------------------------------------------------------------------
20 ! rce 2003-jul-01 - obtained from s.pandis
21 ! rce 2005-jan-21 - converted to a module.  Renamed svod01&2 common blocks 
22 !       to svode_cmn_01&2 to make their name more unique.
23 !       In xerrwv, changed msg to char*(nmes) to eliminate compiler warnings.
24 !-----------------------------------------------------------------------
26 !* ======================================================================
27 !* nist guide to available math software.
28 !* fullsource for module svode from package ode.
29 !* retrieved from netlib on tue jun 16 20:57:19 1998.
30 !* ======================================================================
32 !*deck svode
33       subroutine svode (f, neq, y, t, tout, itol, rtol, atol, itask,   &
34                   istate, iopt, rwork, lrw, iwork, liw, jac, mf,   &
35                   rpar, ipar)
36       external f, jac
37       real y, t, tout, rtol, atol, rwork, rpar
38       integer neq, itol, itask, istate, iopt, lrw, iwork, liw,   &
39               mf, ipar
40       dimension y(*), rtol(*), atol(*), rwork(lrw), iwork(liw),   &
41                 rpar(*), ipar(*)
42 !-----------------------------------------------------------------------
43 ! svode.. variable-coefficient ordinary differential equation solver,
44 ! with fixed-leading-coefficient implementation.
45 ! this version is in single precision.
47 ! svode solves the initial value problem for stiff or nonstiff
48 ! systems of first order odes,
49 !     dy/dt = f(t,y) ,  or, in component form,
50 !     dy(i)/dt = f(i) = f(i,t,y(1),y(2),...,y(neq)) (i = 1,...,neq).
51 ! svode is a package based on the episode and episodeb packages, and
52 ! on the odepack user interface standard, with minor modifications.
53 !-----------------------------------------------------------------------
54 ! revision history (yymmdd)
55 !   890615  date written
56 !   890922  added interrupt/restart ability, minor changes throughout.
57 !   910228  minor revisions in line format,  prologue, etc.
58 !   920227  modifications by d. pang:
59 !           (1) applied subgennam to get generic intrinsic names.
60 !           (2) changed intrinsic names to generic in comments.
61 !           (3) added *deck lines before each routine.
62 !   920721  names of routines and labeled common blocks changed, so as
63 !           to be unique in combined single/double precision code (ach).
64 !   920722  minor revisions to prologue (ach).
65 !   921106  fixed minor bug: etaq,etaqm1 in svstep save statement (ach).
66 !   921118  changed lunsav/mflgsv to ixsav (ach).
67 !   941222  removed mf overwrite; attached sign to h in estimated second
68 !           derivative in svhin; misc. comment corrections throughout.
69 !   970515  minor corrections to comments in prologue, svjac.
70 !-----------------------------------------------------------------------
71 ! references..
73 ! 1. p. n. brown, g. d. byrne, and a. c. hindmarsh, "vode: a variable
74 !    coefficient ode solver," siam j. sci. stat. comput., 10 (1989),
75 !    pp. 1038-1051.  also, llnl report ucrl-98412, june 1988.
76 ! 2. g. d. byrne and a. c. hindmarsh, "a polyalgorithm for the
77 !    numerical solution of ordinary differential equations,"
78 !    acm trans. math. software, 1 (1975), pp. 71-96.
79 ! 3. a. c. hindmarsh and g. d. byrne, "episode: an effective package
80 !    for the integration of systems of ordinary differential
81 !    equations," llnl report ucid-30112, rev. 1, april 1977.
82 ! 4. g. d. byrne and a. c. hindmarsh, "episodeb: an experimental
83 !    package for the integration of systems of ordinary differential
84 !    equations with banded jacobians," llnl report ucid-30132, april
85 !    1976.
86 ! 5. a. c. hindmarsh, "odepack, a systematized collection of ode
87 !    solvers," in scientific computing, r. s. stepleman et al., eds.,
88 !    north-holland, amsterdam, 1983, pp. 55-64.
89 ! 6. k. r. jackson and r. sacks-davis, "an alternative implementation
90 !    of variable step-size multistep formulas for stiff odes," acm
91 !    trans. math. software, 6 (1980), pp. 295-318.
92 !-----------------------------------------------------------------------
93 ! authors..
95 !               peter n. brown and alan c. hindmarsh
96 !               center for applied scientific computing, l-561
97 !               lawrence livermore national laboratory
98 !               livermore, ca 94551
99 ! and
100 !               george d. byrne
101 !               illinois institute of technology
102 !               chicago, il 60616
103 !-----------------------------------------------------------------------
104 ! summary of usage.
106 ! communication between the user and the svode package, for normal
107 ! situations, is summarized here.  this summary describes only a subset
108 ! of the full set of options available.  see the full description for
109 ! details, including optional communication, nonstandard options,
110 ! and instructions for special situations.  see also the example
111 ! problem (with program and output) following this summary.
113 ! a. first provide a subroutine of the form..
115 !           subroutine f (neq, t, y, ydot, rpar, ipar)
116 !           real t, y, ydot, rpar
117 !           dimension y(neq), ydot(neq)
119 ! which supplies the vector function f by loading ydot(i) with f(i).
121 ! b. next determine (or guess) whether or not the problem is stiff.
122 ! stiffness occurs when the jacobian matrix df/dy has an eigenvalue
123 ! whose real part is negative and large in magnitude, compared to the
124 ! reciprocal of the t span of interest.  if the problem is nonstiff,
125 ! use a method flag mf = 10.  if it is stiff, there are four standard
126 ! choices for mf (21, 22, 24, 25), and svode requires the jacobian
127 ! matrix in some form.  in these cases (mf .gt. 0), svode will use a
128 ! saved copy of the jacobian matrix.  if this is undesirable because of
129 ! storage limitations, set mf to the corresponding negative value
130 ! (-21, -22, -24, -25).  (see full description of mf below.)
131 ! the jacobian matrix is regarded either as full (mf = 21 or 22),
132 ! or banded (mf = 24 or 25).  in the banded case, svode requires two
133 ! half-bandwidth parameters ml and mu.  these are, respectively, the
134 ! widths of the lower and upper parts of the band, excluding the main
135 ! diagonal.  thus the band consists of the locations (i,j) with
136 ! i-ml .le. j .le. i+mu, and the full bandwidth is ml+mu+1.
138 ! c. if the problem is stiff, you are encouraged to supply the jacobian
139 ! directly (mf = 21 or 24), but if this is not feasible, svode will
140 ! compute it internally by difference quotients (mf = 22 or 25).
141 ! if you are supplying the jacobian, provide a subroutine of the form..
143 !           subroutine jac (neq, t, y, ml, mu, pd, nrowpd, rpar, ipar)
144 !           real t, y, pd, rpar
145 !           dimension y(neq), pd(nrowpd,neq)
147 ! which supplies df/dy by loading pd as follows..
148 !     for a full jacobian (mf = 21), load pd(i,j) with df(i)/dy(j),
149 ! the partial derivative of f(i) with respect to y(j).  (ignore the
150 ! ml and mu arguments in this case.)
151 !     for a banded jacobian (mf = 24), load pd(i-j+mu+1,j) with
152 ! df(i)/dy(j), i.e. load the diagonal lines of df/dy into the rows of
153 ! pd from the top down.
154 !     in either case, only nonzero elements need be loaded.
156 ! d. write a main program which calls subroutine svode once for
157 ! each point at which answers are desired.  this should also provide
158 ! for possible use of logical unit 6 for output of error messages
159 ! by svode.  on the first call to svode, supply arguments as follows..
160 ! f      = name of subroutine for right-hand side vector f.
161 !          this name must be declared external in calling program.
162 ! neq    = number of first order ode-s.
163 ! y      = array of initial values, of length neq.
164 ! t      = the initial value of the independent variable.
165 ! tout   = first point where output is desired (.ne. t).
166 ! itol   = 1 or 2 according as atol (below) is a scalar or array.
167 ! rtol   = relative tolerance parameter (scalar).
168 ! atol   = absolute tolerance parameter (scalar or array).
169 !          the estimated local error in y(i) will be controlled so as
170 !          to be roughly less (in magnitude) than
171 !             ewt(i) = rtol*abs(y(i)) + atol     if itol = 1, or
172 !             ewt(i) = rtol*abs(y(i)) + atol(i)  if itol = 2.
173 !          thus the local error test passes if, in each component,
174 !          either the absolute error is less than atol (or atol(i)),
175 !          or the relative error is less than rtol.
176 !          use rtol = 0.0 for pure absolute error control, and
177 !          use atol = 0.0 (or atol(i) = 0.0) for pure relative error
178 !          control.  caution.. actual (global) errors may exceed these
179 !          local tolerances, so choose them conservatively.
180 ! itask  = 1 for normal computation of output values of y at t = tout.
181 ! istate = integer flag (input and output).  set istate = 1.
182 ! iopt   = 0 to indicate no optional input used.
183 ! rwork  = real work array of length at least..
184 !             20 + 16*neq                      for mf = 10,
185 !             22 +  9*neq + 2*neq**2           for mf = 21 or 22,
186 !             22 + 11*neq + (3*ml + 2*mu)*neq  for mf = 24 or 25.
187 ! lrw    = declared length of rwork (in user's dimension statement).
188 ! iwork  = integer work array of length at least..
189 !             30        for mf = 10,
190 !             30 + neq  for mf = 21, 22, 24, or 25.
191 !          if mf = 24 or 25, input in iwork(1),iwork(2) the lower
192 !          and upper half-bandwidths ml,mu.
193 ! liw    = declared length of iwork (in user's dimension statement).
194 ! jac    = name of subroutine for jacobian matrix (mf = 21 or 24).
195 !          if used, this name must be declared external in calling
196 !          program.  if not used, pass a dummy name.
197 ! mf     = method flag.  standard values are..
198 !          10 for nonstiff (adams) method, no jacobian used.
199 !          21 for stiff (bdf) method, user-supplied full jacobian.
200 !          22 for stiff method, internally generated full jacobian.
201 !          24 for stiff method, user-supplied banded jacobian.
202 !          25 for stiff method, internally generated banded jacobian.
203 ! rpar,ipar = user-defined real and integer arrays passed to f and jac.
204 ! note that the main program must declare arrays y, rwork, iwork,
205 ! and possibly atol, rpar, and ipar.
207 ! e. the output from the first call (or any call) is..
208 !      y = array of computed values of y(t) vector.
209 !      t = corresponding value of independent variable (normally tout).
210 ! istate = 2  if svode was successful, negative otherwise.
211 !          -1 means excess work done on this call. (perhaps wrong mf.)
212 !          -2 means excess accuracy requested. (tolerances too small.)
213 !          -3 means illegal input detected. (see printed message.)
214 !          -4 means repeated error test failures. (check all input.)
215 !          -5 means repeated convergence failures. (perhaps bad
216 !             jacobian supplied or wrong choice of mf or tolerances.)
217 !          -6 means error weight became zero during problem. (solution
218 !             component i vanished, and atol or atol(i) = 0.)
220 ! f. to continue the integration after a successful return, simply
221 ! reset tout and call svode again.  no other parameters need be reset.
223 !-----------------------------------------------------------------------
224 ! example problem
226 ! the following is a simple example problem, with the coding
227 ! needed for its solution by svode.  the problem is from chemical
228 ! kinetics, and consists of the following three rate equations..
229 !     dy1/dt = -.04*y1 + 1.e4*y2*y3
230 !     dy2/dt = .04*y1 - 1.e4*y2*y3 - 3.e7*y2**2
231 !     dy3/dt = 3.e7*y2**2
232 ! on the interval from t = 0.0 to t = 4.e10, with initial conditions
233 ! y1 = 1.0, y2 = y3 = 0.  the problem is stiff.
235 ! the following coding solves this problem with svode, using mf = 21
236 ! and printing results at t = .4, 4., ..., 4.e10.  it uses
237 ! itol = 2 and atol much smaller for y2 than y1 or y3 because
238 ! y2 has much smaller values.
239 ! at the end of the run, statistical quantities of interest are
240 ! printed. (see optional output in the full description below.)
241 ! to generate fortran source code, replace c in column 1 with a blank
242 ! in the coding below.
244 !     external fex, jex
245 !     real atol, rpar, rtol, rwork, t, tout, y
246 !     dimension y(3), atol(3), rwork(67), iwork(33)
247 !     neq = 3
248 !     y(1) = 1.0e0
249 !     y(2) = 0.0e0
250 !     y(3) = 0.0e0
251 !     t = 0.0e0
252 !     tout = 0.4e0
253 !     itol = 2
254 !     rtol = 1.e-4
255 !     atol(1) = 1.e-8
256 !     atol(2) = 1.e-14
257 !     atol(3) = 1.e-6
258 !     itask = 1
259 !     istate = 1
260 !     iopt = 0
261 !     lrw = 67
262 !     liw = 33
263 !     mf = 21
264 !     do 40 iout = 1,12
265 !       call svode(fex,neq,y,t,tout,itol,rtol,atol,itask,istate,
266 !    1            iopt,rwork,lrw,iwork,liw,jex,mf,rpar,ipar)
267 !       write(6,20)t,y(1),y(2),y(3)
268 ! 20    format(' at t =',e12.4,'   y =',3e14.6)
269 !       if (istate .lt. 0) go to 80
270 ! 40    tout = tout*10.
271 !     write(6,60) iwork(11),iwork(12),iwork(13),iwork(19),
272 !    1            iwork(20),iwork(21),iwork(22)
273 ! 60  format(/' no. steps =',i4,'   no. f-s =',i4,
274 !    1       '   no. j-s =',i4,'   no. lu-s =',i4/
275 !    2       '  no. nonlinear iterations =',i4/
276 !    3       '  no. nonlinear convergence failures =',i4/
277 !    4       '  no. error test failures =',i4/)
278 !     stop
279 ! 80  write(6,90)istate
280 ! 90  format(///' error halt.. istate =',i3)
281 !     stop
282 !     end
284 !     subroutine fex (neq, t, y, ydot, rpar, ipar)
285 !     real rpar, t, y, ydot
286 !     dimension y(neq), ydot(neq)
287 !     ydot(1) = -.04e0*y(1) + 1.e4*y(2)*y(3)
288 !     ydot(3) = 3.e7*y(2)*y(2)
289 !     ydot(2) = -ydot(1) - ydot(3)
290 !     return
291 !     end
293 !     subroutine jex (neq, t, y, ml, mu, pd, nrpd, rpar, ipar)
294 !     real pd, rpar, t, y
295 !     dimension y(neq), pd(nrpd,neq)
296 !     pd(1,1) = -.04e0
297 !     pd(1,2) = 1.e4*y(3)
298 !     pd(1,3) = 1.e4*y(2)
299 !     pd(2,1) = .04e0
300 !     pd(2,3) = -pd(1,3)
301 !     pd(3,2) = 6.e7*y(2)
302 !     pd(2,2) = -pd(1,2) - pd(3,2)
303 !     return
304 !     end
306 ! the following output was obtained from the above program on a
307 ! cray-1 computer with the cft compiler.
309 ! at t =  4.0000e-01   y =  9.851680e-01  3.386314e-05  1.479817e-02
310 ! at t =  4.0000e+00   y =  9.055255e-01  2.240539e-05  9.445214e-02
311 ! at t =  4.0000e+01   y =  7.158108e-01  9.184883e-06  2.841800e-01
312 ! at t =  4.0000e+02   y =  4.505032e-01  3.222940e-06  5.494936e-01
313 ! at t =  4.0000e+03   y =  1.832053e-01  8.942690e-07  8.167938e-01
314 ! at t =  4.0000e+04   y =  3.898560e-02  1.621875e-07  9.610142e-01
315 ! at t =  4.0000e+05   y =  4.935882e-03  1.984013e-08  9.950641e-01
316 ! at t =  4.0000e+06   y =  5.166183e-04  2.067528e-09  9.994834e-01
317 ! at t =  4.0000e+07   y =  5.201214e-05  2.080593e-10  9.999480e-01
318 ! at t =  4.0000e+08   y =  5.213149e-06  2.085271e-11  9.999948e-01
319 ! at t =  4.0000e+09   y =  5.183495e-07  2.073399e-12  9.999995e-01
320 ! at t =  4.0000e+10   y =  5.450996e-08  2.180399e-13  9.999999e-01
322 ! no. steps = 595   no. f-s = 832   no. j-s =  13   no. lu-s = 112
323 !  no. nonlinear iterations = 831
324 !  no. nonlinear convergence failures =   0
325 !  no. error test failures =  22
326 !-----------------------------------------------------------------------
327 ! full description of user interface to svode.
329 ! the user interface to svode consists of the following parts.
331 ! i.   the call sequence to subroutine svode, which is a driver
332 !      routine for the solver.  this includes descriptions of both
333 !      the call sequence arguments and of user-supplied routines.
334 !      following these descriptions is
335 !        * a description of optional input available through the
336 !          call sequence,
337 !        * a description of optional output (in the work arrays), and
338 !        * instructions for interrupting and restarting a solution.
340 ! ii.  descriptions of other routines in the svode package that may be
341 !      (optionally) called by the user.  these provide the ability to
342 !      alter error message handling, save and restore the internal
343 !      common, and obtain specified derivatives of the solution y(t).
345 ! iii. descriptions of common blocks to be declared in overlay
346 !      or similar environments.
348 ! iv.  description of two routines in the svode package, either of
349 !      which the user may replace with his own version, if desired.
350 !      these relate to the measurement of errors.
352 !-----------------------------------------------------------------------
353 ! part i.  call sequence.
355 ! the call sequence parameters used for input only are
356 !     f, neq, tout, itol, rtol, atol, itask, iopt, lrw, liw, jac, mf,
357 ! and those used for both input and output are
358 !     y, t, istate.
359 ! the work arrays rwork and iwork are also used for conditional and
360 ! optional input and optional output.  (the term output here refers
361 ! to the return from subroutine svode to the user's calling program.)
363 ! the legality of input parameters will be thoroughly checked on the
364 ! initial call for the problem, but not checked thereafter unless a
365 ! change in input parameters is flagged by istate = 3 in the input.
367 ! the descriptions of the call arguments are as follows.
369 ! f      = the name of the user-supplied subroutine defining the
370 !          ode system.  the system must be put in the first-order
371 !          form dy/dt = f(t,y), where f is a vector-valued function
372 !          of the scalar t and the vector y.  subroutine f is to
373 !          compute the function f.  it is to have the form
374 !               subroutine f (neq, t, y, ydot, rpar, ipar)
375 !               real t, y, ydot, rpar
376 !               dimension y(neq), ydot(neq)
377 !          where neq, t, and y are input, and the array ydot = f(t,y)
378 !          is output.  y and ydot are arrays of length neq.
379 !          (in the dimension statement above, neq  can be replaced by
380 !          *  to make  y  and  ydot  assumed size arrays.)
381 !          subroutine f should not alter y(1),...,y(neq).
382 !          f must be declared external in the calling program.
384 !          subroutine f may access user-defined real and integer
385 !          work arrays rpar and ipar, which are to be dimensioned
386 !          in the main program.
388 !          if quantities computed in the f routine are needed
389 !          externally to svode, an extra call to f should be made
390 !          for this purpose, for consistent and accurate results.
391 !          if only the derivative dy/dt is needed, use svindy instead.
393 ! neq    = the size of the ode system (number of first order
394 !          ordinary differential equations).  used only for input.
395 !          neq may not be increased during the problem, but
396 !          can be decreased (with istate = 3 in the input).
398 ! y      = a real array for the vector of dependent variables, of
399 !          length neq or more.  used for both input and output on the
400 !          first call (istate = 1), and only for output on other calls.
401 !          on the first call, y must contain the vector of initial
402 !          values.  in the output, y contains the computed solution
403 !          evaluated at t.  if desired, the y array may be used
404 !          for other purposes between calls to the solver.
406 !          this array is passed as the y argument in all calls to
407 !          f and jac.
409 ! t      = the independent variable.  in the input, t is used only on
410 !          the first call, as the initial point of the integration.
411 !          in the output, after each call, t is the value at which a
412 !          computed solution y is evaluated (usually the same as tout).
413 !          on an error return, t is the farthest point reached.
415 ! tout   = the next value of t at which a computed solution is desired.
416 !          used only for input.
418 !          when starting the problem (istate = 1), tout may be equal
419 !          to t for one call, then should .ne. t for the next call.
420 !          for the initial t, an input value of tout .ne. t is used
421 !          in order to determine the direction of the integration
422 !          (i.e. the algebraic sign of the step sizes) and the rough
423 !          scale of the problem.  integration in either direction
424 !          (forward or backward in t) is permitted.
426 !          if itask = 2 or 5 (one-step modes), tout is ignored after
427 !          the first call (i.e. the first call with tout .ne. t).
428 !          otherwise, tout is required on every call.
430 !          if itask = 1, 3, or 4, the values of tout need not be
431 !          monotone, but a value of tout which backs up is limited
432 !          to the current internal t interval, whose endpoints are
433 !          tcur - hu and tcur.  (see optional output, below, for
434 !          tcur and hu.)
436 ! itol   = an indicator for the type of error control.  see
437 !          description below under atol.  used only for input.
439 ! rtol   = a relative error tolerance parameter, either a scalar or
440 !          an array of length neq.  see description below under atol.
441 !          input only.
443 ! atol   = an absolute error tolerance parameter, either a scalar or
444 !          an array of length neq.  input only.
446 !          the input parameters itol, rtol, and atol determine
447 !          the error control performed by the solver.  the solver will
448 !          control the vector e = (e(i)) of estimated local errors
449 !          in y, according to an inequality of the form
450 !                      rms-norm of ( e(i)/ewt(i) )   .le.   1,
451 !          where       ewt(i) = rtol(i)*abs(y(i)) + atol(i),
452 !          and the rms-norm (root-mean-square norm) here is
453 !          rms-norm(v) = sqrt(sum v(i)**2 / neq).  here ewt = (ewt(i))
454 !          is a vector of weights which must always be positive, and
455 !          the values of rtol and atol should all be non-negative.
456 !          the following table gives the types (scalar/array) of
457 !          rtol and atol, and the corresponding form of ewt(i).
459 !             itol    rtol       atol          ewt(i)
460 !              1     scalar     scalar     rtol*abs(y(i)) + atol
461 !              2     scalar     array      rtol*abs(y(i)) + atol(i)
462 !              3     array      scalar     rtol(i)*abs(y(i)) + atol
463 !              4     array      array      rtol(i)*abs(y(i)) + atol(i)
465 !          when either of these parameters is a scalar, it need not
466 !          be dimensioned in the user's calling program.
468 !          if none of the above choices (with itol, rtol, and atol
469 !          fixed throughout the problem) is suitable, more general
470 !          error controls can be obtained by substituting
471 !          user-supplied routines for the setting of ewt and/or for
472 !          the norm calculation.  see part iv below.
474 !          if global errors are to be estimated by making a repeated
475 !          run on the same problem with smaller tolerances, then all
476 !          components of rtol and atol (i.e. of ewt) should be scaled
477 !          down uniformly.
479 ! itask  = an index specifying the task to be performed.
480 !          input only.  itask has the following values and meanings.
481 !          1  means normal computation of output values of y(t) at
482 !             t = tout (by overshooting and interpolating).
483 !          2  means take one step only and return.
484 !          3  means stop at the first internal mesh point at or
485 !             beyond t = tout and return.
486 !          4  means normal computation of output values of y(t) at
487 !             t = tout but without overshooting t = tcrit.
488 !             tcrit must be input as rwork(1).  tcrit may be equal to
489 !             or beyond tout, but not behind it in the direction of
490 !             integration.  this option is useful if the problem
491 !             has a singularity at or beyond t = tcrit.
492 !          5  means take one step, without passing tcrit, and return.
493 !             tcrit must be input as rwork(1).
495 !          note..  if itask = 4 or 5 and the solver reaches tcrit
496 !          (within roundoff), it will return t = tcrit (exactly) to
497 !          indicate this (unless itask = 4 and tout comes before tcrit,
498 !          in which case answers at t = tout are returned first).
500 ! istate = an index used for input and output to specify the
501 !          the state of the calculation.
503 !          in the input, the values of istate are as follows.
504 !          1  means this is the first call for the problem
505 !             (initializations will be done).  see note below.
506 !          2  means this is not the first call, and the calculation
507 !             is to continue normally, with no change in any input
508 !             parameters except possibly tout and itask.
509 !             (if itol, rtol, and/or atol are changed between calls
510 !             with istate = 2, the new values will be used but not
511 !             tested for legality.)
512 !          3  means this is not the first call, and the
513 !             calculation is to continue normally, but with
514 !             a change in input parameters other than
515 !             tout and itask.  changes are allowed in
516 !             neq, itol, rtol, atol, iopt, lrw, liw, mf, ml, mu,
517 !             and any of the optional input except h0.
518 !             (see iwork description for ml and mu.)
519 !          note..  a preliminary call with tout = t is not counted
520 !          as a first call here, as no initialization or checking of
521 !          input is done.  (such a call is sometimes useful to include
522 !          the initial conditions in the output.)
523 !          thus the first call for which tout .ne. t requires
524 !          istate = 1 in the input.
526 !          in the output, istate has the following values and meanings.
527 !           1  means nothing was done, as tout was equal to t with
528 !              istate = 1 in the input.
529 !           2  means the integration was performed successfully.
530 !          -1  means an excessive amount of work (more than mxstep
531 !              steps) was done on this call, before completing the
532 !              requested task, but the integration was otherwise
533 !              successful as far as t.  (mxstep is an optional input
534 !              and is normally 500.)  to continue, the user may
535 !              simply reset istate to a value .gt. 1 and call again.
536 !              (the excess work step counter will be reset to 0.)
537 !              in addition, the user may increase mxstep to avoid
538 !              this error return.  (see optional input below.)
539 !          -2  means too much accuracy was requested for the precision
540 !              of the machine being used.  this was detected before
541 !              completing the requested task, but the integration
542 !              was successful as far as t.  to continue, the tolerance
543 !              parameters must be reset, and istate must be set
544 !              to 3.  the optional output tolsf may be used for this
545 !              purpose.  (note.. if this condition is detected before
546 !              taking any steps, then an illegal input return
547 !              (istate = -3) occurs instead.)
548 !          -3  means illegal input was detected, before taking any
549 !              integration steps.  see written message for details.
550 !              note..  if the solver detects an infinite loop of calls
551 !              to the solver with illegal input, it will cause
552 !              the run to stop.
553 !          -4  means there were repeated error test failures on
554 !              one attempted step, before completing the requested
555 !              task, but the integration was successful as far as t.
556 !              the problem may have a singularity, or the input
557 !              may be inappropriate.
558 !          -5  means there were repeated convergence test failures on
559 !              one attempted step, before completing the requested
560 !              task, but the integration was successful as far as t.
561 !              this may be caused by an inaccurate jacobian matrix,
562 !              if one is being used.
563 !          -6  means ewt(i) became zero for some i during the
564 !              integration.  pure relative error control (atol(i)=0.0)
565 !              was requested on a variable which has now vanished.
566 !              the integration was successful as far as t.
568 !          note..  since the normal output value of istate is 2,
569 !          it does not need to be reset for normal continuation.
570 !          also, since a negative input value of istate will be
571 !          regarded as illegal, a negative output value requires the
572 !          user to change it, and possibly other input, before
573 !          calling the solver again.
575 ! iopt   = an integer flag to specify whether or not any optional
576 !          input is being used on this call.  input only.
577 !          the optional input is listed separately below.
578 !          iopt = 0 means no optional input is being used.
579 !                   default values will be used in all cases.
580 !          iopt = 1 means optional input is being used.
582 ! rwork  = a real working array (single precision).
583 !          the length of rwork must be at least
584 !             20 + nyh*(maxord + 1) + 3*neq + lwm    where
585 !          nyh    = the initial value of neq,
586 !          maxord = 12 (if meth = 1) or 5 (if meth = 2) (unless a
587 !                   smaller value is given as an optional input),
588 !          lwm = length of work space for matrix-related data..
589 !          lwm = 0             if miter = 0,
590 !          lwm = 2*neq**2 + 2  if miter = 1 or 2, and mf.gt.0,
591 !          lwm = neq**2 + 2    if miter = 1 or 2, and mf.lt.0,
592 !          lwm = neq + 2       if miter = 3,
593 !          lwm = (3*ml+2*mu+2)*neq + 2 if miter = 4 or 5, and mf.gt.0,
594 !          lwm = (2*ml+mu+1)*neq + 2   if miter = 4 or 5, and mf.lt.0.
595 !          (see the mf description for meth and miter.)
596 !          thus if maxord has its default value and neq is constant,
597 !          this length is..
598 !             20 + 16*neq                    for mf = 10,
599 !             22 + 16*neq + 2*neq**2         for mf = 11 or 12,
600 !             22 + 16*neq + neq**2           for mf = -11 or -12,
601 !             22 + 17*neq                    for mf = 13,
602 !             22 + 18*neq + (3*ml+2*mu)*neq  for mf = 14 or 15,
603 !             22 + 17*neq + (2*ml+mu)*neq    for mf = -14 or -15,
604 !             20 +  9*neq                    for mf = 20,
605 !             22 +  9*neq + 2*neq**2         for mf = 21 or 22,
606 !             22 +  9*neq + neq**2           for mf = -21 or -22,
607 !             22 + 10*neq                    for mf = 23,
608 !             22 + 11*neq + (3*ml+2*mu)*neq  for mf = 24 or 25.
609 !             22 + 10*neq + (2*ml+mu)*neq    for mf = -24 or -25.
610 !          the first 20 words of rwork are reserved for conditional
611 !          and optional input and optional output.
613 !          the following word in rwork is a conditional input..
614 !            rwork(1) = tcrit = critical value of t which the solver
615 !                       is not to overshoot.  required if itask is
616 !                       4 or 5, and ignored otherwise.  (see itask.)
618 ! lrw    = the length of the array rwork, as declared by the user.
619 !          (this will be checked by the solver.)
621 ! iwork  = an integer work array.  the length of iwork must be at least
622 !             30        if miter = 0 or 3 (mf = 10, 13, 20, 23), or
623 !             30 + neq  otherwise (abs(mf) = 11,12,14,15,21,22,24,25).
624 !          the first 30 words of iwork are reserved for conditional and
625 !          optional input and optional output.
627 !          the following 2 words in iwork are conditional input..
628 !            iwork(1) = ml     these are the lower and upper
629 !            iwork(2) = mu     half-bandwidths, respectively, of the
630 !                       banded jacobian, excluding the main diagonal.
631 !                       the band is defined by the matrix locations
632 !                       (i,j) with i-ml .le. j .le. i+mu.  ml and mu
633 !                       must satisfy  0 .le.  ml,mu  .le. neq-1.
634 !                       these are required if miter is 4 or 5, and
635 !                       ignored otherwise.  ml and mu may in fact be
636 !                       the band parameters for a matrix to which
637 !                       df/dy is only approximately equal.
639 ! liw    = the length of the array iwork, as declared by the user.
640 !          (this will be checked by the solver.)
642 ! note..  the work arrays must not be altered between calls to svode
643 ! for the same problem, except possibly for the conditional and
644 ! optional input, and except for the last 3*neq words of rwork.
645 ! the latter space is used for internal scratch space, and so is
646 ! available for use by the user outside svode between calls, if
647 ! desired (but not for use by f or jac).
649 ! jac    = the name of the user-supplied routine (miter = 1 or 4) to
650 !          compute the jacobian matrix, df/dy, as a function of
651 !          the scalar t and the vector y.  it is to have the form
652 !               subroutine jac (neq, t, y, ml, mu, pd, nrowpd,
653 !                               rpar, ipar)
654 !               real t, y, pd, rpar
655 !               dimension y(neq), pd(nrowpd, neq)
656 !          where neq, t, y, ml, mu, and nrowpd are input and the array
657 !          pd is to be loaded with partial derivatives (elements of the
658 !          jacobian matrix) in the output.  pd must be given a first
659 !          dimension of nrowpd.  t and y have the same meaning as in
660 !          subroutine f.  (in the dimension statement above, neq can
661 !          be replaced by  *  to make y and pd assumed size arrays.)
662 !               in the full matrix case (miter = 1), ml and mu are
663 !          ignored, and the jacobian is to be loaded into pd in
664 !          columnwise manner, with df(i)/dy(j) loaded into pd(i,j).
665 !               in the band matrix case (miter = 4), the elements
666 !          within the band are to be loaded into pd in columnwise
667 !          manner, with diagonal lines of df/dy loaded into the rows
668 !          of pd. thus df(i)/dy(j) is to be loaded into pd(i-j+mu+1,j).
669 !          ml and mu are the half-bandwidth parameters. (see iwork).
670 !          the locations in pd in the two triangular areas which
671 !          correspond to nonexistent matrix elements can be ignored
672 !          or loaded arbitrarily, as they are overwritten by svode.
673 !               jac need not provide df/dy exactly.  a crude
674 !          approximation (possibly with a smaller bandwidth) will do.
675 !               in either case, pd is preset to zero by the solver,
676 !          so that only the nonzero elements need be loaded by jac.
677 !          each call to jac is preceded by a call to f with the same
678 !          arguments neq, t, and y.  thus to gain some efficiency,
679 !          intermediate quantities shared by both calculations may be
680 !          saved in a user common block by f and not recomputed by jac,
681 !          if desired.  also, jac may alter the y array, if desired.
682 !          jac must be declared external in the calling program.
683 !               subroutine jac may access user-defined real and integer
684 !          work arrays, rpar and ipar, whose dimensions are set by the
685 !          user in the main program.
687 ! mf     = the method flag.  used only for input.  the legal values of
688 !          mf are 10, 11, 12, 13, 14, 15, 20, 21, 22, 23, 24, 25,
689 !          -11, -12, -14, -15, -21, -22, -24, -25.
690 !          mf is a signed two-digit integer, mf = jsv*(10*meth + miter).
691 !          jsv = sign(mf) indicates the jacobian-saving strategy..
692 !            jsv =  1 means a copy of the jacobian is saved for reuse
693 !                     in the corrector iteration algorithm.
694 !            jsv = -1 means a copy of the jacobian is not saved
695 !                     (valid only for miter = 1, 2, 4, or 5).
696 !          meth indicates the basic linear multistep method..
697 !            meth = 1 means the implicit adams method.
698 !            meth = 2 means the method based on backward
699 !                     differentiation formulas (bdf-s).
700 !          miter indicates the corrector iteration method..
701 !            miter = 0 means functional iteration (no jacobian matrix
702 !                      is involved).
703 !            miter = 1 means chord iteration with a user-supplied
704 !                      full (neq by neq) jacobian.
705 !            miter = 2 means chord iteration with an internally
706 !                      generated (difference quotient) full jacobian
707 !                      (using neq extra calls to f per df/dy value).
708 !            miter = 3 means chord iteration with an internally
709 !                      generated diagonal jacobian approximation
710 !                      (using 1 extra call to f per df/dy evaluation).
711 !            miter = 4 means chord iteration with a user-supplied
712 !                      banded jacobian.
713 !            miter = 5 means chord iteration with an internally
714 !                      generated banded jacobian (using ml+mu+1 extra
715 !                      calls to f per df/dy evaluation).
716 !          if miter = 1 or 4, the user must supply a subroutine jac
717 !          (the name is arbitrary) as described above under jac.
718 !          for other values of miter, a dummy argument can be used.
720 ! rpar     user-specified array used to communicate real parameters
721 !          to user-supplied subroutines.  if rpar is a vector, then
722 !          it must be dimensioned in the user's main program.  if it
723 !          is unused or it is a scalar, then it need not be
724 !          dimensioned.
726 ! ipar     user-specified array used to communicate integer parameter
727 !          to user-supplied subroutines.  the comments on dimensioning
728 !          rpar apply to ipar.
729 !-----------------------------------------------------------------------
730 ! optional input.
732 ! the following is a list of the optional input provided for in the
733 ! call sequence.  (see also part ii.)  for each such input variable,
734 ! this table lists its name as used in this documentation, its
735 ! location in the call sequence, its meaning, and the default value.
736 ! the use of any of this input requires iopt = 1, and in that
737 ! case all of this input is examined.  a value of zero for any
738 ! of these optional input variables will cause the default value to be
739 ! used.  thus to use a subset of the optional input, simply preload
740 ! locations 5 to 10 in rwork and iwork to 0.0 and 0 respectively, and
741 ! then set those of interest to nonzero values.
743 ! name    location      meaning and default value
745 ! h0      rwork(5)  the step size to be attempted on the first step.
746 !                   the default value is determined by the solver.
748 ! hmax    rwork(6)  the maximum absolute step size allowed.
749 !                   the default value is infinite.
751 ! hmin    rwork(7)  the minimum absolute step size allowed.
752 !                   the default value is 0.  (this lower bound is not
753 !                   enforced on the final step before reaching tcrit
754 !                   when itask = 4 or 5.)
756 ! maxord  iwork(5)  the maximum order to be allowed.  the default
757 !                   value is 12 if meth = 1, and 5 if meth = 2.
758 !                   if maxord exceeds the default value, it will
759 !                   be reduced to the default value.
760 !                   if maxord is changed during the problem, it may
761 !                   cause the current order to be reduced.
763 ! mxstep  iwork(6)  maximum number of (internally defined) steps
764 !                   allowed during one call to the solver.
765 !                   the default value is 500.
767 ! mxhnil  iwork(7)  maximum number of messages printed (per problem)
768 !                   warning that t + h = t on a step (h = step size).
769 !                   this must be positive to result in a non-default
770 !                   value.  the default value is 10.
772 !-----------------------------------------------------------------------
773 ! optional output.
775 ! as optional additional output from svode, the variables listed
776 ! below are quantities related to the performance of svode
777 ! which are available to the user.  these are communicated by way of
778 ! the work arrays, but also have internal mnemonic names as shown.
779 ! except where stated otherwise, all of this output is defined
780 ! on any successful return from svode, and on any return with
781 ! istate = -1, -2, -4, -5, or -6.  on an illegal input return
782 ! (istate = -3), they will be unchanged from their existing values
783 ! (if any), except possibly for tolsf, lenrw, and leniw.
784 ! on any error return, output relevant to the error will be defined,
785 ! as noted below.
787 ! name    location      meaning
789 ! hu      rwork(11) the step size in t last used (successfully).
791 ! hcur    rwork(12) the step size to be attempted on the next step.
793 ! tcur    rwork(13) the current value of the independent variable
794 !                   which the solver has actually reached, i.e. the
795 !                   current internal mesh point in t.  in the output,
796 !                   tcur will always be at least as far from the
797 !                   initial value of t as the current argument t,
798 !                   but may be farther (if interpolation was done).
800 ! tolsf   rwork(14) a tolerance scale factor, greater than 1.0,
801 !                   computed when a request for too much accuracy was
802 !                   detected (istate = -3 if detected at the start of
803 !                   the problem, istate = -2 otherwise).  if itol is
804 !                   left unaltered but rtol and atol are uniformly
805 !                   scaled up by a factor of tolsf for the next call,
806 !                   then the solver is deemed likely to succeed.
807 !                   (the user may also ignore tolsf and alter the
808 !                   tolerance parameters in any other way appropriate.)
810 ! nst     iwork(11) the number of steps taken for the problem so far.
812 ! nfe     iwork(12) the number of f evaluations for the problem so far.
814 ! nje     iwork(13) the number of jacobian evaluations so far.
816 ! nqu     iwork(14) the method order last used (successfully).
818 ! nqcur   iwork(15) the order to be attempted on the next step.
820 ! imxer   iwork(16) the index of the component of largest magnitude in
821 !                   the weighted local error vector ( e(i)/ewt(i) ),
822 !                   on an error return with istate = -4 or -5.
824 ! lenrw   iwork(17) the length of rwork actually required.
825 !                   this is defined on normal returns and on an illegal
826 !                   input return for insufficient storage.
828 ! leniw   iwork(18) the length of iwork actually required.
829 !                   this is defined on normal returns and on an illegal
830 !                   input return for insufficient storage.
832 ! nlu     iwork(19) the number of matrix lu decompositions so far.
834 ! nni     iwork(20) the number of nonlinear (newton) iterations so far.
836 ! ncfn    iwork(21) the number of convergence failures of the nonlinear
837 !                   solver so far.
839 ! netf    iwork(22) the number of error test failures of the integrator
840 !                   so far.
842 ! the following two arrays are segments of the rwork array which
843 ! may also be of interest to the user as optional output.
844 ! for each array, the table below gives its internal name,
845 ! its base address in rwork, and its description.
847 ! name    base address      description
849 ! yh      21             the nordsieck history array, of size nyh by
850 !                        (nqcur + 1), where nyh is the initial value
851 !                        of neq.  for j = 0,1,...,nqcur, column j+1
852 !                        of yh contains hcur**j/factorial(j) times
853 !                        the j-th derivative of the interpolating
854 !                        polynomial currently representing the
855 !                        solution, evaluated at t = tcur.
857 ! acor     lenrw-neq+1   array of size neq used for the accumulated
858 !                        corrections on each step, scaled in the output
859 !                        to represent the estimated local error in y
860 !                        on the last step.  this is the vector e in
861 !                        the description of the error control.  it is
862 !                        defined only on a successful return from svode.
864 !-----------------------------------------------------------------------
865 ! interrupting and restarting
867 ! if the integration of a given problem by svode is to be
868 ! interrrupted and then later continued, such as when restarting
869 ! an interrupted run or alternating between two or more ode problems,
870 ! the user should save, following the return from the last svode call
871 ! prior to the interruption, the contents of the call sequence
872 ! variables and internal common blocks, and later restore these
873 ! values before the next svode call for that problem.  to save
874 ! and restore the common blocks, use subroutine svsrco, as
875 ! described below in part ii.
877 ! in addition, if non-default values for either lun or mflag are
878 ! desired, an extra call to xsetun and/or xsetf should be made just
879 ! before continuing the integration.  see part ii below for details.
881 !-----------------------------------------------------------------------
882 ! part ii.  other routines callable.
884 ! the following are optional calls which the user may make to
885 ! gain additional capabilities in conjunction with svode.
886 ! (the routines xsetun and xsetf are designed to conform to the
887 ! slatec error handling package.)
889 !     form of call                  function
890 !  call xsetun(lun)           set the logical unit number, lun, for
891 !                             output of messages from svode, if
892 !                             the default is not desired.
893 !                             the default value of lun is 6.
895 !  call xsetf(mflag)          set a flag to control the printing of
896 !                             messages by svode.
897 !                             mflag = 0 means do not print. (danger..
898 !                             this risks losing valuable information.)
899 !                             mflag = 1 means print (the default).
901 !                             either of the above calls may be made at
902 !                             any time and will take effect immediately.
904 !  call svsrco(rsav,isav,job) saves and restores the contents of
905 !                             the internal common blocks used by
906 !                             svode. (see part iii below.)
907 !                             rsav must be a real array of length 49
908 !                             or more, and isav must be an integer
909 !                             array of length 40 or more.
910 !                             job=1 means save common into rsav/isav.
911 !                             job=2 means restore common from rsav/isav.
912 !                                svsrco is useful if one is
913 !                             interrupting a run and restarting
914 !                             later, or alternating between two or
915 !                             more problems solved with svode.
917 !  call svindy(,,,,,)         provide derivatives of y, of various
918 !        (see below.)         orders, at a specified point t, if
919 !                             desired.  it may be called only after
920 !                             a successful return from svode.
922 ! the detailed instructions for using svindy are as follows.
923 ! the form of the call is..
925 !  call svindy (t, k, rwork(21), nyh, dky, iflag)
927 ! the input parameters are..
929 ! t         = value of independent variable where answers are desired
930 !             (normally the same as the t last returned by svode).
931 !             for valid results, t must lie between tcur - hu and tcur.
932 !             (see optional output for tcur and hu.)
933 ! k         = integer order of the derivative desired.  k must satisfy
934 !             0 .le. k .le. nqcur, where nqcur is the current order
935 !             (see optional output).  the capability corresponding
936 !             to k = 0, i.e. computing y(t), is already provided
937 !             by svode directly.  since nqcur .ge. 1, the first
938 !             derivative dy/dt is always available with svindy.
939 ! rwork(21) = the base address of the history array yh.
940 ! nyh       = column length of yh, equal to the initial value of neq.
942 ! the output parameters are..
944 ! dky       = a real array of length neq containing the computed value
945 !             of the k-th derivative of y(t).
946 ! iflag     = integer flag, returned as 0 if k and t were legal,
947 !             -1 if k was illegal, and -2 if t was illegal.
948 !             on an error return, a message is also written.
949 !-----------------------------------------------------------------------
950 ! part iii.  common blocks.
951 ! if svode is to be used in an overlay situation, the user
952 ! must declare, in the primary overlay, the variables in..
953 !   (1) the call sequence to svode,
954 !   (2) the two internal common blocks
955 !         /svode_cmn_01/  of length  81  (48 single precision words
956 !                         followed by 33 integer words),
957 !         /svode_cmn_02/  of length  9  (1 single precision word
958 !                         followed by 8 integer words),
960 ! if svode is used on a system in which the contents of internal
961 ! common blocks are not preserved between calls, the user should
962 ! declare the above two common blocks in his main program to insure
963 ! that their contents are preserved.
965 !-----------------------------------------------------------------------
966 ! part iv.  optionally replaceable solver routines.
968 ! below are descriptions of two routines in the svode package which
969 ! relate to the measurement of errors.  either routine can be
970 ! replaced by a user-supplied version, if desired.  however, since such
971 ! a replacement may have a major impact on performance, it should be
972 ! done only when absolutely necessary, and only with great caution.
973 ! (note.. the means by which the package version of a routine is
974 ! superseded by the user's version may be system-dependent.)
976 ! (a) sewset.
977 ! the following subroutine is called just before each internal
978 ! integration step, and sets the array of error weights, ewt, as
979 ! described under itol/rtol/atol above..
980 !     subroutine sewset (neq, itol, rtol, atol, ycur, ewt)
981 ! where neq, itol, rtol, and atol are as in the svode call sequence,
982 ! ycur contains the current dependent variable vector, and
983 ! ewt is the array of weights set by sewset.
985 ! if the user supplies this subroutine, it must return in ewt(i)
986 ! (i = 1,...,neq) a positive quantity suitable for comparison with
987 ! errors in y(i).  the ewt array returned by sewset is passed to the
988 ! svnorm routine (see below.), and also used by svode in the computation
989 ! of the optional output imxer, the diagonal jacobian approximation,
990 ! and the increments for difference quotient jacobians.
992 ! in the user-supplied version of sewset, it may be desirable to use
993 ! the current values of derivatives of y.  derivatives up to order nq
994 ! are available from the history array yh, described above under
995 ! optional output.  in sewset, yh is identical to the ycur array,
996 ! extended to nq + 1 columns with a column length of nyh and scale
997 ! factors of h**j/factorial(j).  on the first call for the problem,
998 ! given by nst = 0, nq is 1 and h is temporarily set to 1.0.
999 ! nyh is the initial value of neq.  the quantities nq, h, and nst
1000 ! can be obtained by including in sewset the statements..
1001 !     real rvod, h, hu
1002 !     common /svode_cmn_01/ rvod(48), ivod(33)
1003 !     common /svode_cmn_02/ hu, ncfn, netf, nfe, nje, nlu, nni, nqu, nst
1004 !     nq = ivod(28)
1005 !     h = rvod(21)
1006 ! thus, for example, the current value of dy/dt can be obtained as
1007 ! ycur(nyh+i)/h  (i=1,...,neq)  (and the division by h is
1008 ! unnecessary when nst = 0).
1010 ! (b) svnorm.
1011 ! the following is a real function routine which computes the weighted
1012 ! root-mean-square norm of a vector v..
1013 !     d = svnorm (n, v, w)
1014 ! where..
1015 !   n = the length of the vector,
1016 !   v = real array of length n containing the vector,
1017 !   w = real array of length n containing weights,
1018 !   d = sqrt( (1/n) * sum(v(i)*w(i))**2 ).
1019 ! svnorm is called with n = neq and with w(i) = 1.0/ewt(i), where
1020 ! ewt is as set by subroutine sewset.
1022 ! if the user supplies this function, it should return a non-negative
1023 ! value of svnorm suitable for use in the error control in svode.
1024 ! none of the arguments should be altered by svnorm.
1025 ! for example, a user-supplied svnorm routine might..
1026 !   -substitute a max-norm of (v(i)*w(i)) for the rms-norm, or
1027 !   -ignore some components of v in the norm, with the effect of
1028 !    suppressing the error control on those components of y.
1029 !-----------------------------------------------------------------------
1030 ! other routines in the svode package.
1032 ! in addition to subroutine svode, the svode package includes the
1033 ! following subroutines and function routines..
1034 !  svhin     computes an approximate step size for the initial step.
1035 !  svindy    computes an interpolated value of the y vector at t = tout.
1036 !  svstep    is the core integrator, which does one step of the
1037 !            integration and the associated error control.
1038 !  svset     sets all method coefficients and test constants.
1039 !  svnlsd    solves the underlying nonlinear system -- the corrector.
1040 !  svjac     computes and preprocesses the jacobian matrix j = df/dy
1041 !            and the newton iteration matrix p = i - (h/l1)*j.
1042 !  svsol     manages solution of linear system in chord iteration.
1043 !  svjust    adjusts the history array on a change of order.
1044 !  sewset    sets the error weight vector ewt before each step.
1045 !  svnorm    computes the weighted r.m.s. norm of a vector.
1046 !  svsrco    is a user-callable routine to save and restore
1047 !            the contents of the internal common blocks.
1048 !  sacopy    is a routine to copy one two-dimensional array to another.
1049 !  sgefa and sgesl   are routines from linpack for solving full
1050 !            systems of linear algebraic equations.
1051 !  sgbfa and sgbsl   are routines from linpack for solving banded
1052 !            linear systems.
1053 !  saxpy, sscal, and scopy are basic linear algebra modules (blas).
1054 !  r1mach    sets the unit roundoff of the machine.
1055 !  xerrwv, xsetun, xsetf, and ixsav handle the printing of all
1056 !            error messages and warnings.  xerrwv is machine-dependent.
1057 ! note..  svnorm, r1mach, and ixsav are function routines.
1058 ! all the others are subroutines.
1060 ! the intrinsic and external routines used by the svode package are..
1061 ! abs, max, min, real, sign, sqrt, and write.
1063 !-----------------------------------------------------------------------
1065 ! type declarations for labeled common block svode_cmn_01 --------------------
1067       real acnrm, ccmxj, conp, crate, drc, el,   &
1068            eta, etamax, h, hmin, hmxi, hnew, hscal, prl1,   &
1069            rc, rl1, tau, tq, tn, uround
1070       integer icf, init, ipup, jcur, jstart, jsv, kflag, kuth,   &
1071               l, lmax, lyh, lewt, lacor, lsavf, lwm, liwm,   &
1072               locjs, maxord, meth, miter, msbj, mxhnil, mxstep,   &
1073               n, newh, newq, nhnil, nq, nqnyh, nqwait, nslj,   &
1074               nslp, nyh
1076 ! type declarations for labeled common block svode_cmn_02 --------------------
1078       real hu
1079       integer ncfn, netf, nfe, nje, nlu, nni, nqu, nst
1081 ! type declarations for local variables --------------------------------
1083 !     external svnlsd      ! rce 2005-jan-21 - module conversion
1084       logical ihit
1085       real atoli, big, ewti, four, h0, hmax, hmx, hun, one,   &
1086          pt2, rh, rtoli, size, tcrit, tnext, tolsf, tp, two, zero
1087       integer i, ier, iflag, imxer, jco, kgo, leniw, lenj, lenp, lenrw,   &
1088          lenwm, lf0, mband, mfa, ml, mord, mu, mxhnl0, mxstp0, niter,   &
1089          nslast
1090       character*80 msg
1092 ! type declaration for function subroutines called ---------------------
1094 !     real r1mach, svnorm      ! rce 2005-jan-21 - module conversion
1096       dimension mord(2)
1097 !-----------------------------------------------------------------------
1098 ! the following fortran-77 declaration is to cause the values of the
1099 ! listed (local) variables to be saved between calls to svode.
1100 !-----------------------------------------------------------------------
1101       save mord, mxhnl0, mxstp0
1102       save zero, one, two, four, pt2, hun
1103 !-----------------------------------------------------------------------
1104 ! the following internal common blocks contain variables which are
1105 ! communicated between subroutines in the svode package, or which are
1106 ! to be saved between calls to svode.
1107 ! in each block, real variables precede integers.
1108 ! the block /svode_cmn_01/ appears in subroutines svode, svindy, svstep,
1109 ! svset, svnlsd, svjac, svsol, svjust and svsrco.
1110 ! the block /svode_cmn_02/ appears in subroutines svode, svindy, svstep,
1111 ! svnlsd, svjac, and svsrco.
1113 ! the variables stored in the internal common blocks are as follows..
1115 ! acnrm  = weighted r.m.s. norm of accumulated correction vectors.
1116 ! ccmxj  = threshhold on drc for updating the jacobian. (see drc.)
1117 ! conp   = the saved value of tq(5).
1118 ! crate  = estimated corrector convergence rate constant.
1119 ! drc    = relative change in h*rl1 since last svjac call.
1120 ! el     = real array of integration coefficients.  see svset.
1121 ! eta    = saved tentative ratio of new to old h.
1122 ! etamax = saved maximum value of eta to be allowed.
1123 ! h      = the step size.
1124 ! hmin   = the minimum absolute value of the step size h to be used.
1125 ! hmxi   = inverse of the maximum absolute value of h to be used.
1126 !          hmxi = 0.0 is allowed and corresponds to an infinite hmax.
1127 ! hnew   = the step size to be attempted on the next step.
1128 ! hscal  = stepsize in scaling of yh array.
1129 ! prl1   = the saved value of rl1.
1130 ! rc     = ratio of current h*rl1 to value on last svjac call.
1131 ! rl1    = the reciprocal of the coefficient el(1).
1132 ! tau    = real vector of past nq step sizes, length 13.
1133 ! tq     = a real vector of length 5 in which svset stores constants
1134 !          used for the convergence test, the error test, and the
1135 !          selection of h at a new order.
1136 ! tn     = the independent variable, updated on each step taken.
1137 ! uround = the machine unit roundoff.  the smallest positive real number
1138 !          such that  1.0 + uround .ne. 1.0
1139 ! icf    = integer flag for convergence failure in svnlsd..
1140 !            0 means no failures.
1141 !            1 means convergence failure with out of date jacobian
1142 !                   (recoverable error).
1143 !            2 means convergence failure with current jacobian or
1144 !                   singular matrix (unrecoverable error).
1145 ! init   = saved integer flag indicating whether initialization of the
1146 !          problem has been done (init = 1) or not.
1147 ! ipup   = saved flag to signal updating of newton matrix.
1148 ! jcur   = output flag from svjac showing jacobian status..
1149 !            jcur = 0 means j is not current.
1150 !            jcur = 1 means j is current.
1151 ! jstart = integer flag used as input to svstep..
1152 !            0  means perform the first step.
1153 !            1  means take a new step continuing from the last.
1154 !            -1 means take the next step with a new value of maxord,
1155 !                  hmin, hmxi, n, meth, miter, and/or matrix parameters.
1156 !          on return, svstep sets jstart = 1.
1157 ! jsv    = integer flag for jacobian saving, = sign(mf).
1158 ! kflag  = a completion code from svstep with the following meanings..
1159 !               0      the step was succesful.
1160 !              -1      the requested error could not be achieved.
1161 !              -2      corrector convergence could not be achieved.
1162 !              -3, -4  fatal error in vnls (can not occur here).
1163 ! kuth   = input flag to svstep showing whether h was reduced by the
1164 !          driver.  kuth = 1 if h was reduced, = 0 otherwise.
1165 ! l      = integer variable, nq + 1, current order plus one.
1166 ! lmax   = maxord + 1 (used for dimensioning).
1167 ! locjs  = a pointer to the saved jacobian, whose storage starts at
1168 !          wm(locjs), if jsv = 1.
1169 ! lyh, lewt, lacor, lsavf, lwm, liwm = saved integer pointers
1170 !          to segments of rwork and iwork.
1171 ! maxord = the maximum order of integration method to be allowed.
1172 ! meth/miter = the method flags.  see mf.
1173 ! msbj   = the maximum number of steps between j evaluations, = 50.
1174 ! mxhnil = saved value of optional input mxhnil.
1175 ! mxstep = saved value of optional input mxstep.
1176 ! n      = the number of first-order odes, = neq.
1177 ! newh   = saved integer to flag change of h.
1178 ! newq   = the method order to be used on the next step.
1179 ! nhnil  = saved counter for occurrences of t + h = t.
1180 ! nq     = integer variable, the current integration method order.
1181 ! nqnyh  = saved value of nq*nyh.
1182 ! nqwait = a counter controlling the frequency of order changes.
1183 !          an order change is about to be considered if nqwait = 1.
1184 ! nslj   = the number of steps taken as of the last jacobian update.
1185 ! nslp   = saved value of nst as of last newton matrix update.
1186 ! nyh    = saved value of the initial value of neq.
1187 ! hu     = the step size in t last used.
1188 ! ncfn   = number of nonlinear convergence failures so far.
1189 ! netf   = the number of error test failures of the integrator so far.
1190 ! nfe    = the number of f evaluations for the problem so far.
1191 ! nje    = the number of jacobian evaluations so far.
1192 ! nlu    = the number of matrix lu decompositions so far.
1193 ! nni    = number of nonlinear iterations so far.
1194 ! nqu    = the method order last used.
1195 ! nst    = the number of steps taken for the problem so far.
1196 !-----------------------------------------------------------------------
1197       common /svode_cmn_01/ acnrm, ccmxj, conp, crate, drc, el(13),   &
1198                       eta, etamax, h, hmin, hmxi, hnew, hscal, prl1,   &
1199                       rc, rl1, tau(13), tq(5), tn, uround,   &
1200                       icf, init, ipup, jcur, jstart, jsv, kflag, kuth,   &
1201                       l, lmax, lyh, lewt, lacor, lsavf, lwm, liwm,   &
1202                       locjs, maxord, meth, miter, msbj, mxhnil, mxstep,   &
1203                       n, newh, newq, nhnil, nq, nqnyh, nqwait, nslj,   &
1204                       nslp, nyh
1205       common /svode_cmn_02/ hu, ncfn, netf, nfe, nje, nlu, nni, nqu, nst
1207       data  mord(1) /12/, mord(2) /5/, mxstp0 /500/, mxhnl0 /10/
1208       data zero /0.0e0/, one /1.0e0/, two /2.0e0/, four /4.0e0/,   &
1209            pt2 /0.2e0/, hun /100.0e0/
1210 !-----------------------------------------------------------------------
1211 ! block a.
1212 ! this code block is executed on every call.
1213 ! it tests istate and itask for legality and branches appropriately.
1214 ! if istate .gt. 1 but the flag init shows that initialization has
1215 ! not yet been done, an error return occurs.
1216 ! if istate = 1 and tout = t, return immediately.
1217 !-----------------------------------------------------------------------
1218       if (istate .lt. 1 .or. istate .gt. 3) go to 601
1219       if (itask .lt. 1 .or. itask .gt. 5) go to 602
1220       if (istate .eq. 1) go to 10
1221       if (init .ne. 1) go to 603
1222       if (istate .eq. 2) go to 200
1223       go to 20
1224  10   init = 0
1225       if (tout .eq. t) return
1226 !-----------------------------------------------------------------------
1227 ! block b.
1228 ! the next code block is executed for the initial call (istate = 1),
1229 ! or for a continuation call with parameter changes (istate = 3).
1230 ! it contains checking of all input and various initializations.
1232 ! first check legality of the non-optional input neq, itol, iopt,
1233 ! mf, ml, and mu.
1234 !-----------------------------------------------------------------------
1235  20   if (neq .le. 0) go to 604
1236       if (istate .eq. 1) go to 25
1237       if (neq .gt. n) go to 605
1238  25   n = neq
1239       if (itol .lt. 1 .or. itol .gt. 4) go to 606
1240       if (iopt .lt. 0 .or. iopt .gt. 1) go to 607
1241       jsv = sign(1,mf)
1242       mfa = abs(mf)
1243       meth = mfa/10
1244       miter = mfa - 10*meth
1245       if (meth .lt. 1 .or. meth .gt. 2) go to 608
1246       if (miter .lt. 0 .or. miter .gt. 5) go to 608
1247       if (miter .le. 3) go to 30
1248       ml = iwork(1)
1249       mu = iwork(2)
1250       if (ml .lt. 0 .or. ml .ge. n) go to 609
1251       if (mu .lt. 0 .or. mu .ge. n) go to 610
1252  30   continue
1253 ! next process and check the optional input. ---------------------------
1254       if (iopt .eq. 1) go to 40
1255       maxord = mord(meth)
1256       mxstep = mxstp0
1257       mxhnil = mxhnl0
1258       if (istate .eq. 1) h0 = zero
1259       hmxi = zero
1260       hmin = zero
1261       go to 60
1262  40   maxord = iwork(5)
1263       if (maxord .lt. 0) go to 611
1264       if (maxord .eq. 0) maxord = 100
1265       maxord = min(maxord,mord(meth))
1266       mxstep = iwork(6)
1267       if (mxstep .lt. 0) go to 612
1268       if (mxstep .eq. 0) mxstep = mxstp0
1269       mxhnil = iwork(7)
1270       if (mxhnil .lt. 0) go to 613
1271       if (mxhnil .eq. 0) mxhnil = mxhnl0
1272       if (istate .ne. 1) go to 50
1273       h0 = rwork(5)
1274       if ((tout - t)*h0 .lt. zero) go to 614
1275  50   hmax = rwork(6)
1276       if (hmax .lt. zero) go to 615
1277       hmxi = zero
1278       if (hmax .gt. zero) hmxi = one/hmax
1279       hmin = rwork(7)
1280       if (hmin .lt. zero) go to 616
1281 !-----------------------------------------------------------------------
1282 ! set work array pointers and check lengths lrw and liw.
1283 ! pointers to segments of rwork and iwork are named by prefixing l to
1284 ! the name of the segment.  e.g., the segment yh starts at rwork(lyh).
1285 ! segments of rwork (in order) are denoted  yh, wm, ewt, savf, acor.
1286 ! within wm, locjs is the location of the saved jacobian (jsv .gt. 0).
1287 !-----------------------------------------------------------------------
1288  60   lyh = 21
1289       if (istate .eq. 1) nyh = n
1290       lwm = lyh + (maxord + 1)*nyh
1291       jco = max(0,jsv)
1292       if (miter .eq. 0) lenwm = 0
1293       if (miter .eq. 1 .or. miter .eq. 2) then
1294         lenwm = 2 + (1 + jco)*n*n
1295         locjs = n*n + 3
1296       endif
1297       if (miter .eq. 3) lenwm = 2 + n
1298       if (miter .eq. 4 .or. miter .eq. 5) then
1299         mband = ml + mu + 1
1300         lenp = (mband + ml)*n
1301         lenj = mband*n
1302         lenwm = 2 + lenp + jco*lenj
1303         locjs = lenp + 3
1304         endif
1305       lewt = lwm + lenwm
1306       lsavf = lewt + n
1307       lacor = lsavf + n
1308       lenrw = lacor + n - 1
1309       iwork(17) = lenrw
1310       liwm = 1
1311       leniw = 30 + n
1312       if (miter .eq. 0 .or. miter .eq. 3) leniw = 30
1313       iwork(18) = leniw
1314       if (lenrw .gt. lrw) go to 617
1315       if (leniw .gt. liw) go to 618
1316 ! check rtol and atol for legality. ------------------------------------
1317       rtoli = rtol(1)
1318       atoli = atol(1)
1319       do 70 i = 1,n
1320         if (itol .ge. 3) rtoli = rtol(i)
1321         if (itol .eq. 2 .or. itol .eq. 4) atoli = atol(i)
1322         if (rtoli .lt. zero) go to 619
1323         if (atoli .lt. zero) go to 620
1324  70     continue
1325       if (istate .eq. 1) go to 100
1326 ! if istate = 3, set flag to signal parameter changes to svstep. -------
1327       jstart = -1
1328       if (nq .le. maxord) go to 90
1329 ! maxord was reduced below nq.  copy yh(*,maxord+2) into savf. ---------
1330       call scopy (n, rwork(lwm), 1, rwork(lsavf), 1)
1331 ! reload wm(1) = rwork(lwm), since lwm may have changed. ---------------
1332  90   if (miter .gt. 0) rwork(lwm) = sqrt(uround)
1333 !-----------------------------------------------------------------------
1334 ! block c.
1335 ! the next block is for the initial call only (istate = 1).
1336 ! it contains all remaining initializations, the initial call to f,
1337 ! and the calculation of the initial step size.
1338 ! the error weights in ewt are inverted after being loaded.
1339 !-----------------------------------------------------------------------
1340  100  uround = r1mach(4)
1341       tn = t
1342       if (itask .ne. 4 .and. itask .ne. 5) go to 110
1343       tcrit = rwork(1)
1344       if ((tcrit - tout)*(tout - t) .lt. zero) go to 625
1345       if (h0 .ne. zero .and. (t + h0 - tcrit)*h0 .gt. zero)   &
1346          h0 = tcrit - t
1347  110  jstart = 0
1348       if (miter .gt. 0) rwork(lwm) = sqrt(uround)
1349       ccmxj = pt2
1350       msbj = 50
1351       nhnil = 0
1352       nst = 0
1353       nje = 0
1354       nni = 0
1355       ncfn = 0
1356       netf = 0
1357       nlu = 0
1358       nslj = 0
1359       nslast = 0
1360       hu = zero
1361       nqu = 0
1362 ! initial call to f.  (lf0 points to yh(*,2).) -------------------------
1363       lf0 = lyh + nyh
1364       call f (n, t, y, rwork(lf0), rpar, ipar)
1365       nfe = 1
1366 ! load the initial value vector in yh. ---------------------------------
1367       call scopy (n, y, 1, rwork(lyh), 1)
1368 ! load and invert the ewt array.  (h is temporarily set to 1.0.) -------
1369       nq = 1
1370       h = one
1371       call sewset (n, itol, rtol, atol, rwork(lyh), rwork(lewt))
1372       do 120 i = 1,n
1373         if (rwork(i+lewt-1) .le. zero) go to 621
1374  120    rwork(i+lewt-1) = one/rwork(i+lewt-1)
1375       if (h0 .ne. zero) go to 180
1376 ! call svhin to set initial step size h0 to be attempted. --------------
1377       call svhin (n, t, rwork(lyh), rwork(lf0), f, rpar, ipar, tout,   &
1378          uround, rwork(lewt), itol, atol, y, rwork(lacor), h0,   &
1379          niter, ier)
1380       nfe = nfe + niter
1381       if (ier .ne. 0) go to 622
1382 ! adjust h0 if necessary to meet hmax bound. ---------------------------
1383  180  rh = abs(h0)*hmxi
1384       if (rh .gt. one) h0 = h0/rh
1385 ! load h with h0 and scale yh(*,2) by h0. ------------------------------
1386       h = h0
1387       call sscal (n, h0, rwork(lf0), 1)
1388       go to 270
1389 !-----------------------------------------------------------------------
1390 ! block d.
1391 ! the next code block is for continuation calls only (istate = 2 or 3)
1392 ! and is to check stop conditions before taking a step.
1393 !-----------------------------------------------------------------------
1394  200  nslast = nst
1395       kuth = 0
1396       go to (210, 250, 220, 230, 240), itask
1397  210  if ((tn - tout)*h .lt. zero) go to 250
1398       call svindy (tout, 0, rwork(lyh), nyh, y, iflag)
1399       if (iflag .ne. 0) go to 627
1400       t = tout
1401       go to 420
1402  220  tp = tn - hu*(one + hun*uround)
1403       if ((tp - tout)*h .gt. zero) go to 623
1404       if ((tn - tout)*h .lt. zero) go to 250
1405       go to 400
1406  230  tcrit = rwork(1)
1407       if ((tn - tcrit)*h .gt. zero) go to 624
1408       if ((tcrit - tout)*h .lt. zero) go to 625
1409       if ((tn - tout)*h .lt. zero) go to 245
1410       call svindy (tout, 0, rwork(lyh), nyh, y, iflag)
1411       if (iflag .ne. 0) go to 627
1412       t = tout
1413       go to 420
1414  240  tcrit = rwork(1)
1415       if ((tn - tcrit)*h .gt. zero) go to 624
1416  245  hmx = abs(tn) + abs(h)
1417       ihit = abs(tn - tcrit) .le. hun*uround*hmx
1418       if (ihit) go to 400
1419       tnext = tn + hnew*(one + four*uround)
1420       if ((tnext - tcrit)*h .le. zero) go to 250
1421       h = (tcrit - tn)*(one - four*uround)
1422       kuth = 1
1423 !-----------------------------------------------------------------------
1424 ! block e.
1425 ! the next block is normally executed for all calls and contains
1426 ! the call to the one-step core integrator svstep.
1428 ! this is a looping point for the integration steps.
1430 ! first check for too many steps being taken, update ewt (if not at
1431 ! start of problem), check for too much accuracy being requested, and
1432 ! check for h below the roundoff level in t.
1433 !-----------------------------------------------------------------------
1434  250  continue
1435       if ((nst-nslast) .ge. mxstep) go to 500
1436       call sewset (n, itol, rtol, atol, rwork(lyh), rwork(lewt))
1437       do 260 i = 1,n
1438         if (rwork(i+lewt-1) .le. zero) go to 510
1439  260    rwork(i+lewt-1) = one/rwork(i+lewt-1)
1440  270  tolsf = uround*svnorm (n, rwork(lyh), rwork(lewt))
1441       if (tolsf .le. one) go to 280
1442       tolsf = tolsf*two
1443       if (nst .eq. 0) go to 626
1444       go to 520
1445  280  if ((tn + h) .ne. tn) go to 290
1446       nhnil = nhnil + 1
1447       if (nhnil .gt. mxhnil) go to 290
1448       msg = 'svode--  warning..internal t (=r1) and h (=r2) are'
1449       call xerrwv (msg, 50, 101, 1, 0, 0, 0, 0, zero, zero)
1450       msg='      such that in the machine, t + h = t on the next step  '
1451       call xerrwv (msg, 60, 101, 1, 0, 0, 0, 0, zero, zero)
1452       msg = '      (h = step size). solver will continue anyway'
1453       call xerrwv (msg, 50, 101, 1, 0, 0, 0, 2, tn, h)
1454       if (nhnil .lt. mxhnil) go to 290
1455       msg = 'svode--  above warning has been issued i1 times.  '
1456       call xerrwv (msg, 50, 102, 1, 0, 0, 0, 0, zero, zero)
1457       msg = '      it will not be issued again for this problem'
1458       call xerrwv (msg, 50, 102, 1, 1, mxhnil, 0, 0, zero, zero)
1459  290  continue
1460 !-----------------------------------------------------------------------
1461 ! call svstep (y, yh, nyh, yh, ewt, savf, vsav, acor,
1462 !              wm, iwm, f, jac, f, svnlsd, rpar, ipar)
1463 !-----------------------------------------------------------------------
1464       call svstep (y, rwork(lyh), nyh, rwork(lyh), rwork(lewt),   &
1465          rwork(lsavf), y, rwork(lacor), rwork(lwm), iwork(liwm),   &
1466          f, jac, f, svnlsd, rpar, ipar)
1467       kgo = 1 - kflag
1468 ! branch on kflag.  note..in this version, kflag can not be set to -3.
1469 !  kflag .eq. 0,   -1,  -2
1470       go to (300, 530, 540), kgo
1471 !-----------------------------------------------------------------------
1472 ! block f.
1473 ! the following block handles the case of a successful return from the
1474 ! core integrator (kflag = 0).  test for stop conditions.
1475 !-----------------------------------------------------------------------
1476  300  init = 1
1477       kuth = 0
1478       go to (310, 400, 330, 340, 350), itask
1479 ! itask = 1.  if tout has been reached, interpolate. -------------------
1480  310  if ((tn - tout)*h .lt. zero) go to 250
1481       call svindy (tout, 0, rwork(lyh), nyh, y, iflag)
1482       t = tout
1483       go to 420
1484 ! itask = 3.  jump to exit if tout was reached. ------------------------
1485  330  if ((tn - tout)*h .ge. zero) go to 400
1486       go to 250
1487 ! itask = 4.  see if tout or tcrit was reached.  adjust h if necessary.
1488  340  if ((tn - tout)*h .lt. zero) go to 345
1489       call svindy (tout, 0, rwork(lyh), nyh, y, iflag)
1490       t = tout
1491       go to 420
1492  345  hmx = abs(tn) + abs(h)
1493       ihit = abs(tn - tcrit) .le. hun*uround*hmx
1494       if (ihit) go to 400
1495       tnext = tn + hnew*(one + four*uround)
1496       if ((tnext - tcrit)*h .le. zero) go to 250
1497       h = (tcrit - tn)*(one - four*uround)
1498       kuth = 1
1499       go to 250
1500 ! itask = 5.  see if tcrit was reached and jump to exit. ---------------
1501  350  hmx = abs(tn) + abs(h)
1502       ihit = abs(tn - tcrit) .le. hun*uround*hmx
1503 !-----------------------------------------------------------------------
1504 ! block g.
1505 ! the following block handles all successful returns from svode.
1506 ! if itask .ne. 1, y is loaded from yh and t is set accordingly.
1507 ! istate is set to 2, and the optional output is loaded into the work
1508 ! arrays before returning.
1509 !-----------------------------------------------------------------------
1510  400  continue
1511       call scopy (n, rwork(lyh), 1, y, 1)
1512       t = tn
1513       if (itask .ne. 4 .and. itask .ne. 5) go to 420
1514       if (ihit) t = tcrit
1515  420  istate = 2
1516       rwork(11) = hu
1517       rwork(12) = hnew
1518       rwork(13) = tn
1519       iwork(11) = nst
1520       iwork(12) = nfe
1521       iwork(13) = nje
1522       iwork(14) = nqu
1523       iwork(15) = newq
1524       iwork(19) = nlu
1525       iwork(20) = nni
1526       iwork(21) = ncfn
1527       iwork(22) = netf
1528       return
1529 !-----------------------------------------------------------------------
1530 ! block h.
1531 ! the following block handles all unsuccessful returns other than
1532 ! those for illegal input.  first the error message routine is called.
1533 ! if there was an error test or convergence test failure, imxer is set.
1534 ! then y is loaded from yh, and t is set to tn.
1535 ! the optional output is loaded into the work arrays before returning.
1536 !-----------------------------------------------------------------------
1537 ! the maximum number of steps was taken before reaching tout. ----------
1538  500  msg = 'svode--  at current t (=r1), mxstep (=i1) steps   '
1539       call xerrwv (msg, 50, 201, 1, 0, 0, 0, 0, zero, zero)
1540       msg = '      taken on this call before reaching tout     '
1541       call xerrwv (msg, 50, 201, 1, 1, mxstep, 0, 1, tn, zero)
1542       istate = -1
1543       go to 580
1544 ! ewt(i) .le. 0.0 for some i (not at start of problem). ----------------
1545  510  ewti = rwork(lewt+i-1)
1546       msg = 'svode--  at t (=r1), ewt(i1) has become r2 .le. 0.'
1547       call xerrwv (msg, 50, 202, 1, 1, i, 0, 2, tn, ewti)
1548       istate = -6
1549       go to 580
1550 ! too much accuracy requested for machine precision. -------------------
1551  520  msg = 'svode--  at t (=r1), too much accuracy requested  '
1552       call xerrwv (msg, 50, 203, 1, 0, 0, 0, 0, zero, zero)
1553       msg = '      for precision of machine..  see tolsf (=r2) '
1554       call xerrwv (msg, 50, 203, 1, 0, 0, 0, 2, tn, tolsf)
1555       rwork(14) = tolsf
1556       istate = -2
1557       go to 580
1558 ! kflag = -1.  error test failed repeatedly or with abs(h) = hmin. -----
1559  530  msg = 'svode--  at t(=r1) and step size h(=r2), the error'
1560       call xerrwv (msg, 50, 204, 1, 0, 0, 0, 0, zero, zero)
1561       msg = '      test failed repeatedly or with abs(h) = hmin'
1562       call xerrwv (msg, 50, 204, 1, 0, 0, 0, 2, tn, h)
1563       istate = -4
1564       go to 560
1565 ! kflag = -2.  convergence failed repeatedly or with abs(h) = hmin. ----
1566  540  msg = 'svode--  at t (=r1) and step size h (=r2), the    '
1567       call xerrwv (msg, 50, 205, 1, 0, 0, 0, 0, zero, zero)
1568       msg = '      corrector convergence failed repeatedly     '
1569       call xerrwv (msg, 50, 205, 1, 0, 0, 0, 0, zero, zero)
1570       msg = '      or with abs(h) = hmin   '
1571       call xerrwv (msg, 30, 205, 1, 0, 0, 0, 2, tn, h)
1572       istate = -5
1573 ! compute imxer if relevant. -------------------------------------------
1574  560  big = zero
1575       imxer = 1
1576       do 570 i = 1,n
1577         size = abs(rwork(i+lacor-1)*rwork(i+lewt-1))
1578         if (big .ge. size) go to 570
1579         big = size
1580         imxer = i
1581  570    continue
1582       iwork(16) = imxer
1583 ! set y vector, t, and optional output. --------------------------------
1584  580  continue
1585       call scopy (n, rwork(lyh), 1, y, 1)
1586       t = tn
1587       rwork(11) = hu
1588       rwork(12) = h
1589       rwork(13) = tn
1590       iwork(11) = nst
1591       iwork(12) = nfe
1592       iwork(13) = nje
1593       iwork(14) = nqu
1594       iwork(15) = nq
1595       iwork(19) = nlu
1596       iwork(20) = nni
1597       iwork(21) = ncfn
1598       iwork(22) = netf
1599       return
1600 !-----------------------------------------------------------------------
1601 ! block i.
1602 ! the following block handles all error returns due to illegal input
1603 ! (istate = -3), as detected before calling the core integrator.
1604 ! first the error message routine is called.   if the illegal input
1605 ! is a negative istate, the run is aborted (apparent infinite loop).
1606 !-----------------------------------------------------------------------
1607  601  msg = 'svode--  istate (=i1) illegal '
1608       call xerrwv (msg, 30, 1, 1, 1, istate, 0, 0, zero, zero)
1609       if (istate .lt. 0) go to 800
1610       go to 700
1611  602  msg = 'svode--  itask (=i1) illegal  '
1612       call xerrwv (msg, 30, 2, 1, 1, itask, 0, 0, zero, zero)
1613       go to 700
1614  603  msg='svode--  istate (=i1) .gt. 1 but svode not initialized      '
1615       call xerrwv (msg, 60, 3, 1, 1, istate, 0, 0, zero, zero)
1616       go to 700
1617  604  msg = 'svode--  neq (=i1) .lt. 1     '
1618       call xerrwv (msg, 30, 4, 1, 1, neq, 0, 0, zero, zero)
1619       go to 700
1620  605  msg = 'svode--  istate = 3 and neq increased (i1 to i2)  '
1621       call xerrwv (msg, 50, 5, 1, 2, n, neq, 0, zero, zero)
1622       go to 700
1623  606  msg = 'svode--  itol (=i1) illegal   '
1624       call xerrwv (msg, 30, 6, 1, 1, itol, 0, 0, zero, zero)
1625       go to 700
1626  607  msg = 'svode--  iopt (=i1) illegal   '
1627       call xerrwv (msg, 30, 7, 1, 1, iopt, 0, 0, zero, zero)
1628       go to 700
1629  608  msg = 'svode--  mf (=i1) illegal     '
1630       call xerrwv (msg, 30, 8, 1, 1, mf, 0, 0, zero, zero)
1631       go to 700
1632  609  msg = 'svode--  ml (=i1) illegal.. .lt.0 or .ge.neq (=i2)'
1633       call xerrwv (msg, 50, 9, 1, 2, ml, neq, 0, zero, zero)
1634       go to 700
1635  610  msg = 'svode--  mu (=i1) illegal.. .lt.0 or .ge.neq (=i2)'
1636       call xerrwv (msg, 50, 10, 1, 2, mu, neq, 0, zero, zero)
1637       go to 700
1638  611  msg = 'svode--  maxord (=i1) .lt. 0  '
1639       call xerrwv (msg, 30, 11, 1, 1, maxord, 0, 0, zero, zero)
1640       go to 700
1641  612  msg = 'svode--  mxstep (=i1) .lt. 0  '
1642       call xerrwv (msg, 30, 12, 1, 1, mxstep, 0, 0, zero, zero)
1643       go to 700
1644  613  msg = 'svode--  mxhnil (=i1) .lt. 0  '
1645       call xerrwv (msg, 30, 13, 1, 1, mxhnil, 0, 0, zero, zero)
1646       go to 700
1647  614  msg = 'svode--  tout (=r1) behind t (=r2)      '
1648       call xerrwv (msg, 40, 14, 1, 0, 0, 0, 2, tout, t)
1649       msg = '      integration direction is given by h0 (=r1)  '
1650       call xerrwv (msg, 50, 14, 1, 0, 0, 0, 1, h0, zero)
1651       go to 700
1652  615  msg = 'svode--  hmax (=r1) .lt. 0.0  '
1653       call xerrwv (msg, 30, 15, 1, 0, 0, 0, 1, hmax, zero)
1654       go to 700
1655  616  msg = 'svode--  hmin (=r1) .lt. 0.0  '
1656       call xerrwv (msg, 30, 16, 1, 0, 0, 0, 1, hmin, zero)
1657       go to 700
1658  617  continue
1659       msg='svode--  rwork length needed, lenrw (=i1), exceeds lrw (=i2)'
1660       call xerrwv (msg, 60, 17, 1, 2, lenrw, lrw, 0, zero, zero)
1661       go to 700
1662  618  continue
1663       msg='svode--  iwork length needed, leniw (=i1), exceeds liw (=i2)'
1664       call xerrwv (msg, 60, 18, 1, 2, leniw, liw, 0, zero, zero)
1665       go to 700
1666  619  msg = 'svode--  rtol(i1) is r1 .lt. 0.0        '
1667       call xerrwv (msg, 40, 19, 1, 1, i, 0, 1, rtoli, zero)
1668       go to 700
1669  620  msg = 'svode--  atol(i1) is r1 .lt. 0.0        '
1670       call xerrwv (msg, 40, 20, 1, 1, i, 0, 1, atoli, zero)
1671       go to 700
1672  621  ewti = rwork(lewt+i-1)
1673       msg = 'svode--  ewt(i1) is r1 .le. 0.0         '
1674       call xerrwv (msg, 40, 21, 1, 1, i, 0, 1, ewti, zero)
1675       go to 700
1676  622  continue
1677       msg='svode--  tout (=r1) too close to t(=r2) to start integration'
1678       call xerrwv (msg, 60, 22, 1, 0, 0, 0, 2, tout, t)
1679       go to 700
1680  623  continue
1681       msg='svode--  itask = i1 and tout (=r1) behind tcur - hu (= r2)  '
1682       call xerrwv (msg, 60, 23, 1, 1, itask, 0, 2, tout, tp)
1683       go to 700
1684  624  continue
1685       msg='svode--  itask = 4 or 5 and tcrit (=r1) behind tcur (=r2)   '
1686       call xerrwv (msg, 60, 24, 1, 0, 0, 0, 2, tcrit, tn)
1687       go to 700
1688  625  continue
1689       msg='svode--  itask = 4 or 5 and tcrit (=r1) behind tout (=r2)   '
1690       call xerrwv (msg, 60, 25, 1, 0, 0, 0, 2, tcrit, tout)
1691       go to 700
1692  626  msg = 'svode--  at start of problem, too much accuracy   '
1693       call xerrwv (msg, 50, 26, 1, 0, 0, 0, 0, zero, zero)
1694       msg='      requested for precision of machine..  see tolsf (=r1) '
1695       call xerrwv (msg, 60, 26, 1, 0, 0, 0, 1, tolsf, zero)
1696       rwork(14) = tolsf
1697       go to 700
1698  627  msg='svode--  trouble from svindy.  itask = i1, tout = r1.       '
1699       call xerrwv (msg, 60, 27, 1, 1, itask, 0, 1, tout, zero)
1701  700  continue
1702       istate = -3
1703       return
1705  800  msg = 'svode--  run aborted.. apparent infinite loop     '
1706       call xerrwv (msg, 50, 303, 2, 0, 0, 0, 0, zero, zero)
1707       return
1708 !----------------------- end of subroutine svode -----------------------
1709       end subroutine svode                                              
1710 !*deck svhin
1711       subroutine svhin (n, t0, y0, ydot, f, rpar, ipar, tout, uround,   &
1712          ewt, itol, atol, y, temp, h0, niter, ier)
1713       external f
1714       real t0, y0, ydot, rpar, tout, uround, ewt, atol, y,   &
1715          temp, h0
1716       integer n, ipar, itol, niter, ier
1717       dimension y0(*), ydot(*), ewt(*), atol(*), y(*),   &
1718          temp(*), rpar(*), ipar(*)
1719 !-----------------------------------------------------------------------
1720 ! call sequence input -- n, t0, y0, ydot, f, rpar, ipar, tout, uround,
1721 !                        ewt, itol, atol, y, temp
1722 ! call sequence output -- h0, niter, ier
1723 ! common block variables accessed -- none
1725 ! subroutines called by svhin.. f
1726 ! function routines called by svhin.. svnorm
1727 !-----------------------------------------------------------------------
1728 ! this routine computes the step size, h0, to be attempted on the
1729 ! first step, when the user has not supplied a value for this.
1731 ! first we check that tout - t0 differs significantly from zero.  then
1732 ! an iteration is done to approximate the initial second derivative
1733 ! and this is used to define h from w.r.m.s.norm(h**2 * yddot / 2) = 1.
1734 ! a bias factor of 1/2 is applied to the resulting h.
1735 ! the sign of h0 is inferred from the initial values of tout and t0.
1737 ! communication with svhin is done with the following variables..
1739 ! n      = size of ode system, input.
1740 ! t0     = initial value of independent variable, input.
1741 ! y0     = vector of initial conditions, input.
1742 ! ydot   = vector of initial first derivatives, input.
1743 ! f      = name of subroutine for right-hand side f(t,y), input.
1744 ! rpar, ipar = dummy names for user's real and integer work arrays.
1745 ! tout   = first output value of independent variable
1746 ! uround = machine unit roundoff
1747 ! ewt, itol, atol = error weights and tolerance parameters
1748 !                   as described in the driver routine, input.
1749 ! y, temp = work arrays of length n.
1750 ! h0     = step size to be attempted, output.
1751 ! niter  = number of iterations (and of f evaluations) to compute h0,
1752 !          output.
1753 ! ier    = the error flag, returned with the value
1754 !          ier = 0  if no trouble occurred, or
1755 !          ier = -1 if tout and t0 are considered too close to proceed.
1756 !-----------------------------------------------------------------------
1758 ! type declarations for local variables --------------------------------
1760       real afi, atoli, delyi, h, half, hg, hlb, hnew, hrat,   &
1761            hub, hun, pt1, t1, tdist, tround, two, yddnrm
1762       integer i, iter
1764 ! type declaration for function subroutines called ---------------------
1766 !     real svnorm      ! rce 2005-jan-21 - module conversion
1767 !-----------------------------------------------------------------------
1768 ! the following fortran-77 declaration is to cause the values of the
1769 ! listed (local) variables to be saved between calls to this integrator.
1770 !-----------------------------------------------------------------------
1771       save half, hun, pt1, two
1772       data half /0.5e0/, hun /100.0e0/, pt1 /0.1e0/, two /2.0e0/
1774       niter = 0
1775       tdist = abs(tout - t0)
1776       tround = uround*max(abs(t0),abs(tout))
1777       if (tdist .lt. two*tround) go to 100
1779 ! set a lower bound on h based on the roundoff level in t0 and tout. ---
1780       hlb = hun*tround
1781 ! set an upper bound on h based on tout-t0 and the initial y and ydot. -
1782       hub = pt1*tdist
1783       atoli = atol(1)
1784       do 10 i = 1, n
1785         if (itol .eq. 2 .or. itol .eq. 4) atoli = atol(i)
1786         delyi = pt1*abs(y0(i)) + atoli
1787         afi = abs(ydot(i))
1788         if (afi*hub .gt. delyi) hub = delyi/afi
1789  10     continue
1791 ! set initial guess for h as geometric mean of upper and lower bounds. -
1792       iter = 0
1793       hg = sqrt(hlb*hub)
1794 ! if the bounds have crossed, exit with the mean value. ----------------
1795       if (hub .lt. hlb) then
1796         h0 = hg
1797         go to 90
1798       endif
1800 ! looping point for iteration. -----------------------------------------
1801  50   continue
1802 ! estimate the second derivative as a difference quotient in f. --------
1803       h = sign (hg, tout - t0)
1804       t1 = t0 + h
1805       do 60 i = 1, n
1806  60     y(i) = y0(i) + h*ydot(i)
1807       call f (n, t1, y, temp, rpar, ipar)
1808       do 70 i = 1, n
1809  70     temp(i) = (temp(i) - ydot(i))/h
1810       yddnrm = svnorm (n, temp, ewt)
1811 ! get the corresponding new value of h. --------------------------------
1812       if (yddnrm*hub*hub .gt. two) then
1813         hnew = sqrt(two/yddnrm)
1814       else
1815         hnew = sqrt(hg*hub)
1816       endif
1817       iter = iter + 1
1818 !-----------------------------------------------------------------------
1819 ! test the stopping conditions.
1820 ! stop if the new and previous h values differ by a factor of .lt. 2.
1821 ! stop if four iterations have been done.  also, stop with previous h
1822 ! if hnew/hg .gt. 2 after first iteration, as this probably means that
1823 ! the second derivative value is bad because of cancellation error.
1824 !-----------------------------------------------------------------------
1825       if (iter .ge. 4) go to 80
1826       hrat = hnew/hg
1827       if ( (hrat .gt. half) .and. (hrat .lt. two) ) go to 80
1828       if ( (iter .ge. 2) .and. (hnew .gt. two*hg) ) then
1829         hnew = hg
1830         go to 80
1831       endif
1832       hg = hnew
1833       go to 50
1835 ! iteration done.  apply bounds, bias factor, and sign.  then exit. ----
1836  80   h0 = hnew*half
1837       if (h0 .lt. hlb) h0 = hlb
1838       if (h0 .gt. hub) h0 = hub
1839  90   h0 = sign(h0, tout - t0)
1840       niter = iter
1841       ier = 0
1842       return
1843 ! error return for tout - t0 too small. --------------------------------
1844  100  ier = -1
1845       return
1846 !----------------------- end of subroutine svhin -----------------------
1847       end subroutine svhin                                               
1848 !*deck svindy
1849       subroutine svindy (t, k, yh, ldyh, dky, iflag)
1850       real t, yh, dky
1851       integer k, ldyh, iflag
1852       dimension yh(ldyh,*), dky(*)
1853 !-----------------------------------------------------------------------
1854 ! call sequence input -- t, k, yh, ldyh
1855 ! call sequence output -- dky, iflag
1856 ! common block variables accessed..
1857 !     /svode_cmn_01/ --  h, tn, uround, l, n, nq
1858 !     /svode_cmn_02/ --  hu
1860 ! subroutines called by svindy.. sscal, xerrwv
1861 ! function routines called by svindy.. none
1862 !-----------------------------------------------------------------------
1863 ! svindy computes interpolated values of the k-th derivative of the
1864 ! dependent variable vector y, and stores it in dky.  this routine
1865 ! is called within the package with k = 0 and t = tout, but may
1866 ! also be called by the user for any k up to the current order.
1867 ! (see detailed instructions in the usage documentation.)
1868 !-----------------------------------------------------------------------
1869 ! the computed values in dky are gotten by interpolation using the
1870 ! nordsieck history array yh.  this array corresponds uniquely to a
1871 ! vector-valued polynomial of degree nqcur or less, and dky is set
1872 ! to the k-th derivative of this polynomial at t.
1873 ! the formula for dky is..
1874 !              q
1875 !  dky(i)  =  sum  c(j,k) * (t - tn)**(j-k) * h**(-j) * yh(i,j+1)
1876 !             j=k
1877 ! where  c(j,k) = j*(j-1)*...*(j-k+1), q = nqcur, tn = tcur, h = hcur.
1878 ! the quantities  nq = nqcur, l = nq+1, n, tn, and h are
1879 ! communicated by common.  the above sum is done in reverse order.
1880 ! iflag is returned negative if either k or t is out of bounds.
1882 ! discussion above and comments in driver explain all variables.
1883 !-----------------------------------------------------------------------
1885 ! type declarations for labeled common block svode_cmn_01 --------------------
1887       real acnrm, ccmxj, conp, crate, drc, el,   &
1888            eta, etamax, h, hmin, hmxi, hnew, hscal, prl1,   &
1889            rc, rl1, tau, tq, tn, uround
1890       integer icf, init, ipup, jcur, jstart, jsv, kflag, kuth,   &
1891               l, lmax, lyh, lewt, lacor, lsavf, lwm, liwm,   &
1892               locjs, maxord, meth, miter, msbj, mxhnil, mxstep,   &
1893               n, newh, newq, nhnil, nq, nqnyh, nqwait, nslj,   &
1894               nslp, nyh
1896 ! type declarations for labeled common block svode_cmn_02 --------------------
1898       real hu
1899       integer ncfn, netf, nfe, nje, nlu, nni, nqu, nst
1901 ! type declarations for local variables --------------------------------
1903       real c, hun, r, s, tfuzz, tn1, tp, zero
1904       integer i, ic, j, jb, jb2, jj, jj1, jp1
1905       character*80 msg
1906 !-----------------------------------------------------------------------
1907 ! the following fortran-77 declaration is to cause the values of the
1908 ! listed (local) variables to be saved between calls to this integrator.
1909 !-----------------------------------------------------------------------
1910       save hun, zero
1912       common /svode_cmn_01/ acnrm, ccmxj, conp, crate, drc, el(13),   &
1913                       eta, etamax, h, hmin, hmxi, hnew, hscal, prl1,   &
1914                       rc, rl1, tau(13), tq(5), tn, uround,   &
1915                       icf, init, ipup, jcur, jstart, jsv, kflag, kuth,   &
1916                       l, lmax, lyh, lewt, lacor, lsavf, lwm, liwm,   &
1917                       locjs, maxord, meth, miter, msbj, mxhnil, mxstep,   &
1918                       n, newh, newq, nhnil, nq, nqnyh, nqwait, nslj,   &
1919                       nslp, nyh
1920       common /svode_cmn_02/ hu, ncfn, netf, nfe, nje, nlu, nni, nqu, nst
1922       data hun /100.0e0/, zero /0.0e0/
1924       iflag = 0
1925       if (k .lt. 0 .or. k .gt. nq) go to 80
1926       tfuzz = hun*uround*(tn + hu)
1927       tp = tn - hu - tfuzz
1928       tn1 = tn + tfuzz
1929       if ((t-tp)*(t-tn1) .gt. zero) go to 90
1931       s = (t - tn)/h
1932       ic = 1
1933       if (k .eq. 0) go to 15
1934       jj1 = l - k
1935       do 10 jj = jj1, nq
1936  10     ic = ic*jj
1937  15   c = real(ic)
1938       do 20 i = 1, n
1939  20     dky(i) = c*yh(i,l)
1940       if (k .eq. nq) go to 55
1941       jb2 = nq - k
1942       do 50 jb = 1, jb2
1943         j = nq - jb
1944         jp1 = j + 1
1945         ic = 1
1946         if (k .eq. 0) go to 35
1947         jj1 = jp1 - k
1948         do 30 jj = jj1, j
1949  30       ic = ic*jj
1950  35     c = real(ic)
1951         do 40 i = 1, n
1952  40       dky(i) = c*yh(i,jp1) + s*dky(i)
1953  50     continue
1954       if (k .eq. 0) return
1955  55   r = h**(-k)
1956       call sscal (n, r, dky, 1)
1957       return
1959  80   msg = 'svindy-- k (=i1) illegal      '
1960       call xerrwv (msg, 30, 51, 1, 1, k, 0, 0, zero, zero)
1961       iflag = -1
1962       return
1963  90   msg = 'svindy-- t (=r1) illegal      '
1964       call xerrwv (msg, 30, 52, 1, 0, 0, 0, 1, t, zero)
1965       msg='      t not in interval tcur - hu (= r1) to tcur (=r2)      '
1966       call xerrwv (msg, 60, 52, 1, 0, 0, 0, 2, tp, tn)
1967       iflag = -2
1968       return
1969 !----------------------- end of subroutine svindy ----------------------
1970       end subroutine svindy                             
1971 !*deck svstep
1972       subroutine svstep (y, yh, ldyh, yh1, ewt, savf, vsav, acor,   &
1973                         wm, iwm, f, jac, psol, vnls, rpar, ipar)
1974       external f, jac, psol, vnls
1975       real y, yh, yh1, ewt, savf, vsav, acor, wm, rpar
1976       integer ldyh, iwm, ipar
1977       dimension y(*), yh(ldyh,*), yh1(*), ewt(*), savf(*), vsav(*),   &
1978          acor(*), wm(*), iwm(*), rpar(*), ipar(*)
1979 !-----------------------------------------------------------------------
1980 ! call sequence input -- y, yh, ldyh, yh1, ewt, savf, vsav,
1981 !                        acor, wm, iwm, f, jac, psol, vnls, rpar, ipar
1982 ! call sequence output -- yh, acor, wm, iwm
1983 ! common block variables accessed..
1984 !     /svode_cmn_01/  acnrm, el(13), h, hmin, hmxi, hnew, hscal, rc, tau(13),
1985 !               tq(5), tn, jcur, jstart, kflag, kuth,
1986 !               l, lmax, maxord, n, newq, nq, nqwait
1987 !     /svode_cmn_02/  hu, ncfn, netf, nfe, nqu, nst
1989 ! subroutines called by svstep.. f, saxpy, scopy, sscal,
1990 !                               svjust, vnls, svset
1991 ! function routines called by svstep.. svnorm
1992 !-----------------------------------------------------------------------
1993 ! svstep performs one step of the integration of an initial value
1994 ! problem for a system of ordinary differential equations.
1995 ! svstep calls subroutine vnls for the solution of the nonlinear system
1996 ! arising in the time step.  thus it is independent of the problem
1997 ! jacobian structure and the type of nonlinear system solution method.
1998 ! svstep returns a completion flag kflag (in common).
1999 ! a return with kflag = -1 or -2 means either abs(h) = hmin or 10
2000 ! consecutive failures occurred.  on a return with kflag negative,
2001 ! the values of tn and the yh array are as of the beginning of the last
2002 ! step, and h is the last step size attempted.
2004 ! communication with svstep is done with the following variables..
2006 ! y      = an array of length n used for the dependent variable vector.
2007 ! yh     = an ldyh by lmax array containing the dependent variables
2008 !          and their approximate scaled derivatives, where
2009 !          lmax = maxord + 1.  yh(i,j+1) contains the approximate
2010 !          j-th derivative of y(i), scaled by h**j/factorial(j)
2011 !          (j = 0,1,...,nq).  on entry for the first step, the first
2012 !          two columns of yh must be set from the initial values.
2013 ! ldyh   = a constant integer .ge. n, the first dimension of yh.
2014 !          n is the number of odes in the system.
2015 ! yh1    = a one-dimensional array occupying the same space as yh.
2016 ! ewt    = an array of length n containing multiplicative weights
2017 !          for local error measurements.  local errors in y(i) are
2018 !          compared to 1.0/ewt(i) in various error tests.
2019 ! savf   = an array of working storage, of length n.
2020 !          also used for input of yh(*,maxord+2) when jstart = -1
2021 !          and maxord .lt. the current order nq.
2022 ! vsav   = a work array of length n passed to subroutine vnls.
2023 ! acor   = a work array of length n, used for the accumulated
2024 !          corrections.  on a successful return, acor(i) contains
2025 !          the estimated one-step local error in y(i).
2026 ! wm,iwm = real and integer work arrays associated with matrix
2027 !          operations in vnls.
2028 ! f      = dummy name for the user supplied subroutine for f.
2029 ! jac    = dummy name for the user supplied jacobian subroutine.
2030 ! psol   = dummy name for the subroutine passed to vnls, for
2031 !          possible use there.
2032 ! vnls   = dummy name for the nonlinear system solving subroutine,
2033 !          whose real name is dependent on the method used.
2034 ! rpar, ipar = dummy names for user's real and integer work arrays.
2035 !-----------------------------------------------------------------------
2037 ! type declarations for labeled common block svode_cmn_01 --------------------
2039       real acnrm, ccmxj, conp, crate, drc, el,   &
2040            eta, etamax, h, hmin, hmxi, hnew, hscal, prl1,   &
2041            rc, rl1, tau, tq, tn, uround
2042       integer icf, init, ipup, jcur, jstart, jsv, kflag, kuth,   &
2043               l, lmax, lyh, lewt, lacor, lsavf, lwm, liwm,   &
2044               locjs, maxord, meth, miter, msbj, mxhnil, mxstep,   &
2045               n, newh, newq, nhnil, nq, nqnyh, nqwait, nslj,   &
2046               nslp, nyh
2048 ! type declarations for labeled common block svode_cmn_02 --------------------
2050       real hu
2051       integer ncfn, netf, nfe, nje, nlu, nni, nqu, nst
2053 ! type declarations for local variables --------------------------------
2055       real addon, bias1,bias2,bias3, cnquot, ddn, dsm, dup,   &
2056            etacf, etamin, etamx1, etamx2, etamx3, etamxf,   &
2057            etaq, etaqm1, etaqp1, flotl, one, onepsm,   &
2058            r, thresh, told, zero
2059       integer i, i1, i2, iback, j, jb, kfc, kfh, mxncf, ncf, nflag
2061 ! type declaration for function subroutines called ---------------------
2063 !     real svnorm      ! rce 2005-jan-21 - module conversion
2064 !-----------------------------------------------------------------------
2065 ! the following fortran-77 declaration is to cause the values of the
2066 ! listed (local) variables to be saved between calls to this integrator.
2067 !-----------------------------------------------------------------------
2068       save addon, bias1, bias2, bias3,   &
2069            etacf, etamin, etamx1, etamx2, etamx3, etamxf, etaq, etaqm1,   &
2070            kfc, kfh, mxncf, onepsm, thresh, one, zero
2071 !-----------------------------------------------------------------------
2072       common /svode_cmn_01/ acnrm, ccmxj, conp, crate, drc, el(13),   &
2073                       eta, etamax, h, hmin, hmxi, hnew, hscal, prl1,   &
2074                       rc, rl1, tau(13), tq(5), tn, uround,   &
2075                       icf, init, ipup, jcur, jstart, jsv, kflag, kuth,   &
2076                       l, lmax, lyh, lewt, lacor, lsavf, lwm, liwm,   &
2077                       locjs, maxord, meth, miter, msbj, mxhnil, mxstep,   &
2078                       n, newh, newq, nhnil, nq, nqnyh, nqwait, nslj,   &
2079                       nslp, nyh
2080       common /svode_cmn_02/ hu, ncfn, netf, nfe, nje, nlu, nni, nqu, nst
2082       data kfc/-3/, kfh/-7/, mxncf/10/
2083       data addon  /1.0e-6/,    bias1  /6.0e0/,     bias2  /6.0e0/,   &
2084            bias3  /10.0e0/,    etacf  /0.25e0/,    etamin /0.1e0/,   &
2085            etamxf /0.2e0/,     etamx1 /1.0e4/,     etamx2 /10.0e0/,   &
2086            etamx3 /10.0e0/,    onepsm /1.00001e0/, thresh /1.5e0/
2087       data one/1.0e0/, zero/0.0e0/
2089       kflag = 0
2090       told = tn
2091       ncf = 0
2092       jcur = 0
2093       nflag = 0
2094       if (jstart .gt. 0) go to 20
2095       if (jstart .eq. -1) go to 100
2096 !-----------------------------------------------------------------------
2097 ! on the first call, the order is set to 1, and other variables are
2098 ! initialized.  etamax is the maximum ratio by which h can be increased
2099 ! in a single step.  it is normally 10, but is larger during the
2100 ! first step to compensate for the small initial h.  if a failure
2101 ! occurs (in corrector convergence or error test), etamax is set to 1
2102 ! for the next increase.
2103 !-----------------------------------------------------------------------
2104       lmax = maxord + 1
2105       nq = 1
2106       l = 2
2107       nqnyh = nq*ldyh
2108       tau(1) = h
2109       prl1 = one
2110       rc = zero
2111       etamax = etamx1
2112       nqwait = 2
2113       hscal = h
2114       go to 200
2115 !-----------------------------------------------------------------------
2116 ! take preliminary actions on a normal continuation step (jstart.gt.0).
2117 ! if the driver changed h, then eta must be reset and newh set to 1.
2118 ! if a change of order was dictated on the previous step, then
2119 ! it is done here and appropriate adjustments in the history are made.
2120 ! on an order decrease, the history array is adjusted by svjust.
2121 ! on an order increase, the history array is augmented by a column.
2122 ! on a change of step size h, the history array yh is rescaled.
2123 !-----------------------------------------------------------------------
2124  20   continue
2125       if (kuth .eq. 1) then
2126         eta = min(eta,h/hscal)
2127         newh = 1
2128         endif
2129  50   if (newh .eq. 0) go to 200
2130       if (newq .eq. nq) go to 150
2131       if (newq .lt. nq) then
2132         call svjust (yh, ldyh, -1)
2133         nq = newq
2134         l = nq + 1
2135         nqwait = l
2136         go to 150
2137         endif
2138       if (newq .gt. nq) then
2139         call svjust (yh, ldyh, 1)
2140         nq = newq
2141         l = nq + 1
2142         nqwait = l
2143         go to 150
2144       endif
2145 !-----------------------------------------------------------------------
2146 ! the following block handles preliminaries needed when jstart = -1.
2147 ! if n was reduced, zero out part of yh to avoid undefined references.
2148 ! if maxord was reduced to a value less than the tentative order newq,
2149 ! then nq is set to maxord, and a new h ratio eta is chosen.
2150 ! otherwise, we take the same preliminary actions as for jstart .gt. 0.
2151 ! in any case, nqwait is reset to l = nq + 1 to prevent further
2152 ! changes in order for that many steps.
2153 ! the new h ratio eta is limited by the input h if kuth = 1,
2154 ! by hmin if kuth = 0, and by hmxi in any case.
2155 ! finally, the history array yh is rescaled.
2156 !-----------------------------------------------------------------------
2157  100  continue
2158       lmax = maxord + 1
2159       if (n .eq. ldyh) go to 120
2160       i1 = 1 + (newq + 1)*ldyh
2161       i2 = (maxord + 1)*ldyh
2162       if (i1 .gt. i2) go to 120
2163       do 110 i = i1, i2
2164  110    yh1(i) = zero
2165  120  if (newq .le. maxord) go to 140
2166       flotl = real(lmax)
2167       if (maxord .lt. nq-1) then
2168         ddn = svnorm (n, savf, ewt)/tq(1)
2169         eta = one/((bias1*ddn)**(one/flotl) + addon)
2170         endif
2171       if (maxord .eq. nq .and. newq .eq. nq+1) eta = etaq
2172       if (maxord .eq. nq-1 .and. newq .eq. nq+1) then
2173         eta = etaqm1
2174         call svjust (yh, ldyh, -1)
2175         endif
2176       if (maxord .eq. nq-1 .and. newq .eq. nq) then
2177         ddn = svnorm (n, savf, ewt)/tq(1)
2178         eta = one/((bias1*ddn)**(one/flotl) + addon)
2179         call svjust (yh, ldyh, -1)
2180         endif
2181       eta = min(eta,one)
2182       nq = maxord
2183       l = lmax
2184  140  if (kuth .eq. 1) eta = min(eta,abs(h/hscal))
2185       if (kuth .eq. 0) eta = max(eta,hmin/abs(hscal))
2186       eta = eta/max(one,abs(hscal)*hmxi*eta)
2187       newh = 1
2188       nqwait = l
2189       if (newq .le. maxord) go to 50
2190 ! rescale the history array for a change in h by a factor of eta. ------
2191  150  r = one
2192       do 180 j = 2, l
2193         r = r*eta
2194         call sscal (n, r, yh(1,j), 1 )
2195  180    continue
2196       h = hscal*eta
2197       hscal = h
2198       rc = rc*eta
2199       nqnyh = nq*ldyh
2200 !-----------------------------------------------------------------------
2201 ! this section computes the predicted values by effectively
2202 ! multiplying the yh array by the pascal triangle matrix.
2203 ! svset is called to calculate all integration coefficients.
2204 ! rc is the ratio of new to old values of the coefficient h/el(2)=h/l1.
2205 !-----------------------------------------------------------------------
2206  200  tn = tn + h
2207       i1 = nqnyh + 1
2208       do 220 jb = 1, nq
2209         i1 = i1 - ldyh
2210         do 210 i = i1, nqnyh
2211  210      yh1(i) = yh1(i) + yh1(i+ldyh)
2212  220  continue
2213       call svset
2214       rl1 = one/el(2)
2215       rc = rc*(rl1/prl1)
2216       prl1 = rl1
2218 ! call the nonlinear system solver. ------------------------------------
2220       call vnls (y, yh, ldyh, vsav, savf, ewt, acor, iwm, wm,   &
2221                  f, jac, psol, nflag, rpar, ipar)
2223       if (nflag .eq. 0) go to 450
2224 !-----------------------------------------------------------------------
2225 ! the vnls routine failed to achieve convergence (nflag .ne. 0).
2226 ! the yh array is retracted to its values before prediction.
2227 ! the step size h is reduced and the step is retried, if possible.
2228 ! otherwise, an error exit is taken.
2229 !-----------------------------------------------------------------------
2230         ncf = ncf + 1
2231         ncfn = ncfn + 1
2232         etamax = one
2233         tn = told
2234         i1 = nqnyh + 1
2235         do 430 jb = 1, nq
2236           i1 = i1 - ldyh
2237           do 420 i = i1, nqnyh
2238  420        yh1(i) = yh1(i) - yh1(i+ldyh)
2239  430      continue
2240         if (nflag .lt. -1) go to 680
2241         if (abs(h) .le. hmin*onepsm) go to 670
2242         if (ncf .eq. mxncf) go to 670
2243         eta = etacf
2244         eta = max(eta,hmin/abs(h))
2245         nflag = -1
2246         go to 150
2247 !-----------------------------------------------------------------------
2248 ! the corrector has converged (nflag = 0).  the local error test is
2249 ! made and control passes to statement 500 if it fails.
2250 !-----------------------------------------------------------------------
2251  450  continue
2252       dsm = acnrm/tq(2)
2253       if (dsm .gt. one) go to 500
2254 !-----------------------------------------------------------------------
2255 ! after a successful step, update the yh and tau arrays and decrement
2256 ! nqwait.  if nqwait is then 1 and nq .lt. maxord, then acor is saved
2257 ! for use in a possible order increase on the next step.
2258 ! if etamax = 1 (a failure occurred this step), keep nqwait .ge. 2.
2259 !-----------------------------------------------------------------------
2260       kflag = 0
2261       nst = nst + 1
2262       hu = h
2263       nqu = nq
2264       do 470 iback = 1, nq
2265         i = l - iback
2266  470    tau(i+1) = tau(i)
2267       tau(1) = h
2268       do 480 j = 1, l
2269         call saxpy (n, el(j), acor, 1, yh(1,j), 1 )
2270  480    continue
2271       nqwait = nqwait - 1
2272       if ((l .eq. lmax) .or. (nqwait .ne. 1)) go to 490
2273       call scopy (n, acor, 1, yh(1,lmax), 1 )
2274       conp = tq(5)
2275  490  if (etamax .ne. one) go to 560
2276       if (nqwait .lt. 2) nqwait = 2
2277       newq = nq
2278       newh = 0
2279       eta = one
2280       hnew = h
2281       go to 690
2282 !-----------------------------------------------------------------------
2283 ! the error test failed.  kflag keeps track of multiple failures.
2284 ! restore tn and the yh array to their previous values, and prepare
2285 ! to try the step again.  compute the optimum step size for the
2286 ! same order.  after repeated failures, h is forced to decrease
2287 ! more rapidly.
2288 !-----------------------------------------------------------------------
2289  500  kflag = kflag - 1
2290       netf = netf + 1
2291       nflag = -2
2292       tn = told
2293       i1 = nqnyh + 1
2294       do 520 jb = 1, nq
2295         i1 = i1 - ldyh
2296         do 510 i = i1, nqnyh
2297  510      yh1(i) = yh1(i) - yh1(i+ldyh)
2298  520  continue
2299       if (abs(h) .le. hmin*onepsm) go to 660
2300       etamax = one
2301       if (kflag .le. kfc) go to 530
2302 ! compute ratio of new h to current h at the current order. ------------
2303       flotl = real(l)
2304       eta = one/((bias2*dsm)**(one/flotl) + addon)
2305       eta = max(eta,hmin/abs(h),etamin)
2306       if ((kflag .le. -2) .and. (eta .gt. etamxf)) eta = etamxf
2307       go to 150
2308 !-----------------------------------------------------------------------
2309 ! control reaches this section if 3 or more consecutive failures
2310 ! have occurred.  it is assumed that the elements of the yh array
2311 ! have accumulated errors of the wrong order.  the order is reduced
2312 ! by one, if possible.  then h is reduced by a factor of 0.1 and
2313 ! the step is retried.  after a total of 7 consecutive failures,
2314 ! an exit is taken with kflag = -1.
2315 !-----------------------------------------------------------------------
2316  530  if (kflag .eq. kfh) go to 660
2317       if (nq .eq. 1) go to 540
2318       eta = max(etamin,hmin/abs(h))
2319       call svjust (yh, ldyh, -1)
2320       l = nq
2321       nq = nq - 1
2322       nqwait = l
2323       go to 150
2324  540  eta = max(etamin,hmin/abs(h))
2325       h = h*eta
2326       hscal = h
2327       tau(1) = h
2328       call f (n, tn, y, savf, rpar, ipar)
2329       nfe = nfe + 1
2330       do 550 i = 1, n
2331  550    yh(i,2) = h*savf(i)
2332       nqwait = 10
2333       go to 200
2334 !-----------------------------------------------------------------------
2335 ! if nqwait = 0, an increase or decrease in order by one is considered.
2336 ! factors etaq, etaqm1, etaqp1 are computed by which h could
2337 ! be multiplied at order q, q-1, or q+1, respectively.
2338 ! the largest of these is determined, and the new order and
2339 ! step size set accordingly.
2340 ! a change of h or nq is made only if h increases by at least a
2341 ! factor of thresh.  if an order change is considered and rejected,
2342 ! then nqwait is set to 2 (reconsider it after 2 steps).
2343 !-----------------------------------------------------------------------
2344 ! compute ratio of new h to current h at the current order. ------------
2345  560  flotl = real(l)
2346       etaq = one/((bias2*dsm)**(one/flotl) + addon)
2347       if (nqwait .ne. 0) go to 600
2348       nqwait = 2
2349       etaqm1 = zero
2350       if (nq .eq. 1) go to 570
2351 ! compute ratio of new h to current h at the current order less one. ---
2352       ddn = svnorm (n, yh(1,l), ewt)/tq(1)
2353       etaqm1 = one/((bias1*ddn)**(one/(flotl - one)) + addon)
2354  570  etaqp1 = zero
2355       if (l .eq. lmax) go to 580
2356 ! compute ratio of new h to current h at current order plus one. -------
2357       cnquot = (tq(5)/conp)*(h/tau(2))**l
2358       do 575 i = 1, n
2359  575    savf(i) = acor(i) - cnquot*yh(i,lmax)
2360       dup = svnorm (n, savf, ewt)/tq(3)
2361       etaqp1 = one/((bias3*dup)**(one/(flotl + one)) + addon)
2362  580  if (etaq .ge. etaqp1) go to 590
2363       if (etaqp1 .gt. etaqm1) go to 620
2364       go to 610
2365  590  if (etaq .lt. etaqm1) go to 610
2366  600  eta = etaq
2367       newq = nq
2368       go to 630
2369  610  eta = etaqm1
2370       newq = nq - 1
2371       go to 630
2372  620  eta = etaqp1
2373       newq = nq + 1
2374       call scopy (n, acor, 1, yh(1,lmax), 1)
2375 ! test tentative new h against thresh, etamax, and hmxi, then exit. ----
2376  630  if (eta .lt. thresh .or. etamax .eq. one) go to 640
2377       eta = min(eta,etamax)
2378       eta = eta/max(one,abs(h)*hmxi*eta)
2379       newh = 1
2380       hnew = h*eta
2381       go to 690
2382  640  newq = nq
2383       newh = 0
2384       eta = one
2385       hnew = h
2386       go to 690
2387 !-----------------------------------------------------------------------
2388 ! all returns are made through this section.
2389 ! on a successful return, etamax is reset and acor is scaled.
2390 !-----------------------------------------------------------------------
2391  660  kflag = -1
2392       go to 720
2393  670  kflag = -2
2394       go to 720
2395  680  if (nflag .eq. -2) kflag = -3
2396       if (nflag .eq. -3) kflag = -4
2397       go to 720
2398  690  etamax = etamx3
2399       if (nst .le. 10) etamax = etamx2
2400  700  r = one/tq(2)
2401       call sscal (n, r, acor, 1)
2402  720  jstart = 1
2403       return
2404 !----------------------- end of subroutine svstep ----------------------
2405       end subroutine svstep                                          
2406 !*deck svset
2407       subroutine svset
2408 !-----------------------------------------------------------------------
2409 ! call sequence communication.. none
2410 ! common block variables accessed..
2411 !     /svode_cmn_01/ -- el(13), h, tau(13), tq(5), l(= nq + 1),
2412 !                 meth, nq, nqwait
2414 ! subroutines called by svset.. none
2415 ! function routines called by svset.. none
2416 !-----------------------------------------------------------------------
2417 ! svset is called by svstep and sets coefficients for use there.
2419 ! for each order nq, the coefficients in el are calculated by use of
2420 !  the generating polynomial lambda(x), with coefficients el(i).
2421 !      lambda(x) = el(1) + el(2)*x + ... + el(nq+1)*(x**nq).
2422 ! for the backward differentiation formulas,
2423 !                                     nq-1
2424 !      lambda(x) = (1 + x/xi*(nq)) * product (1 + x/xi(i) ) .
2425 !                                     i = 1
2426 ! for the adams formulas,
2427 !                              nq-1
2428 !      (d/dx) lambda(x) = c * product (1 + x/xi(i) ) ,
2429 !                              i = 1
2430 !      lambda(-1) = 0,    lambda(0) = 1,
2431 ! where c is a normalization constant.
2432 ! in both cases, xi(i) is defined by
2433 !      h*xi(i) = t sub n  -  t sub (n-i)
2434 !              = h + tau(1) + tau(2) + ... tau(i-1).
2437 ! in addition to variables described previously, communication
2438 ! with svset uses the following..
2439 !   tau    = a vector of length 13 containing the past nq values
2440 !            of h.
2441 !   el     = a vector of length 13 in which vset stores the
2442 !            coefficients for the corrector formula.
2443 !   tq     = a vector of length 5 in which vset stores constants
2444 !            used for the convergence test, the error test, and the
2445 !            selection of h at a new order.
2446 !   meth   = the basic method indicator.
2447 !   nq     = the current order.
2448 !   l      = nq + 1, the length of the vector stored in el, and
2449 !            the number of columns of the yh array being used.
2450 !   nqwait = a counter controlling the frequency of order changes.
2451 !            an order change is about to be considered if nqwait = 1.
2452 !-----------------------------------------------------------------------
2454 ! type declarations for labeled common block svode_cmn_01 --------------------
2456       real acnrm, ccmxj, conp, crate, drc, el,   &
2457            eta, etamax, h, hmin, hmxi, hnew, hscal, prl1,   &
2458            rc, rl1, tau, tq, tn, uround
2459       integer icf, init, ipup, jcur, jstart, jsv, kflag, kuth,   &
2460               l, lmax, lyh, lewt, lacor, lsavf, lwm, liwm,   &
2461               locjs, maxord, meth, miter, msbj, mxhnil, mxstep,   &
2462               n, newh, newq, nhnil, nq, nqnyh, nqwait, nslj,   &
2463               nslp, nyh
2465 ! type declarations for local variables --------------------------------
2467       real ahatn0, alph0, cnqm1, cortes, csum, elp, em,   &
2468            em0, floti, flotl, flotnq, hsum, one, rxi, rxis, s, six,   &
2469            t1, t2, t3, t4, t5, t6, two, xi, zero
2470       integer i, iback, j, jp1, nqm1, nqm2
2472       dimension em(13)
2473 !-----------------------------------------------------------------------
2474 ! the following fortran-77 declaration is to cause the values of the
2475 ! listed (local) variables to be saved between calls to this integrator.
2476 !-----------------------------------------------------------------------
2477       save cortes, one, six, two, zero
2479       common /svode_cmn_01/ acnrm, ccmxj, conp, crate, drc, el(13),   &
2480                       eta, etamax, h, hmin, hmxi, hnew, hscal, prl1,   &
2481                       rc, rl1, tau(13), tq(5), tn, uround,   &
2482                       icf, init, ipup, jcur, jstart, jsv, kflag, kuth,   &
2483                       l, lmax, lyh, lewt, lacor, lsavf, lwm, liwm,   &
2484                       locjs, maxord, meth, miter, msbj, mxhnil, mxstep,   &
2485                       n, newh, newq, nhnil, nq, nqnyh, nqwait, nslj,   &
2486                       nslp, nyh
2488       data cortes /0.1e0/
2489       data one  /1.0e0/, six /6.0e0/, two /2.0e0/, zero /0.0e0/
2491       flotl = real(l)
2492       nqm1 = nq - 1
2493       nqm2 = nq - 2
2494       go to (100, 200), meth
2496 ! set coefficients for adams methods. ----------------------------------
2497  100  if (nq .ne. 1) go to 110
2498       el(1) = one
2499       el(2) = one
2500       tq(1) = one
2501       tq(2) = two
2502       tq(3) = six*tq(2)
2503       tq(5) = one
2504       go to 300
2505  110  hsum = h
2506       em(1) = one
2507       flotnq = flotl - one
2508       do 115 i = 2, l
2509  115    em(i) = zero
2510       do 150 j = 1, nqm1
2511         if ((j .ne. nqm1) .or. (nqwait .ne. 1)) go to 130
2512         s = one
2513         csum = zero
2514         do 120 i = 1, nqm1
2515           csum = csum + s*em(i)/real(i+1)
2516  120      s = -s
2517         tq(1) = em(nqm1)/(flotnq*csum)
2518  130    rxi = h/hsum
2519         do 140 iback = 1, j
2520           i = (j + 2) - iback
2521  140      em(i) = em(i) + em(i-1)*rxi
2522         hsum = hsum + tau(j)
2523  150    continue
2524 ! compute integral from -1 to 0 of polynomial and of x times it. -------
2525       s = one
2526       em0 = zero
2527       csum = zero
2528       do 160 i = 1, nq
2529         floti = real(i)
2530         em0 = em0 + s*em(i)/floti
2531         csum = csum + s*em(i)/(floti+one)
2532  160    s = -s
2533 ! in el, form coefficients of normalized integrated polynomial. --------
2534       s = one/em0
2535       el(1) = one
2536       do 170 i = 1, nq
2537  170    el(i+1) = s*em(i)/real(i)
2538       xi = hsum/h
2539       tq(2) = xi*em0/csum
2540       tq(5) = xi/el(l)
2541       if (nqwait .ne. 1) go to 300
2542 ! for higher order control constant, multiply polynomial by 1+x/xi(q). -
2543       rxi = one/xi
2544       do 180 iback = 1, nq
2545         i = (l + 1) - iback
2546  180    em(i) = em(i) + em(i-1)*rxi
2547 ! compute integral of polynomial. --------------------------------------
2548       s = one
2549       csum = zero
2550       do 190 i = 1, l
2551         csum = csum + s*em(i)/real(i+1)
2552  190    s = -s
2553       tq(3) = flotl*em0/csum
2554       go to 300
2556 ! set coefficients for bdf methods. ------------------------------------
2557  200  do 210 i = 3, l
2558  210    el(i) = zero
2559       el(1) = one
2560       el(2) = one
2561       alph0 = -one
2562       ahatn0 = -one
2563       hsum = h
2564       rxi = one
2565       rxis = one
2566       if (nq .eq. 1) go to 240
2567       do 230 j = 1, nqm2
2568 ! in el, construct coefficients of (1+x/xi(1))*...*(1+x/xi(j+1)). ------
2569         hsum = hsum + tau(j)
2570         rxi = h/hsum
2571         jp1 = j + 1
2572         alph0 = alph0 - one/real(jp1)
2573         do 220 iback = 1, jp1
2574           i = (j + 3) - iback
2575  220      el(i) = el(i) + el(i-1)*rxi
2576  230    continue
2577       alph0 = alph0 - one/real(nq)
2578       rxis = -el(2) - alph0
2579       hsum = hsum + tau(nqm1)
2580       rxi = h/hsum
2581       ahatn0 = -el(2) - rxi
2582       do 235 iback = 1, nq
2583         i = (nq + 2) - iback
2584  235    el(i) = el(i) + el(i-1)*rxis
2585  240  t1 = one - ahatn0 + alph0
2586       t2 = one + real(nq)*t1
2587       tq(2) = abs(alph0*t2/t1)
2588       tq(5) = abs(t2/(el(l)*rxi/rxis))
2589       if (nqwait .ne. 1) go to 300
2590       cnqm1 = rxis/el(l)
2591       t3 = alph0 + one/real(nq)
2592       t4 = ahatn0 + rxi
2593       elp = t3/(one - t4 + t3)
2594       tq(1) = abs(elp/cnqm1)
2595       hsum = hsum + tau(nq)
2596       rxi = h/hsum
2597       t5 = alph0 - one/real(nq+1)
2598       t6 = ahatn0 - rxi
2599       elp = t2/(one - t6 + t5)
2600       tq(3) = abs(elp*rxi*(flotl + one)*t5)
2601  300  tq(4) = cortes*tq(2)
2602       return
2603 !----------------------- end of subroutine svset -----------------------
2604       end subroutine svset
2605 !*deck svjust
2606       subroutine svjust (yh, ldyh, iord)
2607       real yh
2608       integer ldyh, iord
2609       dimension yh(ldyh,*)
2610 !-----------------------------------------------------------------------
2611 ! call sequence input -- yh, ldyh, iord
2612 ! call sequence output -- yh
2613 ! common block input -- nq, meth, lmax, hscal, tau(13), n
2614 ! common block variables accessed..
2615 !     /svode_cmn_01/ -- hscal, tau(13), lmax, meth, n, nq,
2617 ! subroutines called by svjust.. saxpy
2618 ! function routines called by svjust.. none
2619 !-----------------------------------------------------------------------
2620 ! this subroutine adjusts the yh array on reduction of order,
2621 ! and also when the order is increased for the stiff option (meth = 2).
2622 ! communication with svjust uses the following..
2623 ! iord  = an integer flag used when meth = 2 to indicate an order
2624 !         increase (iord = +1) or an order decrease (iord = -1).
2625 ! hscal = step size h used in scaling of nordsieck array yh.
2626 !         (if iord = +1, svjust assumes that hscal = tau(1).)
2627 ! see references 1 and 2 for details.
2628 !-----------------------------------------------------------------------
2630 ! type declarations for labeled common block svode_cmn_01 --------------------
2632       real acnrm, ccmxj, conp, crate, drc, el,   &
2633            eta, etamax, h, hmin, hmxi, hnew, hscal, prl1,   &
2634            rc, rl1, tau, tq, tn, uround
2635       integer icf, init, ipup, jcur, jstart, jsv, kflag, kuth,   &
2636               l, lmax, lyh, lewt, lacor, lsavf, lwm, liwm,   &
2637               locjs, maxord, meth, miter, msbj, mxhnil, mxstep,   &
2638               n, newh, newq, nhnil, nq, nqnyh, nqwait, nslj,   &
2639               nslp, nyh
2641 ! type declarations for local variables --------------------------------
2643       real alph0, alph1, hsum, one, prod, t1, xi,xiold, zero
2644       integer i, iback, j, jp1, lp1, nqm1, nqm2, nqp1
2645 !-----------------------------------------------------------------------
2646 ! the following fortran-77 declaration is to cause the values of the
2647 ! listed (local) variables to be saved between calls to this integrator.
2648 !-----------------------------------------------------------------------
2649       save one, zero
2651       common /svode_cmn_01/ acnrm, ccmxj, conp, crate, drc, el(13),   &
2652                       eta, etamax, h, hmin, hmxi, hnew, hscal, prl1,   &
2653                       rc, rl1, tau(13), tq(5), tn, uround,   &
2654                       icf, init, ipup, jcur, jstart, jsv, kflag, kuth,   &
2655                       l, lmax, lyh, lewt, lacor, lsavf, lwm, liwm,   &
2656                       locjs, maxord, meth, miter, msbj, mxhnil, mxstep,   &
2657                       n, newh, newq, nhnil, nq, nqnyh, nqwait, nslj,   &
2658                       nslp, nyh
2660       data one /1.0e0/, zero /0.0e0/
2662       if ((nq .eq. 2) .and. (iord .ne. 1)) return
2663       nqm1 = nq - 1
2664       nqm2 = nq - 2
2665       go to (100, 200), meth
2666 !-----------------------------------------------------------------------
2667 ! nonstiff option...
2668 ! check to see if the order is being increased or decreased.
2669 !-----------------------------------------------------------------------
2670  100  continue
2671       if (iord .eq. 1) go to 180
2672 ! order decrease. ------------------------------------------------------
2673       do 110 j = 1, lmax
2674  110    el(j) = zero
2675       el(2) = one
2676       hsum = zero
2677       do 130 j = 1, nqm2
2678 ! construct coefficients of x*(x+xi(1))*...*(x+xi(j)). -----------------
2679         hsum = hsum + tau(j)
2680         xi = hsum/hscal
2681         jp1 = j + 1
2682         do 120 iback = 1, jp1
2683           i = (j + 3) - iback
2684  120      el(i) = el(i)*xi + el(i-1)
2685  130    continue
2686 ! construct coefficients of integrated polynomial. ---------------------
2687       do 140 j = 2, nqm1
2688  140    el(j+1) = real(nq)*el(j)/real(j)
2689 ! subtract correction terms from yh array. -----------------------------
2690       do 170 j = 3, nq
2691         do 160 i = 1, n
2692  160      yh(i,j) = yh(i,j) - yh(i,l)*el(j)
2693  170    continue
2694       return
2695 ! order increase. ------------------------------------------------------
2696 ! zero out next column in yh array. ------------------------------------
2697  180  continue
2698       lp1 = l + 1
2699       do 190 i = 1, n
2700  190    yh(i,lp1) = zero
2701       return
2702 !-----------------------------------------------------------------------
2703 ! stiff option...
2704 ! check to see if the order is being increased or decreased.
2705 !-----------------------------------------------------------------------
2706  200  continue
2707       if (iord .eq. 1) go to 300
2708 ! order decrease. ------------------------------------------------------
2709       do 210 j = 1, lmax
2710  210    el(j) = zero
2711       el(3) = one
2712       hsum = zero
2713       do 230 j = 1,nqm2
2714 ! construct coefficients of x*x*(x+xi(1))*...*(x+xi(j)). ---------------
2715         hsum = hsum + tau(j)
2716         xi = hsum/hscal
2717         jp1 = j + 1
2718         do 220 iback = 1, jp1
2719           i = (j + 4) - iback
2720  220      el(i) = el(i)*xi + el(i-1)
2721  230    continue
2722 ! subtract correction terms from yh array. -----------------------------
2723       do 250 j = 3,nq
2724         do 240 i = 1, n
2725  240      yh(i,j) = yh(i,j) - yh(i,l)*el(j)
2726  250    continue
2727       return
2728 ! order increase. ------------------------------------------------------
2729  300  do 310 j = 1, lmax
2730  310    el(j) = zero
2731       el(3) = one
2732       alph0 = -one
2733       alph1 = one
2734       prod = one
2735       xiold = one
2736       hsum = hscal
2737       if (nq .eq. 1) go to 340
2738       do 330 j = 1, nqm1
2739 ! construct coefficients of x*x*(x+xi(1))*...*(x+xi(j)). ---------------
2740         jp1 = j + 1
2741         hsum = hsum + tau(jp1)
2742         xi = hsum/hscal
2743         prod = prod*xi
2744         alph0 = alph0 - one/real(jp1)
2745         alph1 = alph1 + one/xi
2746         do 320 iback = 1, jp1
2747           i = (j + 4) - iback
2748  320      el(i) = el(i)*xiold + el(i-1)
2749         xiold = xi
2750  330    continue
2751  340  continue
2752       t1 = (-alph0 - alph1)/prod
2753 ! load column l + 1 in yh array. ---------------------------------------
2754       lp1 = l + 1
2755       do 350 i = 1, n
2756  350    yh(i,lp1) = t1*yh(i,lmax)
2757 ! add correction terms to yh array. ------------------------------------
2758       nqp1 = nq + 1
2759       do 370 j = 3, nqp1
2760         call saxpy (n, el(j), yh(1,lp1), 1, yh(1,j), 1 )
2761  370  continue
2762       return
2763 !----------------------- end of subroutine svjust ----------------------
2764       end subroutine svjust                 
2765 !*deck svnlsd
2766       subroutine svnlsd (y, yh, ldyh, vsav, savf, ewt, acor, iwm, wm,   &
2767                        f, jac, pdum, nflag, rpar, ipar)
2768       external f, jac, pdum
2769       real y, yh, vsav, savf, ewt, acor, wm, rpar
2770       integer ldyh, iwm, nflag, ipar
2771       dimension y(*), yh(ldyh,*), vsav(*), savf(*), ewt(*), acor(*),   &
2772                 iwm(*), wm(*), rpar(*), ipar(*)
2773 !-----------------------------------------------------------------------
2774 ! call sequence input -- y, yh, ldyh, savf, ewt, acor, iwm, wm,
2775 !                        f, jac, nflag, rpar, ipar
2776 ! call sequence output -- yh, acor, wm, iwm, nflag
2777 ! common block variables accessed..
2778 !     /svode_cmn_01/ acnrm, crate, drc, h, rc, rl1, tq(5), tn, icf,
2779 !                jcur, meth, miter, n, nslp
2780 !     /svode_cmn_02/ hu, ncfn, netf, nfe, nje, nlu, nni, nqu, nst
2782 ! subroutines called by svnlsd.. f, saxpy, scopy, sscal, svjac, svsol
2783 ! function routines called by svnlsd.. svnorm
2784 !-----------------------------------------------------------------------
2785 ! subroutine svnlsd is a nonlinear system solver, which uses functional
2786 ! iteration or a chord (modified newton) method.  for the chord method
2787 ! direct linear algebraic system solvers are used.  subroutine svnlsd
2788 ! then handles the corrector phase of this integration package.
2790 ! communication with svnlsd is done with the following variables. (for
2791 ! more details, please see the comments in the driver subroutine.)
2793 ! y          = the dependent variable, a vector of length n, input.
2794 ! yh         = the nordsieck (taylor) array, ldyh by lmax, input
2795 !              and output.  on input, it contains predicted values.
2796 ! ldyh       = a constant .ge. n, the first dimension of yh, input.
2797 ! vsav       = unused work array.
2798 ! savf       = a work array of length n.
2799 ! ewt        = an error weight vector of length n, input.
2800 ! acor       = a work array of length n, used for the accumulated
2801 !              corrections to the predicted y vector.
2802 ! wm,iwm     = real and integer work arrays associated with matrix
2803 !              operations in chord iteration (miter .ne. 0).
2804 ! f          = dummy name for user supplied routine for f.
2805 ! jac        = dummy name for user supplied jacobian routine.
2806 ! pdum       = unused dummy subroutine name.  included for uniformity
2807 !              over collection of integrators.
2808 ! nflag      = input/output flag, with values and meanings as follows..
2809 !              input
2810 !                  0 first call for this time step.
2811 !                 -1 convergence failure in previous call to svnlsd.
2812 !                 -2 error test failure in svstep.
2813 !              output
2814 !                  0 successful completion of nonlinear solver.
2815 !                 -1 convergence failure or singular matrix.
2816 !                 -2 unrecoverable error in matrix preprocessing
2817 !                    (cannot occur here).
2818 !                 -3 unrecoverable error in solution (cannot occur
2819 !                    here).
2820 ! rpar, ipar = dummy names for user's real and integer work arrays.
2822 ! ipup       = own variable flag with values and meanings as follows..
2823 !              0,            do not update the newton matrix.
2824 !              miter .ne. 0, update newton matrix, because it is the
2825 !                            initial step, order was changed, the error
2826 !                            test failed, or an update is indicated by
2827 !                            the scalar rc or step counter nst.
2829 ! for more details, see comments in driver subroutine.
2830 !-----------------------------------------------------------------------
2831 ! type declarations for labeled common block svode_cmn_01 --------------------
2833       real acnrm, ccmxj, conp, crate, drc, el,   &
2834            eta, etamax, h, hmin, hmxi, hnew, hscal, prl1,   &
2835            rc, rl1, tau, tq, tn, uround
2836       integer icf, init, ipup, jcur, jstart, jsv, kflag, kuth,   &
2837               l, lmax, lyh, lewt, lacor, lsavf, lwm, liwm,   &
2838               locjs, maxord, meth, miter, msbj, mxhnil, mxstep,   &
2839               n, newh, newq, nhnil, nq, nqnyh, nqwait, nslj,   &
2840               nslp, nyh
2842 ! type declarations for labeled common block svode_cmn_02 --------------------
2844       real hu
2845       integer ncfn, netf, nfe, nje, nlu, nni, nqu, nst
2847 ! type declarations for local variables --------------------------------
2849       real ccmax, crdown, cscale, dcon, del, delp, one,   &
2850            rdiv, two, zero
2851       integer i, ierpj, iersl, m, maxcor, msbp
2853 ! type declaration for function subroutines called ---------------------
2855 !     real svnorm      ! rce 2005-jan-21 - module conversion
2856 !-----------------------------------------------------------------------
2857 ! the following fortran-77 declaration is to cause the values of the
2858 ! listed (local) variables to be saved between calls to this integrator.
2859 !-----------------------------------------------------------------------
2860       save ccmax, crdown, maxcor, msbp, rdiv, one, two, zero
2862       common /svode_cmn_01/ acnrm, ccmxj, conp, crate, drc, el(13),   &
2863                       eta, etamax, h, hmin, hmxi, hnew, hscal, prl1,   &
2864                       rc, rl1, tau(13), tq(5), tn, uround,   &
2865                       icf, init, ipup, jcur, jstart, jsv, kflag, kuth,   &
2866                       l, lmax, lyh, lewt, lacor, lsavf, lwm, liwm,   &
2867                       locjs, maxord, meth, miter, msbj, mxhnil, mxstep,   &
2868                       n, newh, newq, nhnil, nq, nqnyh, nqwait, nslj,   &
2869                       nslp, nyh
2870       common /svode_cmn_02/ hu, ncfn, netf, nfe, nje, nlu, nni, nqu, nst
2872       data ccmax /0.3e0/, crdown /0.3e0/, maxcor /3/, msbp /20/,   &
2873            rdiv  /2.0e0/
2874       data one /1.0e0/, two /2.0e0/, zero /0.0e0/
2875 !-----------------------------------------------------------------------
2876 ! on the first step, on a change of method order, or after a
2877 ! nonlinear convergence failure with nflag = -2, set ipup = miter
2878 ! to force a jacobian update when miter .ne. 0.
2879 !-----------------------------------------------------------------------
2880       if (jstart .eq. 0) nslp = 0
2881       if (nflag .eq. 0) icf = 0
2882       if (nflag .eq. -2) ipup = miter
2883       if ( (jstart .eq. 0) .or. (jstart .eq. -1) ) ipup = miter
2884 ! if this is functional iteration, set crate .eq. 1 and drop to 220
2885       if (miter .eq. 0) then
2886         crate = one
2887         go to 220
2888       endif
2889 !-----------------------------------------------------------------------
2890 ! rc is the ratio of new to old values of the coefficient h/el(2)=h/l1.
2891 ! when rc differs from 1 by more than ccmax, ipup is set to miter
2892 ! to force svjac to be called, if a jacobian is involved.
2893 ! in any case, svjac is called at least every msbp steps.
2894 !-----------------------------------------------------------------------
2895       drc = abs(rc-one)
2896       if (drc .gt. ccmax .or. nst .ge. nslp+msbp) ipup = miter
2897 !-----------------------------------------------------------------------
2898 ! up to maxcor corrector iterations are taken.  a convergence test is
2899 ! made on the r.m.s. norm of each correction, weighted by the error
2900 ! weight vector ewt.  the sum of the corrections is accumulated in the
2901 ! vector acor(i).  the yh array is not altered in the corrector loop.
2902 !-----------------------------------------------------------------------
2903  220  m = 0
2904       delp = zero
2905       call scopy (n, yh(1,1), 1, y, 1 )
2906       call f (n, tn, y, savf, rpar, ipar)
2907       nfe = nfe + 1
2908       if (ipup .le. 0) go to 250
2909 !-----------------------------------------------------------------------
2910 ! if indicated, the matrix p = i - h*rl1*j is reevaluated and
2911 ! preprocessed before starting the corrector iteration.  ipup is set
2912 ! to 0 as an indicator that this has been done.
2913 !-----------------------------------------------------------------------
2914       call svjac (y, yh, ldyh, ewt, acor, savf, wm, iwm, f, jac, ierpj,   &
2915                  rpar, ipar)
2916       ipup = 0
2917       rc = one
2918       drc = zero
2919       crate = one
2920       nslp = nst
2921 ! if matrix is singular, take error return to force cut in step size. --
2922       if (ierpj .ne. 0) go to 430
2923  250  do 260 i = 1,n
2924  260    acor(i) = zero
2925 ! this is a looping point for the corrector iteration. -----------------
2926  270  if (miter .ne. 0) go to 350
2927 !-----------------------------------------------------------------------
2928 ! in the case of functional iteration, update y directly from
2929 ! the result of the last function evaluation.
2930 !-----------------------------------------------------------------------
2931       do 280 i = 1,n
2932  280    savf(i) = rl1*(h*savf(i) - yh(i,2))
2933       do 290 i = 1,n
2934  290    y(i) = savf(i) - acor(i)
2935       del = svnorm (n, y, ewt)
2936       do 300 i = 1,n
2937  300    y(i) = yh(i,1) + savf(i)
2938       call scopy (n, savf, 1, acor, 1)
2939       go to 400
2940 !-----------------------------------------------------------------------
2941 ! in the case of the chord method, compute the corrector error,
2942 ! and solve the linear system with that as right-hand side and
2943 ! p as coefficient matrix.  the correction is scaled by the factor
2944 ! 2/(1+rc) to account for changes in h*rl1 since the last svjac call.
2945 !-----------------------------------------------------------------------
2946  350  do 360 i = 1,n
2947  360    y(i) = (rl1*h)*savf(i) - (rl1*yh(i,2) + acor(i))
2948       call svsol (wm, iwm, y, iersl)
2949       nni = nni + 1
2950       if (iersl .gt. 0) go to 410
2951       if (meth .eq. 2 .and. rc .ne. one) then
2952         cscale = two/(one + rc)
2953         call sscal (n, cscale, y, 1)
2954       endif
2955       del = svnorm (n, y, ewt)
2956       call saxpy (n, one, y, 1, acor, 1)
2957       do 380 i = 1,n
2958  380    y(i) = yh(i,1) + acor(i)
2959 !-----------------------------------------------------------------------
2960 ! test for convergence.  if m .gt. 0, an estimate of the convergence
2961 ! rate constant is stored in crate, and this is used in the test.
2962 !-----------------------------------------------------------------------
2963  400  if (m .ne. 0) crate = max(crdown*crate,del/delp)
2964       dcon = del*min(one,crate)/tq(4)
2965       if (dcon .le. one) go to 450
2966       m = m + 1
2967       if (m .eq. maxcor) go to 410
2968       if (m .ge. 2 .and. del .gt. rdiv*delp) go to 410
2969       delp = del
2970       call f (n, tn, y, savf, rpar, ipar)
2971       nfe = nfe + 1
2972       go to 270
2974  410  if (miter .eq. 0 .or. jcur .eq. 1) go to 430
2975       icf = 1
2976       ipup = miter
2977       go to 220
2979  430  continue
2980       nflag = -1
2981       icf = 2
2982       ipup = miter
2983       return
2985 ! return for successful step. ------------------------------------------
2986  450  nflag = 0
2987       jcur = 0
2988       icf = 0
2989       if (m .eq. 0) acnrm = del
2990       if (m .gt. 0) acnrm = svnorm (n, acor, ewt)
2991       return
2992 !----------------------- end of subroutine svnlsd ----------------------
2993       end subroutine svnlsd                                              
2994 !*deck svjac
2995       subroutine svjac (y, yh, ldyh, ewt, ftem, savf, wm, iwm, f, jac,   &
2996                        ierpj, rpar, ipar)
2997       external f, jac
2998       real y, yh, ewt, ftem, savf, wm, rpar
2999       integer ldyh, iwm, ierpj, ipar
3000       dimension y(*), yh(ldyh,*), ewt(*), ftem(*), savf(*),   &
3001          wm(*), iwm(*), rpar(*), ipar(*)
3002 !-----------------------------------------------------------------------
3003 ! call sequence input -- y, yh, ldyh, ewt, ftem, savf, wm, iwm,
3004 !                        f, jac, rpar, ipar
3005 ! call sequence output -- wm, iwm, ierpj
3006 ! common block variables accessed..
3007 !     /svode_cmn_01/  ccmxj, drc, h, rl1, tn, uround, icf, jcur, locjs,
3008 !               miter, msbj, n, nslj
3009 !     /svode_cmn_02/  nfe, nst, nje, nlu
3011 ! subroutines called by svjac.. f, jac, sacopy, scopy, sgbfa, sgefa,
3012 !                              sscal
3013 ! function routines called by svjac.. svnorm
3014 !-----------------------------------------------------------------------
3015 ! svjac is called by svnlsd to compute and process the matrix
3016 ! p = i - h*rl1*j , where j is an approximation to the jacobian.
3017 ! here j is computed by the user-supplied routine jac if
3018 ! miter = 1 or 4, or by finite differencing if miter = 2, 3, or 5.
3019 ! if miter = 3, a diagonal approximation to j is used.
3020 ! if jsv = -1, j is computed from scratch in all cases.
3021 ! if jsv = 1 and miter = 1, 2, 4, or 5, and if the saved value of j is
3022 ! considered acceptable, then p is constructed from the saved j.
3023 ! j is stored in wm and replaced by p.  if miter .ne. 3, p is then
3024 ! subjected to lu decomposition in preparation for later solution
3025 ! of linear systems with p as coefficient matrix. this is done
3026 ! by sgefa if miter = 1 or 2, and by sgbfa if miter = 4 or 5.
3028 ! communication with svjac is done with the following variables.  (for
3029 ! more details, please see the comments in the driver subroutine.)
3030 ! y          = vector containing predicted values on entry.
3031 ! yh         = the nordsieck array, an ldyh by lmax array, input.
3032 ! ldyh       = a constant .ge. n, the first dimension of yh, input.
3033 ! ewt        = an error weight vector of length n.
3034 ! savf       = array containing f evaluated at predicted y, input.
3035 ! wm         = real work space for matrices.  in the output, it contains
3036 !              the inverse diagonal matrix if miter = 3 and the lu
3037 !              decomposition of p if miter is 1, 2 , 4, or 5.
3038 !              storage of matrix elements starts at wm(3).
3039 !              storage of the saved jacobian starts at wm(locjs).
3040 !              wm also contains the following matrix-related data..
3041 !              wm(1) = sqrt(uround), used in numerical jacobian step.
3042 !              wm(2) = h*rl1, saved for later use if miter = 3.
3043 ! iwm        = integer work space containing pivot information,
3044 !              starting at iwm(31), if miter is 1, 2, 4, or 5.
3045 !              iwm also contains band parameters ml = iwm(1) and
3046 !              mu = iwm(2) if miter is 4 or 5.
3047 ! f          = dummy name for the user supplied subroutine for f.
3048 ! jac        = dummy name for the user supplied jacobian subroutine.
3049 ! rpar, ipar = dummy names for user's real and integer work arrays.
3050 ! rl1        = 1/el(2) (input).
3051 ! ierpj      = output error flag,  = 0 if no trouble, 1 if the p
3052 !              matrix is found to be singular.
3053 ! jcur       = output flag to indicate whether the jacobian matrix
3054 !              (or approximation) is now current.
3055 !              jcur = 0 means j is not current.
3056 !              jcur = 1 means j is current.
3057 !-----------------------------------------------------------------------
3059 ! type declarations for labeled common block svode_cmn_01 --------------------
3061       real acnrm, ccmxj, conp, crate, drc, el,   &
3062            eta, etamax, h, hmin, hmxi, hnew, hscal, prl1,   &
3063            rc, rl1, tau, tq, tn, uround
3064       integer icf, init, ipup, jcur, jstart, jsv, kflag, kuth,   &
3065               l, lmax, lyh, lewt, lacor, lsavf, lwm, liwm,   &
3066               locjs, maxord, meth, miter, msbj, mxhnil, mxstep,   &
3067               n, newh, newq, nhnil, nq, nqnyh, nqwait, nslj,   &
3068               nslp, nyh
3070 ! type declarations for labeled common block svode_cmn_02 --------------------
3072       real hu
3073       integer ncfn, netf, nfe, nje, nlu, nni, nqu, nst
3075 ! type declarations for local variables --------------------------------
3077       real con, di, fac, hrl1, one, pt1, r, r0, srur, thou,   &
3078            yi, yj, yjj, zero
3079       integer i, i1, i2, ier, ii, j, j1, jj, jok, lenp, mba, mband,   &
3080               meb1, meband, ml, ml3, mu, np1
3082 ! type declaration for function subroutines called ---------------------
3084 !     real svnorm      ! rce 2005-jan-21 - module conversion
3085 !-----------------------------------------------------------------------
3086 ! the following fortran-77 declaration is to cause the values of the
3087 ! listed (local) variables to be saved between calls to this subroutine.
3088 !-----------------------------------------------------------------------
3089       save one, pt1, thou, zero
3090 !-----------------------------------------------------------------------
3091       common /svode_cmn_01/ acnrm, ccmxj, conp, crate, drc, el(13),   &
3092                       eta, etamax, h, hmin, hmxi, hnew, hscal, prl1,   &
3093                       rc, rl1, tau(13), tq(5), tn, uround,   &
3094                       icf, init, ipup, jcur, jstart, jsv, kflag, kuth,   &
3095                       l, lmax, lyh, lewt, lacor, lsavf, lwm, liwm,   &
3096                       locjs, maxord, meth, miter, msbj, mxhnil, mxstep,   &
3097                       n, newh, newq, nhnil, nq, nqnyh, nqwait, nslj,   &
3098                       nslp, nyh
3099       common /svode_cmn_02/ hu, ncfn, netf, nfe, nje, nlu, nni, nqu, nst
3101       data one /1.0e0/, thou /1000.0e0/, zero /0.0e0/, pt1 /0.1e0/
3103       ierpj = 0
3104       hrl1 = h*rl1
3105 ! see whether j should be evaluated (jok = -1) or not (jok = 1). -------
3106       jok = jsv
3107       if (jsv .eq. 1) then
3108         if (nst .eq. 0 .or. nst .gt. nslj+msbj) jok = -1
3109         if (icf .eq. 1 .and. drc .lt. ccmxj) jok = -1
3110         if (icf .eq. 2) jok = -1
3111       endif
3112 ! end of setting jok. --------------------------------------------------
3114       if (jok .eq. -1 .and. miter .eq. 1) then
3115 ! if jok = -1 and miter = 1, call jac to evaluate jacobian. ------------
3116       nje = nje + 1
3117       nslj = nst
3118       jcur = 1
3119       lenp = n*n
3120       do 110 i = 1,lenp
3121  110    wm(i+2) = zero
3122       call jac (n, tn, y, 0, 0, wm(3), n, rpar, ipar)
3123       if (jsv .eq. 1) call scopy (lenp, wm(3), 1, wm(locjs), 1)
3124       endif
3126       if (jok .eq. -1 .and. miter .eq. 2) then
3127 ! if miter = 2, make n calls to f to approximate the jacobian. ---------
3128       nje = nje + 1
3129       nslj = nst
3130       jcur = 1
3131       fac = svnorm (n, savf, ewt)
3132       r0 = thou*abs(h)*uround*real(n)*fac
3133       if (r0 .eq. zero) r0 = one
3134       srur = wm(1)
3135       j1 = 2
3136       do 230 j = 1,n
3137         yj = y(j)
3138         r = max(srur*abs(yj),r0/ewt(j))
3139         y(j) = y(j) + r
3140         fac = one/r
3141         call f (n, tn, y, ftem, rpar, ipar)
3142         do 220 i = 1,n
3143  220      wm(i+j1) = (ftem(i) - savf(i))*fac
3144         y(j) = yj
3145         j1 = j1 + n
3146  230    continue
3147       nfe = nfe + n
3148       lenp = n*n
3149       if (jsv .eq. 1) call scopy (lenp, wm(3), 1, wm(locjs), 1)
3150       endif
3152       if (jok .eq. 1 .and. (miter .eq. 1 .or. miter .eq. 2)) then
3153       jcur = 0
3154       lenp = n*n
3155       call scopy (lenp, wm(locjs), 1, wm(3), 1)
3156       endif
3158       if (miter .eq. 1 .or. miter .eq. 2) then
3159 ! multiply jacobian by scalar, add identity, and do lu decomposition. --
3160       con = -hrl1
3161       call sscal (lenp, con, wm(3), 1)
3162       j = 3
3163       np1 = n + 1
3164       do 250 i = 1,n
3165         wm(j) = wm(j) + one
3166  250    j = j + np1
3167       nlu = nlu + 1
3168       call sgefa (wm(3), n, n, iwm(31), ier)
3169       if (ier .ne. 0) ierpj = 1
3170       return
3171       endif
3172 ! end of code block for miter = 1 or 2. --------------------------------
3174       if (miter .eq. 3) then
3175 ! if miter = 3, construct a diagonal approximation to j and p. ---------
3176       nje = nje + 1
3177       jcur = 1
3178       wm(2) = hrl1
3179       r = rl1*pt1
3180       do 310 i = 1,n
3181  310    y(i) = y(i) + r*(h*savf(i) - yh(i,2))
3182       call f (n, tn, y, wm(3), rpar, ipar)
3183       nfe = nfe + 1
3184       do 320 i = 1,n
3185         r0 = h*savf(i) - yh(i,2)
3186         di = pt1*r0 - h*(wm(i+2) - savf(i))
3187         wm(i+2) = one
3188         if (abs(r0) .lt. uround/ewt(i)) go to 320
3189         if (abs(di) .eq. zero) go to 330
3190         wm(i+2) = pt1*r0/di
3191  320    continue
3192       return
3193  330  ierpj = 1
3194       return
3195       endif
3196 ! end of code block for miter = 3. -------------------------------------
3198 ! set constants for miter = 4 or 5. ------------------------------------
3199       ml = iwm(1)
3200       mu = iwm(2)
3201       ml3 = ml + 3
3202       mband = ml + mu + 1
3203       meband = mband + ml
3204       lenp = meband*n
3206       if (jok .eq. -1 .and. miter .eq. 4) then
3207 ! if jok = -1 and miter = 4, call jac to evaluate jacobian. ------------
3208       nje = nje + 1
3209       nslj = nst
3210       jcur = 1
3211       do 410 i = 1,lenp
3212  410    wm(i+2) = zero
3213       call jac (n, tn, y, ml, mu, wm(ml3), meband, rpar, ipar)
3214       if (jsv .eq. 1)   &
3215          call sacopy (mband, n, wm(ml3), meband, wm(locjs), mband)
3216       endif
3218       if (jok .eq. -1 .and. miter .eq. 5) then
3219 ! if miter = 5, make ml+mu+1 calls to f to approximate the jacobian. ---
3220       nje = nje + 1
3221       nslj = nst
3222       jcur = 1
3223       mba = min(mband,n)
3224       meb1 = meband - 1
3225       srur = wm(1)
3226       fac = svnorm (n, savf, ewt)
3227       r0 = thou*abs(h)*uround*real(n)*fac
3228       if (r0 .eq. zero) r0 = one
3229       do 560 j = 1,mba
3230         do 530 i = j,n,mband
3231           yi = y(i)
3232           r = max(srur*abs(yi),r0/ewt(i))
3233  530      y(i) = y(i) + r
3234         call f (n, tn, y, ftem, rpar, ipar)
3235         do 550 jj = j,n,mband
3236           y(jj) = yh(jj,1)
3237           yjj = y(jj)
3238           r = max(srur*abs(yjj),r0/ewt(jj))
3239           fac = one/r
3240           i1 = max(jj-mu,1)
3241           i2 = min(jj+ml,n)
3242           ii = jj*meb1 - ml + 2
3243           do 540 i = i1,i2
3244  540        wm(ii+i) = (ftem(i) - savf(i))*fac
3245  550      continue
3246  560    continue
3247       nfe = nfe + mba
3248       if (jsv .eq. 1)   &
3249          call sacopy (mband, n, wm(ml3), meband, wm(locjs), mband)
3250       endif
3252       if (jok .eq. 1) then
3253       jcur = 0
3254       call sacopy (mband, n, wm(locjs), mband, wm(ml3), meband)
3255       endif
3257 ! multiply jacobian by scalar, add identity, and do lu decomposition.
3258       con = -hrl1
3259       call sscal (lenp, con, wm(3), 1 )
3260       ii = mband + 2
3261       do 580 i = 1,n
3262         wm(ii) = wm(ii) + one
3263  580    ii = ii + meband
3264       nlu = nlu + 1
3265       call sgbfa (wm(3), meband, n, ml, mu, iwm(31), ier)
3266       if (ier .ne. 0) ierpj = 1
3267       return
3268 ! end of code block for miter = 4 or 5. --------------------------------
3270 !----------------------- end of subroutine svjac -----------------------
3271       end subroutine svjac                                                
3272 !*deck sacopy
3273       subroutine sacopy (nrow, ncol, a, nrowa, b, nrowb)
3274       real a, b
3275       integer nrow, ncol, nrowa, nrowb
3276       dimension a(nrowa,ncol), b(nrowb,ncol)
3277 !-----------------------------------------------------------------------
3278 ! call sequence input -- nrow, ncol, a, nrowa, nrowb
3279 ! call sequence output -- b
3280 ! common block variables accessed -- none
3282 ! subroutines called by sacopy.. scopy
3283 ! function routines called by sacopy.. none
3284 !-----------------------------------------------------------------------
3285 ! this routine copies one rectangular array, a, to another, b,
3286 ! where a and b may have different row dimensions, nrowa and nrowb.
3287 ! the data copied consists of nrow rows and ncol columns.
3288 !-----------------------------------------------------------------------
3289       integer ic
3291       do 20 ic = 1,ncol
3292         call scopy (nrow, a(1,ic), 1, b(1,ic), 1)
3293  20     continue
3295       return
3296 !----------------------- end of subroutine sacopy ----------------------
3297       end subroutine sacopy                                 
3298 !*deck svsol
3299       subroutine svsol (wm, iwm, x, iersl)
3300       real wm, x
3301       integer iwm, iersl
3302       dimension wm(*), iwm(*), x(*)
3303 !-----------------------------------------------------------------------
3304 ! call sequence input -- wm, iwm, x
3305 ! call sequence output -- x, iersl
3306 ! common block variables accessed..
3307 !     /svode_cmn_01/ -- h, rl1, miter, n
3309 ! subroutines called by svsol.. sgesl, sgbsl
3310 ! function routines called by svsol.. none
3311 !-----------------------------------------------------------------------
3312 ! this routine manages the solution of the linear system arising from
3313 ! a chord iteration.  it is called if miter .ne. 0.
3314 ! if miter is 1 or 2, it calls sgesl to accomplish this.
3315 ! if miter = 3 it updates the coefficient h*rl1 in the diagonal
3316 ! matrix, and then computes the solution.
3317 ! if miter is 4 or 5, it calls sgbsl.
3318 ! communication with svsol uses the following variables..
3319 ! wm    = real work space containing the inverse diagonal matrix if
3320 !         miter = 3 and the lu decomposition of the matrix otherwise.
3321 !         storage of matrix elements starts at wm(3).
3322 !         wm also contains the following matrix-related data..
3323 !         wm(1) = sqrt(uround) (not used here),
3324 !         wm(2) = hrl1, the previous value of h*rl1, used if miter = 3.
3325 ! iwm   = integer work space containing pivot information, starting at
3326 !         iwm(31), if miter is 1, 2, 4, or 5.  iwm also contains band
3327 !         parameters ml = iwm(1) and mu = iwm(2) if miter is 4 or 5.
3328 ! x     = the right-hand side vector on input, and the solution vector
3329 !         on output, of length n.
3330 ! iersl = output flag.  iersl = 0 if no trouble occurred.
3331 !         iersl = 1 if a singular matrix arose with miter = 3.
3332 !-----------------------------------------------------------------------
3334 ! type declarations for labeled common block svode_cmn_01 --------------------
3336       real acnrm, ccmxj, conp, crate, drc, el,   &
3337            eta, etamax, h, hmin, hmxi, hnew, hscal, prl1,   &
3338            rc, rl1, tau, tq, tn, uround
3339       integer icf, init, ipup, jcur, jstart, jsv, kflag, kuth,   &
3340               l, lmax, lyh, lewt, lacor, lsavf, lwm, liwm,   &
3341               locjs, maxord, meth, miter, msbj, mxhnil, mxstep,   &
3342               n, newh, newq, nhnil, nq, nqnyh, nqwait, nslj,   &
3343               nslp, nyh
3345 ! type declarations for local variables --------------------------------
3347       integer i, meband, ml, mu
3348       real di, hrl1, one, phrl1, r, zero
3349 !-----------------------------------------------------------------------
3350 ! the following fortran-77 declaration is to cause the values of the
3351 ! listed (local) variables to be saved between calls to this integrator.
3352 !-----------------------------------------------------------------------
3353       save one, zero
3355       common /svode_cmn_01/ acnrm, ccmxj, conp, crate, drc, el(13),   &
3356                       eta, etamax, h, hmin, hmxi, hnew, hscal, prl1,   &
3357                       rc, rl1, tau(13), tq(5), tn, uround,   &
3358                       icf, init, ipup, jcur, jstart, jsv, kflag, kuth,   &
3359                       l, lmax, lyh, lewt, lacor, lsavf, lwm, liwm,   &
3360                       locjs, maxord, meth, miter, msbj, mxhnil, mxstep,   &
3361                       n, newh, newq, nhnil, nq, nqnyh, nqwait, nslj,   &
3362                       nslp, nyh
3364       data one /1.0e0/, zero /0.0e0/
3366       iersl = 0
3367       go to (100, 100, 300, 400, 400), miter
3368  100  call sgesl (wm(3), n, n, iwm(31), x, 0)
3369       return
3371  300  phrl1 = wm(2)
3372       hrl1 = h*rl1
3373       wm(2) = hrl1
3374       if (hrl1 .eq. phrl1) go to 330
3375       r = hrl1/phrl1
3376       do 320 i = 1,n
3377         di = one - r*(one - one/wm(i+2))
3378         if (abs(di) .eq. zero) go to 390
3379  320    wm(i+2) = one/di
3381  330  do 340 i = 1,n
3382  340    x(i) = wm(i+2)*x(i)
3383       return
3384  390  iersl = 1
3385       return
3387  400  ml = iwm(1)
3388       mu = iwm(2)
3389       meband = 2*ml + mu + 1
3390       call sgbsl (wm(3), meband, n, ml, mu, iwm(31), x, 0)
3391       return
3392 !----------------------- end of subroutine svsol -----------------------
3393       end subroutine svsol                    
3394 !*deck svsrco
3395       subroutine svsrco (rsav, isav, job)
3396       real rsav
3397       integer isav, job
3398       dimension rsav(*), isav(*)
3399 !-----------------------------------------------------------------------
3400 ! call sequence input -- rsav, isav, job
3401 ! call sequence output -- rsav, isav
3402 ! common block variables accessed -- all of /svode_cmn_01/ and /svode_cmn_02/
3404 ! subroutines/functions called by svsrco.. none
3405 !-----------------------------------------------------------------------
3406 ! this routine saves or restores (depending on job) the contents of the
3407 ! common blocks svode_cmn_01 and svode_cmn_02, which are used internally by svode.
3409 ! rsav = real array of length 49 or more.
3410 ! isav = integer array of length 41 or more.
3411 ! job  = flag indicating to save or restore the common blocks..
3412 !        job  = 1 if common is to be saved (written to rsav/isav).
3413 !        job  = 2 if common is to be restored (read from rsav/isav).
3414 !        a call with job = 2 presumes a prior call with job = 1.
3415 !-----------------------------------------------------------------------
3416       real rvod1, rvod2
3417       integer ivod1, ivod2
3418       integer i, leniv1, leniv2, lenrv1, lenrv2
3419 !-----------------------------------------------------------------------
3420 ! the following fortran-77 declaration is to cause the values of the
3421 ! listed (local) variables to be saved between calls to this integrator.
3422 !-----------------------------------------------------------------------
3423       save lenrv1, leniv1, lenrv2, leniv2
3425       common /svode_cmn_01/ rvod1(48), ivod1(33)
3426       common /svode_cmn_02/ rvod2(1), ivod2(8)
3427       data lenrv1/48/, leniv1/33/, lenrv2/1/, leniv2/8/
3429       if (job .eq. 2) go to 100
3430       do 10 i = 1,lenrv1
3431  10     rsav(i) = rvod1(i)
3432       do 15 i = 1,lenrv2
3433  15     rsav(lenrv1+i) = rvod2(i)
3435       do 20 i = 1,leniv1
3436  20     isav(i) = ivod1(i)
3437       do 25 i = 1,leniv2
3438  25     isav(leniv1+i) = ivod2(i)
3440       return
3442  100  continue
3443       do 110 i = 1,lenrv1
3444  110     rvod1(i) = rsav(i)
3445       do 115 i = 1,lenrv2
3446  115     rvod2(i) = rsav(lenrv1+i)
3448       do 120 i = 1,leniv1
3449  120     ivod1(i) = isav(i)
3450       do 125 i = 1,leniv2
3451  125     ivod2(i) = isav(leniv1+i)
3453       return
3454 !----------------------- end of subroutine svsrco ----------------------
3455       end subroutine svsrco                  
3456 !*deck sewset
3457       subroutine sewset (n, itol, rtol, atol, ycur, ewt)
3458       real rtol, atol, ycur, ewt
3459       integer n, itol
3460       dimension rtol(*), atol(*), ycur(n), ewt(n)
3461 !-----------------------------------------------------------------------
3462 ! call sequence input -- n, itol, rtol, atol, ycur
3463 ! call sequence output -- ewt
3464 ! common block variables accessed -- none
3466 ! subroutines/functions called by sewset.. none
3467 !-----------------------------------------------------------------------
3468 ! this subroutine sets the error weight vector ewt according to
3469 !     ewt(i) = rtol(i)*abs(ycur(i)) + atol(i),  i = 1,...,n,
3470 ! with the subscript on rtol and/or atol possibly replaced by 1 above,
3471 ! depending on the value of itol.
3472 !-----------------------------------------------------------------------
3473       integer i
3475       go to (10, 20, 30, 40), itol
3476  10   continue
3477       do 15 i = 1, n
3478  15     ewt(i) = rtol(1)*abs(ycur(i)) + atol(1)
3479       return
3480  20   continue
3481       do 25 i = 1, n
3482  25     ewt(i) = rtol(1)*abs(ycur(i)) + atol(i)
3483       return
3484  30   continue
3485       do 35 i = 1, n
3486  35     ewt(i) = rtol(i)*abs(ycur(i)) + atol(1)
3487       return
3488  40   continue
3489       do 45 i = 1, n
3490  45     ewt(i) = rtol(i)*abs(ycur(i)) + atol(i)
3491       return
3492 !----------------------- end of subroutine sewset ----------------------
3493       end subroutine sewset                                 
3494 !*deck svnorm
3495       real function svnorm (n, v, w)
3496       real v, w
3497       integer n
3498       dimension v(n), w(n)
3499 !-----------------------------------------------------------------------
3500 ! call sequence input -- n, v, w
3501 ! call sequence output -- none
3502 ! common block variables accessed -- none
3504 ! subroutines/functions called by svnorm.. none
3505 !-----------------------------------------------------------------------
3506 ! this function routine computes the weighted root-mean-square norm
3507 ! of the vector of length n contained in the array v, with weights
3508 ! contained in the array w of length n..
3509 !   svnorm = sqrt( (1/n) * sum( v(i)*w(i) )**2 )
3510 !-----------------------------------------------------------------------
3511       real sum
3512       integer i
3514       sum = 0.0e0
3515       do 10 i = 1, n
3516  10     sum = sum + (v(i)*w(i))**2
3517       svnorm = sqrt(sum/real(n))
3518       return
3519 !----------------------- end of function svnorm ------------------------
3520       end function svnorm          
3521 !*deck r1mach
3522       real function r1mach (idum)
3523       integer idum
3524 !-----------------------------------------------------------------------
3525 ! this routine computes the unit roundoff of the machine.
3526 ! this is defined as the smallest positive machine number
3527 ! u such that  1.0 + u .ne. 1.0
3529 ! subroutines/functions called by r1mach.. none
3530 !-----------------------------------------------------------------------
3531       real u, comp
3532       u = 1.0e0
3533  10   u = u*0.5e0
3534       comp = 1.0e0 + u
3535       if (comp .ne. 1.0e0) go to 10
3536       r1mach = u*2.0e0
3537       return
3538 !----------------------- end of function r1mach ------------------------
3539       end function r1mach       
3540 !*deck xerrwv
3541       subroutine xerrwv (msg, nmes, nerr, level, ni, i1, i2, nr, r1, r2)
3542       real r1, r2
3543       integer nmes, nerr, level, ni, i1, i2, nr
3544 !     character*1 msg(nmes)      ! rce 2005-jan-21 - module conversion
3545       character*(nmes) msg
3546 !-----------------------------------------------------------------------
3547 ! subroutines xerrwv, xsetf, xsetun, and the function routine ixsav,
3548 ! as given here, constitute a simplified version of the slatec error
3549 ! handling package.
3550 ! written by a. c. hindmarsh and p. n. brown at llnl.
3551 ! version of 18 november, 1992.
3552 ! this version is in single precision.
3554 ! all arguments are input arguments.
3556 ! msg    = the message (character array).
3557 ! nmes   = the length of msg (number of characters).
3558 ! nerr   = the error number (not used).
3559 ! level  = the error level..
3560 !          0 or 1 means recoverable (control returns to caller).
3561 !          2 means fatal (run is aborted--see note below).
3562 ! ni     = number of integers (0, 1, or 2) to be printed with message.
3563 ! i1,i2  = integers to be printed, depending on ni.
3564 ! nr     = number of reals (0, 1, or 2) to be printed with message.
3565 ! r1,r2  = reals to be printed, depending on nr.
3567 ! note..  this routine is machine-dependent and specialized for use
3568 ! in limited context, in the following ways..
3569 ! 1. the argument msg is assumed to be of type character, and
3570 !    the message is printed with a format of (1x,80a1).
3571 ! 2. the message is assumed to take only one line.
3572 !    multi-line messages are generated by repeated calls.
3573 ! 3. if level = 2, control passes to the statement   stop
3574 !    to abort the run.  this statement may be machine-dependent.
3575 ! 4. r1 and r2 are assumed to be in single precision and are printed
3576 !    in e21.13 format.
3578 ! for a different default logical unit number, change the data
3579 ! statement in function routine ixsav.
3580 ! for a different run-abort command, change the statement following
3581 ! statement 100 at the end.
3582 !-----------------------------------------------------------------------
3583 ! subroutines called by xerrwv.. none
3584 ! function routine called by xerrwv.. ixsav
3585 !-----------------------------------------------------------------------
3587 !     integer i, lunit, ixsav, mesflg      ! rce 2005-jan-21 - module conversion
3588       integer i, lunit,        mesflg
3590 ! get logical unit number and message print flag. ----------------------
3591       lunit = ixsav (1, 0, .false.)
3592       mesflg = ixsav (2, 0, .false.)
3593       if (mesflg .eq. 0) go to 100
3594 ! write the message. ---------------------------------------------------
3595 !     write (lunit,10) (msg(i),i=1,nmes)      ! rce 2005-jan-21 - module conversion
3596       write (lunit,10)  msg
3597  10   format(1x,80a1)
3598       if (ni .eq. 1) write (lunit, 20) i1
3599  20   format(6x,'in above message,  i1 =',i10)
3600       if (ni .eq. 2) write (lunit, 30) i1,i2
3601  30   format(6x,'in above message,  i1 =',i10,3x,'i2 =',i10)
3602       if (nr .eq. 1) write (lunit, 40) r1
3603  40   format(6x,'in above message,  r1 =',e21.13)
3604       if (nr .eq. 2) write (lunit, 50) r1,r2
3605  50   format(6x,'in above,  r1 =',e21.13,3x,'r2 =',e21.13)
3607 ! rce 2005-may-05 - do not write to unit 2
3609 ! added by sp (writing of error messages in unit 2)
3611 !     write (2,10) (msg(i),i=1,nmes)      ! rce 2005-jan-21 - module conversion
3612 !     write (2,10)  msg
3613 !     if (ni .eq. 1) write (2, 20) i1
3614 !     if (ni .eq. 2) write (2, 30) i1,i2
3615 !     if (nr .eq. 1) write (2, 40) r1
3616 !     if (nr .eq. 2) write (2, 50) r1,r2
3618 ! abort the run if level = 2. ------------------------------------------
3619  100  if (level .ne. 2) return
3620 ! rce 2005-may-05 - use peg_error_fatal
3621 !     stop
3622       call peg_error_fatal( -1,   &
3623           '*** fatal error - module_svode_solver, subr xerrwv' )
3624 !----------------------- end of subroutine xerrwv ----------------------
3625       end subroutine xerrwv                                                 
3626 !*deck xsetun
3627       subroutine xsetun (lun)
3628 !-----------------------------------------------------------------------
3629 ! this routine resets the logical unit number for messages.
3631 ! subroutines called by xsetun.. none
3632 ! function routine called by xsetun.. ixsav
3633 !-----------------------------------------------------------------------
3634 !     integer lun, junk, ixsav      ! rce 2005-jan-21 - module conversion
3635       integer lun, junk
3637       if (lun .gt. 0) junk = ixsav (1,lun,.true.)
3638       return
3639 !----------------------- end of subroutine xsetun ----------------------
3640       end subroutine xsetun      
3641 !*deck xsetf
3642       subroutine xsetf (mflag)
3643 !-----------------------------------------------------------------------
3644 ! this routine resets the print control flag mflag.
3646 ! subroutines called by xsetf.. none
3647 ! function routine called by xsetf.. ixsav
3648 !-----------------------------------------------------------------------
3649 !     integer mflag, junk, ixsav      ! rce 2005-jan-21 - module conversion
3650       integer mflag, junk
3652       if (mflag .eq. 0 .or. mflag .eq. 1) junk = ixsav (2,mflag,.true.)
3653       return
3654 !----------------------- end of subroutine xsetf -----------------------
3655       end subroutine xsetf        
3656 !*deck ixsav
3657       integer function ixsav (ipar, ivalue, iset)
3658       logical iset
3659       integer ipar, ivalue
3660 !-----------------------------------------------------------------------
3661 ! ixsav saves and recalls one of two error message parameters:
3662 !   lunit, the logical unit number to which messages are printed, and
3663 !   mesflg, the message print flag.
3664 ! this is a modification of the slatec library routine j4save.
3666 ! saved local variables..
3667 !  lunit  = logical unit number for messages.
3668 !           the default is 6 (machine-dependent).
3669 !  mesflg = print control flag..
3670 !           1 means print all messages (the default).
3671 !           0 means no printing.
3673 ! on input..
3674 !   ipar   = parameter indicator (1 for lunit, 2 for mesflg).
3675 !   ivalue = the value to be set for the parameter, if iset = .true.
3676 !   iset   = logical flag to indicate whether to read or write.
3677 !            if iset = .true., the parameter will be given
3678 !            the value ivalue.  if iset = .false., the parameter
3679 !            will be unchanged, and ivalue is a dummy argument.
3681 ! on return..
3682 !   ixsav = the (old) value of the parameter.
3684 ! subroutines/functions called by ixsav.. none
3685 !-----------------------------------------------------------------------
3686       integer lunit, mesflg
3687 !-----------------------------------------------------------------------
3688 ! the following fortran-77 declaration is to cause the values of the
3689 ! listed (local) variables to be saved between calls to this routine.
3690 !-----------------------------------------------------------------------
3691       save lunit, mesflg
3692       data lunit/6/, mesflg/1/
3694       if (ipar .eq. 1) then
3695         ixsav = lunit
3696         if (iset) lunit = ivalue
3697         endif
3699       if (ipar .eq. 2) then
3700         ixsav = mesflg
3701         if (iset) mesflg = ivalue
3702         endif
3704       return
3705 !----------------------- end of function ixsav -------------------------
3706       end function ixsav                     
3707 ! routines selected from blas and used by vode
3709       subroutine scopy(n,sx,incx,sy,incy)
3711 !     copies a vector, x, to a vector, y.
3712 !     uses unrolled loops for increments equal to 1.
3713 !     jack dongarra, linpack, 3/11/78.
3714 !     modified 12/3/93, array(1) declarations changed to array(*)
3716       real sx(*),sy(*)
3717       integer i,incx,incy,ix,iy,m,mp1,n
3719       if(n.le.0)return
3720       if(incx.eq.1.and.incy.eq.1)go to 20
3722 !        code for unequal increments or equal increments
3723 !          not equal to 1
3725       ix = 1
3726       iy = 1
3727       if(incx.lt.0)ix = (-n+1)*incx + 1
3728       if(incy.lt.0)iy = (-n+1)*incy + 1
3729       do 10 i = 1,n
3730         sy(iy) = sx(ix)
3731         ix = ix + incx
3732         iy = iy + incy
3733    10 continue
3734       return
3736 !        code for both increments equal to 1
3739 !        clean-up loop
3741    20 m = mod(n,7)
3742       if( m .eq. 0 ) go to 40
3743       do 30 i = 1,m
3744         sy(i) = sx(i)
3745    30 continue
3746       if( n .lt. 7 ) return
3747    40 mp1 = m + 1
3748       do 50 i = mp1,n,7
3749         sy(i) = sx(i)
3750         sy(i + 1) = sx(i + 1)
3751         sy(i + 2) = sx(i + 2)
3752         sy(i + 3) = sx(i + 3)
3753         sy(i + 4) = sx(i + 4)
3754         sy(i + 5) = sx(i + 5)
3755         sy(i + 6) = sx(i + 6)
3756    50 continue
3757       return
3758       end subroutine scopy                   
3763       subroutine sscal(n,sa,sx,incx)
3765 !     scales a vector by a constant.
3766 !     uses unrolled loops for increment equal to 1.
3767 !     jack dongarra, linpack, 3/11/78.
3768 !     modified 3/93 to return if incx .le. 0.
3769 !     modified 12/3/93, array(1) declarations changed to array(*)
3771       real sa,sx(*)
3772       integer i,incx,m,mp1,n,nincx
3774       if( n.le.0 .or. incx.le.0 )return
3775       if(incx.eq.1)go to 20
3777 !        code for increment not equal to 1
3779       nincx = n*incx
3780       do 10 i = 1,nincx,incx
3781         sx(i) = sa*sx(i)
3782    10 continue
3783       return
3785 !        code for increment equal to 1
3788 !        clean-up loop
3790    20 m = mod(n,5)
3791       if( m .eq. 0 ) go to 40
3792       do 30 i = 1,m
3793         sx(i) = sa*sx(i)
3794    30 continue
3795       if( n .lt. 5 ) return
3796    40 mp1 = m + 1
3797       do 50 i = mp1,n,5
3798         sx(i) = sa*sx(i)
3799         sx(i + 1) = sa*sx(i + 1)
3800         sx(i + 2) = sa*sx(i + 2)
3801         sx(i + 3) = sa*sx(i + 3)
3802         sx(i + 4) = sa*sx(i + 4)
3803    50 continue
3804       return
3805       end subroutine sscal              
3809       subroutine saxpy(n,sa,sx,incx,sy,incy)
3811 !     constant times a vector plus a vector.
3812 !     uses unrolled loop for increments equal to one.
3813 !     jack dongarra, linpack, 3/11/78.
3814 !     modified 12/3/93, array(1) declarations changed to array(*)
3816       real sx(*),sy(*),sa
3817       integer i,incx,incy,ix,iy,m,mp1,n
3819       if(n.le.0)return
3820       if (sa .eq. 0.0) return
3821       if(incx.eq.1.and.incy.eq.1)go to 20
3823 !        code for unequal increments or equal increments
3824 !          not equal to 1
3826       ix = 1
3827       iy = 1
3828       if(incx.lt.0)ix = (-n+1)*incx + 1
3829       if(incy.lt.0)iy = (-n+1)*incy + 1
3830       do 10 i = 1,n
3831         sy(iy) = sy(iy) + sa*sx(ix)
3832         ix = ix + incx
3833         iy = iy + incy
3834    10 continue
3835       return
3837 !        code for both increments equal to 1
3840 !        clean-up loop
3842    20 m = mod(n,4)
3843       if( m .eq. 0 ) go to 40
3844       do 30 i = 1,m
3845         sy(i) = sy(i) + sa*sx(i)
3846    30 continue
3847       if( n .lt. 4 ) return
3848    40 mp1 = m + 1
3849       do 50 i = mp1,n,4
3850         sy(i) = sy(i) + sa*sx(i)
3851         sy(i + 1) = sy(i + 1) + sa*sx(i + 1)
3852         sy(i + 2) = sy(i + 2) + sa*sx(i + 2)
3853         sy(i + 3) = sy(i + 3) + sa*sx(i + 3)
3854    50 continue
3855       return
3856       end subroutine saxpy                      
3859 ! linpack routines used by svode
3862 !* ======================================================================
3863 !* nist guide to available math software.
3864 !* fullsource for module sgefa from package linpack.
3865 !* retrieved from netlib on wed jun 17 08:00:19 1998.
3866 !* ======================================================================
3867       subroutine sgefa(a,lda,n,ipvt,info)
3868       integer lda,n,ipvt(*),info
3869       real a(lda,*)
3871 !     sgefa factors a real matrix by gaussian elimination.
3873 !     sgefa is usually called by sgeco, but it can be called
3874 !     directly with a saving in time if  rcond  is not needed.
3875 !     (time for sgeco) = (1 + 9/n)*(time for sgefa) .
3877 !     on entry
3879 !        a       real(lda, n)
3880 !                the matrix to be factored.
3882 !        lda     integer
3883 !                the leading dimension of the array  a .
3885 !        n       integer
3886 !                the order of the matrix  a .
3888 !     on return
3890 !        a       an upper triangular matrix and the multipliers
3891 !                which were used to obtain it.
3892 !                the factorization can be written  a = l*u  where
3893 !                l  is a product of permutation and unit lower
3894 !                triangular matrices and  u  is upper triangular.
3896 !        ipvt    integer(n)
3897 !                an integer vector of pivot indices.
3899 !        info    integer
3900 !                = 0  normal value.
3901 !                = k  if  u(k,k) .eq. 0.0 .  this is not an error
3902 !                     condition for this subroutine, but it does
3903 !                     indicate that sgesl or sgedi will divide by zero
3904 !                     if called.  use  rcond  in sgeco for a reliable
3905 !                     indication of singularity.
3907 !     linpack. this version dated 08/14/78 .
3908 !     cleve moler, university of new mexico, argonne national lab.
3910 !     subroutines and functions
3912 !     blas saxpy,sscal,isamax
3914 !     internal variables
3916       real t
3917 !     integer isamax,j,k,kp1,l,nm1      ! rce 2005-jan-21 - module conversion
3918       integer        j,k,kp1,l,nm1
3921 !     gaussian elimination with partial pivoting
3923       info = 0
3924       nm1 = n - 1
3925       if (nm1 .lt. 1) go to 70
3926       do 60 k = 1, nm1
3927          kp1 = k + 1
3929 !        find l = pivot index
3931          l = isamax(n-k+1,a(k,k),1) + k - 1
3932          ipvt(k) = l
3934 !        zero pivot implies this column already triangularized
3936          if (a(l,k) .eq. 0.0e0) go to 40
3938 !           interchange if necessary
3940             if (l .eq. k) go to 10
3941                t = a(l,k)
3942                a(l,k) = a(k,k)
3943                a(k,k) = t
3944    10       continue
3946 !           compute multipliers
3948             t = -1.0e0/a(k,k)
3949             call sscal(n-k,t,a(k+1,k),1)
3951 !           row elimination with column indexing
3953             do 30 j = kp1, n
3954                t = a(l,j)
3955                if (l .eq. k) go to 20
3956                   a(l,j) = a(k,j)
3957                   a(k,j) = t
3958    20          continue
3959                call saxpy(n-k,t,a(k+1,k),1,a(k+1,j),1)
3960    30       continue
3961          go to 50
3962    40    continue
3963             info = k
3964    50    continue
3965    60 continue
3966    70 continue
3967       ipvt(n) = n
3968       if (a(n,n) .eq. 0.0e0) info = n
3969       return
3970       end subroutine sgefa                   
3973       integer function isamax(n,sx,incx)
3975 !     finds the index of element having max. absolute value.
3976 !     jack dongarra, linpack, 3/11/78.
3977 !     modified 3/93 to return if incx .le. 0.
3978 !     modified 12/3/93, array(1) declarations changed to array(*)
3980       real sx(*),smax
3981       integer i,incx,ix,n
3983       isamax = 0
3984       if( n.lt.1 .or. incx.le.0 ) return
3985       isamax = 1
3986       if(n.eq.1)return
3987       if(incx.eq.1)go to 20
3989 !        code for increment not equal to 1
3991       ix = 1
3992       smax = abs(sx(1))
3993       ix = ix + incx
3994       do 10 i = 2,n
3995          if(abs(sx(ix)).le.smax) go to 5
3996          isamax = i
3997          smax = abs(sx(ix))
3998     5    ix = ix + incx
3999    10 continue
4000       return
4002 !        code for increment equal to 1
4004    20 smax = abs(sx(1))
4005       do 30 i = 2,n
4006          if(abs(sx(i)).le.smax) go to 30
4007          isamax = i
4008          smax = abs(sx(i))
4009    30 continue
4010       return
4011       end function isamax           
4016 !* ======================================================================
4017 !* nist guide to available math software.
4018 !* fullsource for module sgbfa from package linpack.
4019 !* retrieved from netlib on wed jun 17 08:01:13 1998.
4020 !* ======================================================================
4021       subroutine sgbfa(abd,lda,n,ml,mu,ipvt,info)
4022       integer lda,n,ml,mu,ipvt(1),info
4023       real abd(lda,1)
4025 !     sgbfa factors a real band matrix by elimination.
4027 !     sgbfa is usually called by sgbco, but it can be called
4028 !     directly with a saving in time if  rcond  is not needed.
4030 !     on entry
4032 !        abd     real(lda, n)
4033 !                contains the matrix in band storage.  the columns
4034 !                of the matrix are stored in the columns of  abd  and
4035 !                the diagonals of the matrix are stored in rows
4036 !                ml+1 through 2*ml+mu+1 of  abd .
4037 !                see the comments below for details.
4039 !        lda     integer
4040 !                the leading dimension of the array  abd .
4041 !                lda must be .ge. 2*ml + mu + 1 .
4043 !        n       integer
4044 !                the order of the original matrix.
4046 !        ml      integer
4047 !                number of diagonals below the main diagonal.
4048 !                0 .le. ml .lt. n .
4050 !        mu      integer
4051 !                number of diagonals above the main diagonal.
4052 !                0 .le. mu .lt. n .
4053 !                more efficient if  ml .le. mu .
4054 !     on return
4056 !        abd     an upper triangular matrix in band storage and
4057 !                the multipliers which were used to obtain it.
4058 !                the factorization can be written  a = l*u  where
4059 !                l  is a product of permutation and unit lower
4060 !                triangular matrices and  u  is upper triangular.
4062 !        ipvt    integer(n)
4063 !                an integer vector of pivot indices.
4065 !        info    integer
4066 !                = 0  normal value.
4067 !                = k  if  u(k,k) .eq. 0.0 .  this is not an error
4068 !                     condition for this subroutine, but it does
4069 !                     indicate that sgbsl will divide by zero if
4070 !                     called.  use  rcond  in sgbco for a reliable
4071 !                     indication of singularity.
4073 !     band storage
4075 !           if  a  is a band matrix, the following program segment
4076 !           will set up the input.
4078 !                   ml = (band width below the diagonal)
4079 !                   mu = (band width above the diagonal)
4080 !                   m = ml + mu + 1
4081 !                   do 20 j = 1, n
4082 !                      i1 = max0(1, j-mu)
4083 !                      i2 = min0(n, j+ml)
4084 !                      do 10 i = i1, i2
4085 !                         k = i - j + m
4086 !                         abd(k,j) = a(i,j)
4087 !                10    continue
4088 !                20 continue
4090 !           this uses rows  ml+1  through  2*ml+mu+1  of  abd .
4091 !           in addition, the first  ml  rows in  abd  are used for
4092 !           elements generated during the triangularization.
4093 !           the total number of rows needed in  abd  is  2*ml+mu+1 .
4094 !           the  ml+mu by ml+mu  upper left triangle and the
4095 !           ml by ml  lower right triangle are not referenced.
4097 !     linpack. this version dated 08/14/78 .
4098 !     cleve moler, university of new mexico, argonne national lab.
4100 !     subroutines and functions
4102 !     blas saxpy,sscal,isamax
4103 !     fortran max0,min0
4105 !     internal variables
4107       real t
4108 !     integer i,isamax,i0,j,ju,jz,j0,j1,k,kp1,l,lm,m,mm,nm1      ! rce 2005-jan-21 - module conversion
4109       integer i,       i0,j,ju,jz,j0,j1,k,kp1,l,lm,m,mm,nm1
4112       m = ml + mu + 1
4113       info = 0
4115 !     zero initial fill-in columns
4117       j0 = mu + 2
4118       j1 = min0(n,m) - 1
4119       if (j1 .lt. j0) go to 30
4120       do 20 jz = j0, j1
4121          i0 = m + 1 - jz
4122          do 10 i = i0, ml
4123             abd(i,jz) = 0.0e0
4124    10    continue
4125    20 continue
4126    30 continue
4127       jz = j1
4128       ju = 0
4130 !     gaussian elimination with partial pivoting
4132       nm1 = n - 1
4133       if (nm1 .lt. 1) go to 130
4134       do 120 k = 1, nm1
4135          kp1 = k + 1
4137 !        zero next fill-in column
4139          jz = jz + 1
4140          if (jz .gt. n) go to 50
4141          if (ml .lt. 1) go to 50
4142             do 40 i = 1, ml
4143                abd(i,jz) = 0.0e0
4144    40       continue
4145    50    continue
4147 !        find l = pivot index
4149          lm = min0(ml,n-k)
4150          l = isamax(lm+1,abd(m,k),1) + m - 1
4151          ipvt(k) = l + k - m
4153 !        zero pivot implies this column already triangularized
4155          if (abd(l,k) .eq. 0.0e0) go to 100
4157 !           interchange if necessary
4159             if (l .eq. m) go to 60
4160                t = abd(l,k)
4161                abd(l,k) = abd(m,k)
4162                abd(m,k) = t
4163    60       continue
4165 !           compute multipliers
4167             t = -1.0e0/abd(m,k)
4168             call sscal(lm,t,abd(m+1,k),1)
4170 !           row elimination with column indexing
4172             ju = min0(max0(ju,mu+ipvt(k)),n)
4173             mm = m
4174             if (ju .lt. kp1) go to 90
4175             do 80 j = kp1, ju
4176                l = l - 1
4177                mm = mm - 1
4178                t = abd(l,j)
4179                if (l .eq. mm) go to 70
4180                   abd(l,j) = abd(mm,j)
4181                   abd(mm,j) = t
4182    70          continue
4183                call saxpy(lm,t,abd(m+1,k),1,abd(mm+1,j),1)
4184    80       continue
4185    90       continue
4186          go to 110
4187   100    continue
4188             info = k
4189   110    continue
4190   120 continue
4191   130 continue
4192       ipvt(n) = n
4193       if (abd(m,n) .eq. 0.0e0) info = n
4194       return
4195       end subroutine sgbfa                           
4201 !* ======================================================================
4202 !* nist guide to available math software.
4203 !* fullsource for module sgesl from package linpack.
4204 !* retrieved from netlib on wed jun 17 08:01:47 1998.
4205 !* ======================================================================
4206       subroutine sgesl(a,lda,n,ipvt,b,job)
4207       integer lda,n,ipvt(*),job
4208       real a(lda,*),b(*)
4210 !     sgesl solves the real system
4211 !     a * x = b  or  trans(a) * x = b
4212 !     using the factors computed by sgeco or sgefa.
4214 !     on entry
4216 !        a       real(lda, n)
4217 !                the output from sgeco or sgefa.
4219 !        lda     integer
4220 !                the leading dimension of the array  a .
4222 !        n       integer
4223 !                the order of the matrix  a .
4225 !        ipvt    integer(n)
4226 !                the pivot vector from sgeco or sgefa.
4228 !        b       real(n)
4229 !                the right hand side vector.
4231 !        job     integer
4232 !                = 0         to solve  a*x = b ,
4233 !                = nonzero   to solve  trans(a)*x = b  where
4234 !                            trans(a)  is the transpose.
4236 !     on return
4238 !        b       the solution vector  x .
4240 !     error condition
4242 !        a division by zero will occur if the input factor contains a
4243 !        zero on the diagonal.  technically this indicates singularity
4244 !        but it is often caused by improper arguments or improper
4245 !        setting of lda .  it will not occur if the subroutines are
4246 !        called correctly and if sgeco has set rcond .gt. 0.0
4247 !        or sgefa has set info .eq. 0 .
4249 !     to compute  inverse(a) * c  where  c  is a matrix
4250 !     with  p  columns
4251 !           call sgeco(a,lda,n,ipvt,rcond,z)
4252 !           if (rcond is too small) go to ...
4253 !           do 10 j = 1, p
4254 !              call sgesl(a,lda,n,ipvt,c(1,j),0)
4255 !        10 continue
4257 !     linpack. this version dated 08/14/78 .
4258 !     cleve moler, university of new mexico, argonne national lab.
4260 !     subroutines and functions
4262 !     blas saxpy,sdot
4264 !     internal variables
4266 !     real sdot,t      ! rce 2005-jan-21 - module conversion
4267       real      t
4268       integer k,kb,l,nm1
4270       nm1 = n - 1
4271       if (job .ne. 0) go to 50
4273 !        job = 0 , solve  a * x = b
4274 !        first solve  l*y = b
4276          if (nm1 .lt. 1) go to 30
4277          do 20 k = 1, nm1
4278             l = ipvt(k)
4279             t = b(l)
4280             if (l .eq. k) go to 10
4281                b(l) = b(k)
4282                b(k) = t
4283    10       continue
4284             call saxpy(n-k,t,a(k+1,k),1,b(k+1),1)
4285    20    continue
4286    30    continue
4288 !        now solve  u*x = y
4290          do 40 kb = 1, n
4291             k = n + 1 - kb
4292             b(k) = b(k)/a(k,k)
4293             t = -b(k)
4294             call saxpy(k-1,t,a(1,k),1,b(1),1)
4295    40    continue
4296       go to 100
4297    50 continue
4299 !        job = nonzero, solve  trans(a) * x = b
4300 !        first solve  trans(u)*y = b
4302          do 60 k = 1, n
4303             t = sdot(k-1,a(1,k),1,b(1),1)
4304             b(k) = (b(k) - t)/a(k,k)
4305    60    continue
4307 !        now solve trans(l)*x = y
4309          if (nm1 .lt. 1) go to 90
4310          do 80 kb = 1, nm1
4311             k = n - kb
4312             b(k) = b(k) + sdot(n-k,a(k+1,k),1,b(k+1),1)
4313             l = ipvt(k)
4314             if (l .eq. k) go to 70
4315                t = b(l)
4316                b(l) = b(k)
4317                b(k) = t
4318    70       continue
4319    80    continue
4320    90    continue
4321   100 continue
4322       return
4323       end subroutine sgesl                    
4325       real function sdot(n,sx,incx,sy,incy)
4327 !     forms the dot product of two vectors.
4328 !     uses unrolled loops for increments equal to one.
4329 !     jack dongarra, linpack, 3/11/78.
4330 !     modified 12/3/93, array(1) declarations changed to array(*)
4332       real sx(*),sy(*),stemp
4333       integer i,incx,incy,ix,iy,m,mp1,n
4335       stemp = 0.0e0
4336       sdot = 0.0e0
4337       if(n.le.0)return
4338       if(incx.eq.1.and.incy.eq.1)go to 20
4340 !        code for unequal increments or equal increments
4341 !          not equal to 1
4343       ix = 1
4344       iy = 1
4345       if(incx.lt.0)ix = (-n+1)*incx + 1
4346       if(incy.lt.0)iy = (-n+1)*incy + 1
4347       do 10 i = 1,n
4348         stemp = stemp + sx(ix)*sy(iy)
4349         ix = ix + incx
4350         iy = iy + incy
4351    10 continue
4352       sdot = stemp
4353       return
4355 !        code for both increments equal to 1
4358 !        clean-up loop
4360    20 m = mod(n,5)
4361       if( m .eq. 0 ) go to 40
4362       do 30 i = 1,m
4363         stemp = stemp + sx(i)*sy(i)
4364    30 continue
4365       if( n .lt. 5 ) go to 60
4366    40 mp1 = m + 1
4367       do 50 i = mp1,n,5
4368         stemp = stemp + sx(i)*sy(i) + sx(i + 1)*sy(i + 1) +   &
4369          sx(i + 2)*sy(i + 2) + sx(i + 3)*sy(i + 3) + sx(i + 4)*sy(i + 4)
4370    50 continue
4371    60 sdot = stemp
4372       return
4373       end function sdot                   
4375       subroutine sgbsl(abd,lda,n,ml,mu,ipvt,b,job)
4376       integer lda,n,ml,mu,ipvt(1),job
4377       real abd(lda,1),b(1)
4379 !     sgbsl solves the real band system
4380 !     a * x = b  or  trans(a) * x = b
4381 !     using the factors computed by sgbco or sgbfa.
4383 !     on entry
4385 !        abd     real(lda, n)
4386 !                the output from sgbco or sgbfa.
4388 !        lda     integer
4389 !                the leading dimension of the array  abd .
4391 !        n       integer
4392 !                the order of the original matrix.
4394 !        ml      integer
4395 !                number of diagonals below the main diagonal.
4397 !        mu      integer
4398 !                number of diagonals above the main diagonal.
4400 !        ipvt    integer(n)
4401 !                the pivot vector from sgbco or sgbfa.
4403 !        b       real(n)
4404 !                the right hand side vector.
4406 !        job     integer
4407 !                = 0         to solve  a*x = b ,
4408 !                = nonzero   to solve  trans(a)*x = b , where
4409 !                            trans(a)  is the transpose.
4411 !     on return
4413 !        b       the solution vector  x .
4415 !     error condition
4417 !        a division by zero will occur if the input factor contains a
4418 !        zero on the diagonal.  technically this indicates singularity
4419 !        but it is often caused by improper arguments or improper
4420 !        setting of lda .  it will not occur if the subroutines are
4421 !        called correctly and if sgbco has set rcond .gt. 0.0
4422 !        or sgbfa has set info .eq. 0 .
4424 !     to compute  inverse(a) * c  where  c  is a matrix
4425 !     with  p  columns
4426 !           call sgbco(abd,lda,n,ml,mu,ipvt,rcond,z)
4427 !           if (rcond is too small) go to ...
4428 !           do 10 j = 1, p
4429 !              call sgbsl(abd,lda,n,ml,mu,ipvt,c(1,j),0)
4430 !        10 continue
4432 !     linpack. this version dated 08/14/78 .
4433 !     cleve moler, university of new mexico, argonne national lab.
4435 !     subroutines and functions
4437 !     blas saxpy,sdot
4438 !     fortran min0
4440 !     internal variables
4442 !     real sdot,t      ! rce 2005-jan-21 - module conversion
4443       real      t
4444       integer k,kb,l,la,lb,lm,m,nm1
4446       m = mu + ml + 1
4447       nm1 = n - 1
4448       if (job .ne. 0) go to 50
4450 !        job = 0 , solve  a * x = b
4451 !        first solve l*y = b
4453          if (ml .eq. 0) go to 30
4454          if (nm1 .lt. 1) go to 30
4455             do 20 k = 1, nm1
4456                lm = min0(ml,n-k)
4457                l = ipvt(k)
4458                t = b(l)
4459                if (l .eq. k) go to 10
4460                   b(l) = b(k)
4461                   b(k) = t
4462    10          continue
4463                call saxpy(lm,t,abd(m+1,k),1,b(k+1),1)
4464    20       continue
4465    30    continue
4467 !        now solve  u*x = y
4469          do 40 kb = 1, n
4470             k = n + 1 - kb
4471             b(k) = b(k)/abd(m,k)
4472             lm = min0(k,m) - 1
4473             la = m - lm
4474             lb = k - lm
4475             t = -b(k)
4476             call saxpy(lm,t,abd(la,k),1,b(lb),1)
4477    40    continue
4478       go to 100
4479    50 continue
4481 !        job = nonzero, solve  trans(a) * x = b
4482 !        first solve  trans(u)*y = b
4484          do 60 k = 1, n
4485             lm = min0(k,m) - 1
4486             la = m - lm
4487             lb = k - lm
4488             t = sdot(lm,abd(la,k),1,b(lb),1)
4489             b(k) = (b(k) - t)/abd(m,k)
4490    60    continue
4492 !        now solve trans(l)*x = y
4494          if (ml .eq. 0) go to 90
4495          if (nm1 .lt. 1) go to 90
4496             do 80 kb = 1, nm1
4497                k = n - kb
4498                lm = min0(ml,n-k)
4499                b(k) = b(k) + sdot(lm,abd(m+1,k),1,b(k+1),1)
4500                l = ipvt(k)
4501                if (l .eq. k) go to 70
4502                   t = b(l)
4503                   b(l) = b(k)
4504                   b(k) = t
4505    70          continue
4506    80       continue
4507    90    continue
4508   100 continue
4509       return
4510       end subroutine sgbsl                            
4513       end module module_cmu_svode_solver