deabm.f 31 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671
  1. *DECK DEABM
  2. SUBROUTINE DEABM (F, NEQ, T, Y, TOUT, INFO, RTOL, ATOL, IDID,
  3. + RWORK, LRW, IWORK, LIW, RPAR, IPAR)
  4. C***BEGIN PROLOGUE DEABM
  5. C***PURPOSE Solve an initial value problem in ordinary differential
  6. C equations using an Adams-Bashforth method.
  7. C***LIBRARY SLATEC (DEPAC)
  8. C***CATEGORY I1A1B
  9. C***TYPE SINGLE PRECISION (DEABM-S, DDEABM-D)
  10. C***KEYWORDS ADAMS-BASHFORTH METHOD, DEPAC, INITIAL VALUE PROBLEMS,
  11. C ODE, ORDINARY DIFFERENTIAL EQUATIONS, PREDICTOR-CORRECTOR
  12. C***AUTHOR Shampine, L. F., (SNLA)
  13. C Watts, H. A., (SNLA)
  14. C***DESCRIPTION
  15. C
  16. C This is the Adams code in the package of differential equation
  17. C solvers DEPAC, consisting of the codes DERKF, DEABM, and DEBDF.
  18. C Design of the package was by L. F. Shampine and H. A. Watts.
  19. C It is documented in
  20. C SAND79-2374 , DEPAC - Design of a User Oriented Package of ODE
  21. C Solvers.
  22. C DEABM is a driver for a modification of the code ODE written by
  23. C L. F. Shampine and M. K. Gordon
  24. C Sandia Laboratories
  25. C Albuquerque, New Mexico 87185
  26. C
  27. C **********************************************************************
  28. C ** DEPAC PACKAGE OVERVIEW **
  29. C **************************************************
  30. C
  31. C You have a choice of three differential equation solvers from
  32. C DEPAC. The following brief descriptions are meant to aid you
  33. C in choosing the most appropriate code for your problem.
  34. C
  35. C DERKF is a fifth order Runge-Kutta code. It is the simplest of
  36. C the three choices, both algorithmically and in the use of the
  37. C code. DERKF is primarily designed to solve non-stiff and mild-
  38. C ly stiff differential equations when derivative evaluations are
  39. C not expensive. It should generally not be used to get high
  40. C accuracy results nor answers at a great many specific points.
  41. C Because DERKF has very low overhead costs, it will usually
  42. C result in the least expensive integration when solving
  43. C problems requiring a modest amount of accuracy and having
  44. C equations that are not costly to evaluate. DERKF attempts to
  45. C discover when it is not suitable for the task posed.
  46. C
  47. C DEABM is a variable order (one through twelve) Adams code.
  48. C Its complexity lies somewhere between that of DERKF and DEBDF.
  49. C DEABM is primarily designed to solve non-stiff and mildly stiff
  50. C differential equations when derivative evaluations are
  51. C expensive, high accuracy results are needed or answers at
  52. C many specific points are required. DEABM attempts to discover
  53. C when it is not suitable for the task posed.
  54. C
  55. C DEBDF is a variable order (one through five) backward
  56. C differentiation formula code. It is the most complicated of
  57. C the three choices. DEBDF is primarily designed to solve stiff
  58. C differential equations at crude to moderate tolerances.
  59. C If the problem is very stiff at all, DERKF and DEABM will be
  60. C quite inefficient compared to DEBDF. However, DEBDF will be
  61. C inefficient compared to DERKF and DEABM on non-stiff problems
  62. C because it uses much more storage, has a much larger overhead,
  63. C and the low order formulas will not give high accuracies
  64. C efficiently.
  65. C
  66. C The concept of stiffness cannot be described in a few words.
  67. C If you do not know the problem to be stiff, try either DERKF
  68. C or DEABM. Both of these codes will inform you of stiffness
  69. C when the cost of solving such problems becomes important.
  70. C
  71. C **********************************************************************
  72. C ** ABSTRACT **
  73. C **************
  74. C
  75. C Subroutine DEABM uses the Adams-Bashforth-Moulton predictor-
  76. C corrector formulas of orders one through twelve to integrate a
  77. C system of NEQ first order ordinary differential equations of the
  78. C form
  79. C DU/DX = F(X,U)
  80. C when the vector Y(*) of initial values for U(*) at X=T is given. The
  81. C subroutine integrates from T to TOUT. It is easy to continue the
  82. C integration to get results at additional TOUT. This is the interval
  83. C mode of operation. It is also easy for the routine to return with
  84. C the solution at each intermediate step on the way to TOUT. This is
  85. C the intermediate-output mode of operation.
  86. C
  87. C DEABM uses subprograms DES, STEPS, SINTRP, HSTART, HVNRM, R1MACH and
  88. C the error handling routine XERMSG. The only machine dependent
  89. C parameters to be assigned appear in R1MACH.
  90. C
  91. C **********************************************************************
  92. C ** DESCRIPTION OF THE ARGUMENTS TO DEABM (AN OVERVIEW) **
  93. C *********************************************************
  94. C
  95. C The parameters are
  96. C
  97. C F -- This is the name of a subroutine which you provide to
  98. C define the differential equations.
  99. C
  100. C NEQ -- This is the number of (first order) differential
  101. C equations to be integrated.
  102. C
  103. C T -- This is a value of the independent variable.
  104. C
  105. C Y(*) -- This array contains the solution components at T.
  106. C
  107. C TOUT -- This is a point at which a solution is desired.
  108. C
  109. C INFO(*) -- The basic task of the code is to integrate the
  110. C differential equations from T to TOUT and return an
  111. C answer at TOUT. INFO(*) is an integer array which is used
  112. C to communicate exactly how you want this task to be
  113. C carried out.
  114. C
  115. C RTOL, ATOL -- These quantities represent relative and absolute
  116. C error tolerances which you provide to indicate how
  117. C accurately you wish the solution to be computed. You may
  118. C choose them to be both scalars or else both vectors.
  119. C
  120. C IDID -- This scalar quantity is an indicator reporting what
  121. C the code did. You must monitor this integer variable to
  122. C decide what action to take next.
  123. C
  124. C RWORK(*), LRW -- RWORK(*) is a real work array of length LRW
  125. C which provides the code with needed storage space.
  126. C
  127. C IWORK(*), LIW -- IWORK(*) is an integer work array of length LIW
  128. C which provides the code with needed storage space and an
  129. C across call flag.
  130. C
  131. C RPAR, IPAR -- These are real and integer parameter arrays which
  132. C you can use for communication between your calling
  133. C program and the F subroutine.
  134. C
  135. C Quantities which are used as input items are
  136. C NEQ, T, Y(*), TOUT, INFO(*),
  137. C RTOL, ATOL, RWORK(1), LRW and LIW.
  138. C
  139. C Quantities which may be altered by the code are
  140. C T, Y(*), INFO(1), RTOL, ATOL,
  141. C IDID, RWORK(*) and IWORK(*).
  142. C
  143. C **********************************************************************
  144. C ** INPUT -- WHAT TO DO ON THE FIRST CALL TO DEABM **
  145. C ****************************************************
  146. C
  147. C The first call of the code is defined to be the start of each new
  148. C problem. Read through the descriptions of all the following items,
  149. C provide sufficient storage space for designated arrays, set
  150. C appropriate variables for the initialization of the problem, and
  151. C give information about how you want the problem to be solved.
  152. C
  153. C
  154. C F -- Provide a subroutine of the form
  155. C F(X,U,UPRIME,RPAR,IPAR)
  156. C to define the system of first order differential equations
  157. C which is to be solved. For the given values of X and the
  158. C vector U(*)=(U(1),U(2),...,U(NEQ)) , the subroutine must
  159. C evaluate the NEQ components of the system of differential
  160. C equations DU/DX = F(X,U) and store the derivatives in
  161. C array UPRIME(*), that is, UPRIME(I) = * DU(I)/DX * for
  162. C equations I=1,...,NEQ.
  163. C
  164. C Subroutine F must not alter X or U(*). You must declare
  165. C the name F in an external statement in your program that
  166. C calls DEABM. You must dimension U and UPRIME in F.
  167. C
  168. C RPAR and IPAR are real and integer parameter arrays which
  169. C you can use for communication between your calling program
  170. C and subroutine F. They are not used or altered by DEABM.
  171. C If you do not need RPAR or IPAR, ignore these parameters
  172. C by treating them as dummy arguments. If you do choose to
  173. C use them, dimension them in your calling program and in F
  174. C as arrays of appropriate length.
  175. C
  176. C NEQ -- Set it to the number of differential equations.
  177. C (NEQ .GE. 1)
  178. C
  179. C T -- Set it to the initial point of the integration.
  180. C You must use a program variable for T because the code
  181. C changes its value.
  182. C
  183. C Y(*) -- Set this vector to the initial values of the NEQ solution
  184. C components at the initial point. You must dimension Y at
  185. C least NEQ in your calling program.
  186. C
  187. C TOUT -- Set it to the first point at which a solution
  188. C is desired. You can take TOUT = T, in which case the code
  189. C will evaluate the derivative of the solution at T and
  190. C return. Integration either forward in T (TOUT .GT. T)
  191. C or backward in T (TOUT .LT. T) is permitted.
  192. C
  193. C The code advances the solution from T to TOUT using
  194. C step sizes which are automatically selected so as to
  195. C achieve the desired accuracy. If you wish, the code will
  196. C return with the solution and its derivative following
  197. C each intermediate step (intermediate-output mode) so that
  198. C you can monitor them, but you still must provide TOUT in
  199. C accord with the basic aim of the code.
  200. C
  201. C The first step taken by the code is a critical one
  202. C because it must reflect how fast the solution changes near
  203. C the initial point. The code automatically selects an
  204. C initial step size which is practically always suitable for
  205. C the problem. By using the fact that the code will not
  206. C step past TOUT in the first step, you could, if necessary,
  207. C restrict the length of the initial step size.
  208. C
  209. C For some problems it may not be permissible to integrate
  210. C past a point TSTOP because a discontinuity occurs there
  211. C or the solution or its derivative is not defined beyond
  212. C TSTOP. When you have declared a TSTOP point (see INFO(4)
  213. C and RWORK(1)), you have told the code not to integrate
  214. C past TSTOP. In this case any TOUT beyond TSTOP is invalid
  215. C input.
  216. C
  217. C INFO(*) -- Use the INFO array to give the code more details about
  218. C how you want your problem solved. This array should be
  219. C dimensioned of length 15 to accommodate other members of
  220. C DEPAC or possible future extensions, though DEABM uses
  221. C only the first four entries. You must respond to all of
  222. C the following items which are arranged as questions. The
  223. C simplest use of the code corresponds to answering all
  224. C questions as YES ,i.e. setting all entries of INFO to 0.
  225. C
  226. C INFO(1) -- This parameter enables the code to initialize
  227. C itself. You must set it to indicate the start of every
  228. C new problem.
  229. C
  230. C **** Is this the first call for this problem ...
  231. C YES -- Set INFO(1) = 0
  232. C NO -- Not applicable here.
  233. C See below for continuation calls. ****
  234. C
  235. C INFO(2) -- How much accuracy you want of your solution
  236. C is specified by the error tolerances RTOL and ATOL.
  237. C The simplest use is to take them both to be scalars.
  238. C To obtain more flexibility, they can both be vectors.
  239. C The code must be told your choice.
  240. C
  241. C **** Are both error tolerances RTOL, ATOL scalars ...
  242. C YES -- Set INFO(2) = 0
  243. C and input scalars for both RTOL and ATOL
  244. C NO -- Set INFO(2) = 1
  245. C and input arrays for both RTOL and ATOL ****
  246. C
  247. C INFO(3) -- The code integrates from T in the direction
  248. C of TOUT by steps. If you wish, it will return the
  249. C computed solution and derivative at the next
  250. C intermediate step (the intermediate-output mode) or
  251. C TOUT, whichever comes first. This is a good way to
  252. C proceed if you want to see the behavior of the solution.
  253. C If you must have solutions at a great many specific
  254. C TOUT points, this code will compute them efficiently.
  255. C
  256. C **** Do you want the solution only at
  257. C TOUT (and not at the next intermediate step) ...
  258. C YES -- Set INFO(3) = 0
  259. C NO -- Set INFO(3) = 1 ****
  260. C
  261. C INFO(4) -- To handle solutions at a great many specific
  262. C values TOUT efficiently, this code may integrate past
  263. C TOUT and interpolate to obtain the result at TOUT.
  264. C Sometimes it is not possible to integrate beyond some
  265. C point TSTOP because the equation changes there or it is
  266. C not defined past TSTOP. Then you must tell the code
  267. C not to go past.
  268. C
  269. C **** Can the integration be carried out without any
  270. C restrictions on the independent variable T ...
  271. C YES -- Set INFO(4)=0
  272. C NO -- Set INFO(4)=1
  273. C and define the stopping point TSTOP by
  274. C setting RWORK(1)=TSTOP ****
  275. C
  276. C RTOL, ATOL -- You must assign relative (RTOL) and absolute (ATOL)
  277. C error tolerances to tell the code how accurately you want
  278. C the solution to be computed. They must be defined as
  279. C program variables because the code may change them. You
  280. C have two choices --
  281. C both RTOL and ATOL are scalars. (INFO(2)=0)
  282. C both RTOL and ATOL are vectors. (INFO(2)=1)
  283. C In either case all components must be non-negative.
  284. C
  285. C The tolerances are used by the code in a local error test
  286. C at each step which requires roughly that
  287. C ABS(LOCAL ERROR) .LE. RTOL*ABS(Y)+ATOL
  288. C for each vector component.
  289. C (More specifically, a Euclidean norm is used to measure
  290. C the size of vectors, and the error test uses the magnitude
  291. C of the solution at the beginning of the step.)
  292. C
  293. C The true (global) error is the difference between the true
  294. C solution of the initial value problem and the computed
  295. C approximation. Practically all present day codes,
  296. C including this one, control the local error at each step
  297. C and do not even attempt to control the global error
  298. C directly. Roughly speaking, they produce a solution Y(T)
  299. C which satisfies the differential equations with a
  300. C residual R(T), DY(T)/DT = F(T,Y(T)) + R(T) ,
  301. C and, almost always, R(T) is bounded by the error
  302. C tolerances. Usually, but not always, the true accuracy of
  303. C the computed Y is comparable to the error tolerances. This
  304. C code will usually, but not always, deliver a more accurate
  305. C solution if you reduce the tolerances and integrate again.
  306. C By comparing two such solutions you can get a fairly
  307. C reliable idea of the true error in the solution at the
  308. C bigger tolerances.
  309. C
  310. C Setting ATOL=0.0 results in a pure relative error test on
  311. C that component. Setting RTOL=0.0 results in a pure abso-
  312. C lute error test on that component. A mixed test with non-
  313. C zero RTOL and ATOL corresponds roughly to a relative error
  314. C test when the solution component is much bigger than ATOL
  315. C and to an absolute error test when the solution component
  316. C is smaller than the threshold ATOL.
  317. C
  318. C Proper selection of the absolute error control parameters
  319. C ATOL requires you to have some idea of the scale of the
  320. C solution components. To acquire this information may mean
  321. C that you will have to solve the problem more than once.
  322. C In the absence of scale information, you should ask for
  323. C some relative accuracy in all the components (by setting
  324. C RTOL values non-zero) and perhaps impose extremely small
  325. C absolute error tolerances to protect against the danger of
  326. C a solution component becoming zero.
  327. C
  328. C The code will not attempt to compute a solution at an
  329. C accuracy unreasonable for the machine being used. It will
  330. C advise you if you ask for too much accuracy and inform
  331. C you as to the maximum accuracy it believes possible.
  332. C
  333. C RWORK(*) -- Dimension this real work array of length LRW in your
  334. C calling program.
  335. C
  336. C RWORK(1) -- If you have set INFO(4)=0, you can ignore this
  337. C optional input parameter. Otherwise you must define a
  338. C stopping point TSTOP by setting RWORK(1) = TSTOP.
  339. C (for some problems it may not be permissible to integrate
  340. C past a point TSTOP because a discontinuity occurs there
  341. C or the solution or its derivative is not defined beyond
  342. C TSTOP.)
  343. C
  344. C LRW -- Set it to the declared length of the RWORK array.
  345. C You must have LRW .GE. 130+21*NEQ
  346. C
  347. C IWORK(*) -- Dimension this integer work array of length LIW in
  348. C your calling program.
  349. C
  350. C LIW -- Set it to the declared length of the IWORK array.
  351. C You must have LIW .GE. 51
  352. C
  353. C RPAR, IPAR -- These are parameter arrays, of real and integer
  354. C type, respectively. You can use them for communication
  355. C between your program that calls DEABM and the F
  356. C subroutine. They are not used or altered by DEABM. If
  357. C you do not need RPAR or IPAR, ignore these parameters by
  358. C treating them as dummy arguments. If you do choose to use
  359. C them, dimension them in your calling program and in F as
  360. C arrays of appropriate length.
  361. C
  362. C **********************************************************************
  363. C ** OUTPUT -- AFTER ANY RETURN FROM DEABM **
  364. C *******************************************
  365. C
  366. C The principal aim of the code is to return a computed solution at
  367. C TOUT, although it is also possible to obtain intermediate results
  368. C along the way. To find out whether the code achieved its goal
  369. C or if the integration process was interrupted before the task was
  370. C completed, you must check the IDID parameter.
  371. C
  372. C
  373. C T -- The solution was successfully advanced to the
  374. C output value of T.
  375. C
  376. C Y(*) -- Contains the computed solution approximation at T.
  377. C You may also be interested in the approximate derivative
  378. C of the solution at T. It is contained in
  379. C RWORK(21),...,RWORK(20+NEQ).
  380. C
  381. C IDID -- Reports what the code did
  382. C
  383. C *** Task Completed ***
  384. C reported by positive values of IDID
  385. C
  386. C IDID = 1 -- A step was successfully taken in the
  387. C intermediate-output mode. The code has not
  388. C yet reached TOUT.
  389. C
  390. C IDID = 2 -- The integration to TOUT was successfully
  391. C completed (T=TOUT) by stepping exactly to TOUT.
  392. C
  393. C IDID = 3 -- The integration to TOUT was successfully
  394. C completed (T=TOUT) by stepping past TOUT.
  395. C Y(*) is obtained by interpolation.
  396. C
  397. C *** Task Interrupted ***
  398. C reported by negative values of IDID
  399. C
  400. C IDID = -1 -- A large amount of work has been expended.
  401. C (500 steps attempted)
  402. C
  403. C IDID = -2 -- The error tolerances are too stringent.
  404. C
  405. C IDID = -3 -- The local error test cannot be satisfied
  406. C because you specified a zero component in ATOL
  407. C and the corresponding computed solution
  408. C component is zero. Thus, a pure relative error
  409. C test is impossible for this component.
  410. C
  411. C IDID = -4 -- The problem appears to be stiff.
  412. C
  413. C IDID = -5,-6,-7,..,-32 -- Not applicable for this code
  414. C but used by other members of DEPAC or possible
  415. C future extensions.
  416. C
  417. C *** Task Terminated ***
  418. C reported by the value of IDID=-33
  419. C
  420. C IDID = -33 -- The code has encountered trouble from which
  421. C it cannot recover. A message is printed
  422. C explaining the trouble and control is returned
  423. C to the calling program. For example, this
  424. C occurs when invalid input is detected.
  425. C
  426. C RTOL, ATOL -- These quantities remain unchanged except when
  427. C IDID = -2. In this case, the error tolerances have been
  428. C increased by the code to values which are estimated to be
  429. C appropriate for continuing the integration. However, the
  430. C reported solution at T was obtained using the input values
  431. C of RTOL and ATOL.
  432. C
  433. C RWORK, IWORK -- Contain information which is usually of no
  434. C interest to the user but necessary for subsequent calls.
  435. C However, you may find use for
  436. C
  437. C RWORK(11)--Which contains the step size H to be
  438. C attempted on the next step.
  439. C
  440. C RWORK(12)--If the tolerances have been increased by the
  441. C code (IDID = -2) , they were multiplied by the
  442. C value in RWORK(12).
  443. C
  444. C RWORK(13)--Which contains the current value of the
  445. C independent variable, i.e. the farthest point
  446. C integration has reached. This will be dif-
  447. C ferent from T only when interpolation has been
  448. C performed (IDID=3).
  449. C
  450. C RWORK(20+I)--Which contains the approximate derivative of
  451. C the solution component Y(I). In DEABM, it is
  452. C obtained by calling subroutine F to evaluate
  453. C the differential equation using T and Y(*) when
  454. C IDID=1 or 2, and by interpolation when IDID=3.
  455. C
  456. C **********************************************************************
  457. C ** INPUT -- WHAT TO DO TO CONTINUE THE INTEGRATION **
  458. C ** (CALLS AFTER THE FIRST) **
  459. C *****************************************************
  460. C
  461. C This code is organized so that subsequent calls to continue the
  462. C integration involve little (if any) additional effort on your
  463. C part. You must monitor the IDID parameter in order to
  464. C determine what to do next.
  465. C
  466. C Recalling that the principal task of the code is to integrate
  467. C from T to TOUT (the interval mode), usually all you will need
  468. C to do is specify a new TOUT upon reaching the current TOUT.
  469. C
  470. C Do not alter any quantity not specifically permitted below,
  471. C in particular do not alter NEQ, T, Y(*), RWORK(*), IWORK(*) or
  472. C the differential equation in subroutine F. Any such alteration
  473. C constitutes a new problem and must be treated as such, i.e.
  474. C you must start afresh.
  475. C
  476. C You cannot change from vector to scalar error control or vice
  477. C versa (INFO(2)) but you can change the size of the entries of
  478. C RTOL, ATOL. Increasing a tolerance makes the equation easier
  479. C to integrate. Decreasing a tolerance will make the equation
  480. C harder to integrate and should generally be avoided.
  481. C
  482. C You can switch from the intermediate-output mode to the
  483. C interval mode (INFO(3)) or vice versa at any time.
  484. C
  485. C If it has been necessary to prevent the integration from going
  486. C past a point TSTOP (INFO(4), RWORK(1)), keep in mind that the
  487. C code will not integrate to any TOUT beyond the currently
  488. C specified TSTOP. Once TSTOP has been reached you must change
  489. C the value of TSTOP or set INFO(4)=0. You may change INFO(4)
  490. C or TSTOP at any time but you must supply the value of TSTOP in
  491. C RWORK(1) whenever you set INFO(4)=1.
  492. C
  493. C The parameter INFO(1) is used by the code to indicate the
  494. C beginning of a new problem and to indicate whether integration
  495. C is to be continued. You must input the value INFO(1) = 0
  496. C when starting a new problem. You must input the value
  497. C INFO(1) = 1 if you wish to continue after an interrupted task.
  498. C Do not set INFO(1) = 0 on a continuation call unless you
  499. C want the code to restart at the current T.
  500. C
  501. C *** Following a Completed Task ***
  502. C If
  503. C IDID = 1, call the code again to continue the integration
  504. C another step in the direction of TOUT.
  505. C
  506. C IDID = 2 or 3, define a new TOUT and call the code again.
  507. C TOUT must be different from T. You cannot change
  508. C the direction of integration without restarting.
  509. C
  510. C *** Following an Interrupted Task ***
  511. C To show the code that you realize the task was
  512. C interrupted and that you want to continue, you
  513. C must take appropriate action and reset INFO(1) = 1
  514. C If
  515. C IDID = -1, the code has attempted 500 steps.
  516. C If you want to continue, set INFO(1) = 1 and
  517. C call the code again. An additional 500 steps
  518. C will be allowed.
  519. C
  520. C IDID = -2, the error tolerances RTOL, ATOL have been
  521. C increased to values the code estimates appropriate
  522. C for continuing. You may want to change them
  523. C yourself. If you are sure you want to continue
  524. C with relaxed error tolerances, set INFO(1)=1 and
  525. C call the code again.
  526. C
  527. C IDID = -3, a solution component is zero and you set the
  528. C corresponding component of ATOL to zero. If you
  529. C are sure you want to continue, you must first
  530. C alter the error criterion to use positive values
  531. C for those components of ATOL corresponding to zero
  532. C solution components, then set INFO(1)=1 and call
  533. C the code again.
  534. C
  535. C IDID = -4, the problem appears to be stiff. It is very
  536. C inefficient to solve such problems with DEABM. The
  537. C code DEBDF in DEPAC handles this task efficiently.
  538. C If you are absolutely sure you want to continue
  539. C with DEABM, set INFO(1)=1 and call the code again.
  540. C
  541. C IDID = -5,-6,-7,..,-32 --- cannot occur with this code
  542. C but used by other members of DEPAC or possible
  543. C future extensions.
  544. C
  545. C *** Following a Terminated Task ***
  546. C If
  547. C IDID = -33, you cannot continue the solution of this
  548. C problem. An attempt to do so will result in your
  549. C run being terminated.
  550. C
  551. C **********************************************************************
  552. C
  553. C***REFERENCES L. F. Shampine and H. A. Watts, DEPAC - design of a user
  554. C oriented package of ODE solvers, Report SAND79-2374,
  555. C Sandia Laboratories, 1979.
  556. C***ROUTINES CALLED DES, XERMSG
  557. C***REVISION HISTORY (YYMMDD)
  558. C 800501 DATE WRITTEN
  559. C 890831 Modified array declarations. (WRB)
  560. C 891024 Changed references from VNORM to HVNRM. (WRB)
  561. C 891024 REVISION DATE from Version 3.2
  562. C 891214 Prologue converted to Version 4.0 format. (BAB)
  563. C 900510 Convert XERRWV calls to XERMSG calls. (RWC)
  564. C 920501 Reformatted the REFERENCES section. (WRB)
  565. C***END PROLOGUE DEABM
  566. C
  567. LOGICAL START,PHASE1,NORND,STIFF,INTOUT
  568. C
  569. DIMENSION Y(*),INFO(15),RTOL(*),ATOL(*),RWORK(*),IWORK(*),
  570. 1 RPAR(*),IPAR(*)
  571. C
  572. CHARACTER*8 XERN1
  573. CHARACTER*16 XERN3
  574. C
  575. EXTERNAL F
  576. C
  577. C CHECK FOR AN APPARENT INFINITE LOOP
  578. C
  579. C***FIRST EXECUTABLE STATEMENT DEABM
  580. IF ( INFO(1) .EQ. 0 ) IWORK(LIW) = 0
  581. IF (IWORK(LIW) .GE. 5) THEN
  582. IF (T .EQ. RWORK(21 + NEQ)) THEN
  583. WRITE (XERN3, '(1PE15.6)') T
  584. CALL XERMSG ('SLATEC', 'DEABM',
  585. * 'AN APPARENT INFINITE LOOP HAS BEEN DETECTED.$$' //
  586. * 'YOU HAVE MADE REPEATED CALLS AT T = ' // XERN3 //
  587. * ' AND THE INTEGRATION HAS NOT ADVANCED. CHECK THE ' //
  588. * 'WAY YOU HAVE SET PARAMETERS FOR THE CALL TO THE ' //
  589. * 'CODE, PARTICULARLY INFO(1).', 13, 2)
  590. RETURN
  591. ENDIF
  592. ENDIF
  593. C
  594. C CHECK LRW AND LIW FOR SUFFICIENT STORAGE ALLOCATION
  595. C
  596. IDID=0
  597. IF (LRW .LT. 130+21*NEQ) THEN
  598. WRITE (XERN1, '(I8)') LRW
  599. CALL XERMSG ('SLATEC', 'DEABM', 'THE LENGTH OF THE RWORK ' //
  600. * 'ARRAY MUST BE AT LEAST 130 + 21*NEQ.$$' //
  601. * 'YOU HAVE CALLED THE CODE WITH LRW = ' // XERN1, 1, 1)
  602. IDID=-33
  603. ENDIF
  604. C
  605. IF (LIW .LT. 51) THEN
  606. WRITE (XERN1, '(I8)') LIW
  607. CALL XERMSG ('SLATEC', 'DEABM', 'THE LENGTH OF THE IWORK ' //
  608. * 'ARRAY MUST BE AT LEAST 51.$$YOU HAVE CALLED THE CODE ' //
  609. * 'WITH LIW = ' // XERN1, 2, 1)
  610. IDID=-33
  611. ENDIF
  612. C
  613. C COMPUTE THE INDICES FOR THE ARRAYS TO BE STORED IN THE WORK ARRAY
  614. C
  615. IYPOUT = 21
  616. ITSTAR = NEQ + 21
  617. IYP = 1 + ITSTAR
  618. IYY = NEQ + IYP
  619. IWT = NEQ + IYY
  620. IP = NEQ + IWT
  621. IPHI = NEQ + IP
  622. IALPHA = (NEQ*16) + IPHI
  623. IBETA = 12 + IALPHA
  624. IPSI = 12 + IBETA
  625. IV = 12 + IPSI
  626. IW = 12 + IV
  627. ISIG = 12 + IW
  628. IG = 13 + ISIG
  629. IGI = 13 + IG
  630. IXOLD = 11 + IGI
  631. IHOLD = 1 + IXOLD
  632. ITOLD = 1 + IHOLD
  633. IDELSN = 1 + ITOLD
  634. ITWOU = 1 + IDELSN
  635. IFOURU = 1 + ITWOU
  636. C
  637. RWORK(ITSTAR) = T
  638. IF (INFO(1) .EQ. 0) GO TO 50
  639. START = IWORK(21) .NE. (-1)
  640. PHASE1 = IWORK(22) .NE. (-1)
  641. NORND = IWORK(23) .NE. (-1)
  642. STIFF = IWORK(24) .NE. (-1)
  643. INTOUT = IWORK(25) .NE. (-1)
  644. C
  645. 50 CALL DES(F,NEQ,T,Y,TOUT,INFO,RTOL,ATOL,IDID,RWORK(IYPOUT),
  646. 1 RWORK(IYP),RWORK(IYY),RWORK(IWT),RWORK(IP),RWORK(IPHI),
  647. 2 RWORK(IALPHA),RWORK(IBETA),RWORK(IPSI),RWORK(IV),
  648. 3 RWORK(IW),RWORK(ISIG),RWORK(IG),RWORK(IGI),RWORK(11),
  649. 4 RWORK(12),RWORK(13),RWORK(IXOLD),RWORK(IHOLD),
  650. 5 RWORK(ITOLD),RWORK(IDELSN),RWORK(1),RWORK(ITWOU),
  651. 5 RWORK(IFOURU),START,PHASE1,NORND,STIFF,INTOUT,IWORK(26),
  652. 6 IWORK(27),IWORK(28),IWORK(29),IWORK(30),IWORK(31),
  653. 7 IWORK(32),IWORK(33),IWORK(34),IWORK(35),IWORK(45),
  654. 8 RPAR,IPAR)
  655. C
  656. IWORK(21) = -1
  657. IF (START) IWORK(21) = 1
  658. IWORK(22) = -1
  659. IF (PHASE1) IWORK(22) = 1
  660. IWORK(23) = -1
  661. IF (NORND) IWORK(23) = 1
  662. IWORK(24) = -1
  663. IF (STIFF) IWORK(24) = 1
  664. IWORK(25) = -1
  665. IF (INTOUT) IWORK(25) = 1
  666. C
  667. IF (IDID .NE. (-2)) IWORK(LIW) = IWORK(LIW) + 1
  668. IF (T .NE. RWORK(ITSTAR)) IWORK(LIW) = 0
  669. C
  670. RETURN
  671. END