debdf.f 44 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925
  1. *DECK DEBDF
  2. SUBROUTINE DEBDF (F, NEQ, T, Y, TOUT, INFO, RTOL, ATOL, IDID,
  3. + RWORK, LRW, IWORK, LIW, RPAR, IPAR, JAC)
  4. C***BEGIN PROLOGUE DEBDF
  5. C***PURPOSE Solve an initial value problem in ordinary differential
  6. C equations using backward differentiation formulas. It is
  7. C intended primarily for stiff problems.
  8. C***LIBRARY SLATEC (DEPAC)
  9. C***CATEGORY I1A2
  10. C***TYPE SINGLE PRECISION (DEBDF-S, DDEBDF-D)
  11. C***KEYWORDS BACKWARD DIFFERENTIATION FORMULAS, DEPAC,
  12. C INITIAL VALUE PROBLEMS, ODE,
  13. C ORDINARY DIFFERENTIAL EQUATIONS, STIFF
  14. C***AUTHOR Shampine, L. F., (SNLA)
  15. C Watts, H. A., (SNLA)
  16. C***DESCRIPTION
  17. C
  18. C This is the backward differentiation code in the package of
  19. C differential equation solvers DEPAC, consisting of the codes
  20. C DERKF, DEABM, and DEBDF. Design of the package was by
  21. C L. F. Shampine and H. A. Watts. It is documented in
  22. C SAND-79-2374 , DEPAC - Design of a User Oriented Package of ODE
  23. C Solvers.
  24. C DEBDF is a driver for a modification of the code LSODE written by
  25. C A. C. Hindmarsh
  26. C Lawrence Livermore Laboratory
  27. C Livermore, California 94550
  28. C
  29. C **********************************************************************
  30. C ** DEPAC PACKAGE OVERVIEW **
  31. C **********************************************************************
  32. C
  33. C You have a choice of three differential equation solvers from
  34. C DEPAC. The following brief descriptions are meant to aid you
  35. C in choosing the most appropriate code for your problem.
  36. C
  37. C DERKF is a fifth order Runge-Kutta code. It is the simplest of
  38. C the three choices, both algorithmically and in the use of the
  39. C code. DERKF is primarily designed to solve non-stiff and mild-
  40. C ly stiff differential equations when derivative evaluations are
  41. C not expensive. It should generally not be used to get high
  42. C accuracy results nor answers at a great many specific points.
  43. C Because DERKF has very low overhead costs, it will usually
  44. C result in the least expensive integration when solving
  45. C problems requiring a modest amount of accuracy and having
  46. C equations that are not costly to evaluate. DERKF attempts to
  47. C discover when it is not suitable for the task posed.
  48. C
  49. C DEABM is a variable order (one through twelve) Adams code.
  50. C Its complexity lies somewhere between that of DERKF and DEBDF.
  51. C DEABM is primarily designed to solve non-stiff and mildly
  52. C stiff differential equations when derivative evaluations are
  53. C expensive, high accuracy results are needed or answers at
  54. C many specific points are required. DEABM attempts to discover
  55. C when it is not suitable for the task posed.
  56. C
  57. C DEBDF is a variable order (one through five) backward
  58. C differentiation formula code. It is the most complicated of
  59. C the three choices. DEBDF is primarily designed to solve stiff
  60. C differential equations at crude to moderate tolerances.
  61. C If the problem is very stiff at all, DERKF and DEABM will be
  62. C quite inefficient compared to DEBDF. However, DEBDF will be
  63. C inefficient compared to DERKF and DEABM on non-stiff problems
  64. C because it uses much more storage, has a much larger overhead,
  65. C and the low order formulas will not give high accuracies
  66. C efficiently.
  67. C
  68. C The concept of stiffness cannot be described in a few words.
  69. C If you do not know the problem to be stiff, try either DERKF
  70. C or DEABM. Both of these codes will inform you of stiffness
  71. C when the cost of solving such problems becomes important.
  72. C
  73. C **********************************************************************
  74. C ** ABSTRACT **
  75. C **********************************************************************
  76. C
  77. C Subroutine DEBDF uses the backward differentiation formulas of
  78. C orders one through five to integrate a system of NEQ first order
  79. C ordinary differential equations of the form
  80. C DU/DX = F(X,U)
  81. C when the vector Y(*) of initial values for U(*) at X=T is given. The
  82. C subroutine integrates from T to TOUT. It is easy to continue the
  83. C integration to get results at additional TOUT. This is the interval
  84. C mode of operation. It is also easy for the routine to return with
  85. C The solution at each intermediate step on the way to TOUT. This is
  86. C the intermediate-output mode of operation.
  87. C
  88. C **********************************************************************
  89. C ** DESCRIPTION OF THE ARGUMENTS TO DEBDF (AN OVERVIEW) **
  90. C **********************************************************************
  91. C
  92. C The Parameters are:
  93. C
  94. C F -- This is the name of a subroutine which you provide to
  95. C define the differential equations.
  96. C
  97. C NEQ -- This is the number of (first order) differential
  98. C equations to be integrated.
  99. C
  100. C T -- This is a value of the independent variable.
  101. C
  102. C Y(*) -- This array contains the solution components at T.
  103. C
  104. C TOUT -- This is a point at which a solution is desired.
  105. C
  106. C INFO(*) -- The basic task of the code is to integrate the
  107. C differential equations from T to TOUT and return an
  108. C answer at TOUT. INFO(*) is an INTEGER array which is used
  109. C to communicate exactly how you want this task to be
  110. C carried out.
  111. C
  112. C RTOL, ATOL -- These quantities
  113. C represent relative and absolute error tolerances which you
  114. C provide to indicate how accurately you wish the solution
  115. C to be computed. You may choose them to be both scalars
  116. C or else both vectors.
  117. C
  118. C IDID -- This scalar quantity is an indicator reporting what
  119. C the code did. You must monitor this INTEGER variable to
  120. C decide what action to take next.
  121. C
  122. C RWORK(*), LRW -- RWORK(*) is a REAL work array of
  123. C length LRW which provides the code with needed storage
  124. C space.
  125. C
  126. C IWORK(*), LIW -- IWORK(*) is an INTEGER work array of length LIW
  127. C which provides the code with needed storage space and an
  128. C across call flag.
  129. C
  130. C RPAR, IPAR -- These are REAL and INTEGER parameter
  131. C arrays which you can use for communication between your
  132. C calling program and the F subroutine (and the JAC
  133. C subroutine).
  134. C
  135. C JAC -- This is the name of a subroutine which you may choose to
  136. C provide for defining the Jacobian matrix of partial
  137. C derivatives DF/DU.
  138. C
  139. C Quantities which are used as input items are
  140. C NEQ, T, Y(*), TOUT, INFO(*),
  141. C RTOL, ATOL, RWORK(1), LRW,
  142. C IWORK(1), IWORK(2), and LIW.
  143. C
  144. C Quantities which may be altered by the code are
  145. C T, Y(*), INFO(1), RTOL, ATOL,
  146. C IDID, RWORK(*) and IWORK(*).
  147. C
  148. C **********************************************************************
  149. C * INPUT -- What To Do On The First Call To DEBDF *
  150. C **********************************************************************
  151. C
  152. C The first call of the code is defined to be the start of each new
  153. C problem. Read through the descriptions of all the following items,
  154. C provide sufficient storage space for designated arrays, set
  155. C appropriate variables for the initialization of the problem, and
  156. C give information about how you want the problem to be solved.
  157. C
  158. C
  159. C F -- provide a subroutine of the form
  160. C F(X,U,UPRIME,RPAR,IPAR)
  161. C to define the system of first order differential equations
  162. C which is to be solved. For the given values of X and the
  163. C vector U(*)=(U(1),U(2),...,U(NEQ)) , the subroutine must
  164. C evaluate the NEQ components of the system of differential
  165. C equations DU/DX=F(X,U) and store the derivatives in the
  166. C array UPRIME(*), that is, UPRIME(I) = * DU(I)/DX * for
  167. C equations I=1,...,NEQ.
  168. C
  169. C Subroutine F must not alter X or U(*). You must declare
  170. C the name F in an external statement in your program that
  171. C calls DEBDF. You must dimension U and UPRIME in F.
  172. C
  173. C RPAR and IPAR are REAL and INTEGER parameter arrays which
  174. C you can use for communication between your calling program
  175. C and subroutine F. They are not used or altered by DEBDF.
  176. C If you do not need RPAR or IPAR, ignore these parameters
  177. C by treating them as dummy arguments. If you do choose to
  178. C use them, dimension them in your calling program and in F
  179. C as arrays of appropriate length.
  180. C
  181. C NEQ -- Set it to the number of differential equations.
  182. C (NEQ .GE. 1)
  183. C
  184. C T -- Set it to the initial point of the integration.
  185. C You must use a program variable for T because the code
  186. C changes its value.
  187. C
  188. C Y(*) -- Set this vector to the initial values of the NEQ solution
  189. C components at the initial point. You must dimension Y at
  190. C least NEQ in your calling program.
  191. C
  192. C TOUT -- Set it to the first point at which a solution is desired.
  193. C You can take TOUT = T, in which case the code
  194. C will evaluate the derivative of the solution at T and
  195. C return. Integration either forward in T (TOUT .GT. T)
  196. C or backward in T (TOUT .LT. T) is permitted.
  197. C
  198. C The code advances the solution from T to TOUT using
  199. C step sizes which are automatically selected so as to
  200. C achieve the desired accuracy. If you wish, the code will
  201. C return with the solution and its derivative following
  202. C each intermediate step (intermediate-output mode) so that
  203. C you can monitor them, but you still must provide TOUT in
  204. C accord with the basic aim of the code.
  205. C
  206. C The first step taken by the code is a critical one
  207. C because it must reflect how fast the solution changes near
  208. C the initial point. The code automatically selects an
  209. C initial step size which is practically always suitable for
  210. C the problem. By using the fact that the code will not
  211. C step past TOUT in the first step, you could, if necessary,
  212. C restrict the length of the initial step size.
  213. C
  214. C For some problems it may not be permissible to integrate
  215. C past a point TSTOP because a discontinuity occurs there
  216. C or the solution or its derivative is not defined beyond
  217. C TSTOP. When you have declared a TSTOP point (see INFO(4)
  218. C and RWORK(1)), you have told the code not to integrate
  219. C past TSTOP. In this case any TOUT beyond TSTOP is invalid
  220. C input.
  221. C
  222. C INFO(*) -- Use the INFO array to give the code more details about
  223. C how you want your problem solved. This array should be
  224. C dimensioned of length 15 to accommodate other members of
  225. C DEPAC or possible future extensions, though DEBDF uses
  226. C only the first six entries. You must respond to all of
  227. C the following items which are arranged as questions. The
  228. C simplest use of the code corresponds to answering all
  229. C questions as YES ,i.e. setting all entries of INFO to 0.
  230. C
  231. C INFO(1) -- This parameter enables the code to initialize
  232. C itself. You must set it to indicate the start of every
  233. C new problem.
  234. C
  235. C **** Is this the first call for this problem ...
  236. C YES -- Set INFO(1) = 0
  237. C NO -- Not applicable here.
  238. C See below for continuation calls. ****
  239. C
  240. C INFO(2) -- How much accuracy you want of your solution
  241. C is specified by the error tolerances RTOL and ATOL.
  242. C The simplest use is to take them both to be scalars.
  243. C To obtain more flexibility, they can both be vectors.
  244. C The code must be told your choice.
  245. C
  246. C **** Are both error tolerances RTOL, ATOL scalars ...
  247. C YES -- Set INFO(2) = 0
  248. C and input scalars for both RTOL and ATOL
  249. C NO -- Set INFO(2) = 1
  250. C and input arrays for both RTOL and ATOL ****
  251. C
  252. C INFO(3) -- The code integrates from T in the direction
  253. C of TOUT by steps. If you wish, it will return the
  254. C computed solution and derivative at the next
  255. C intermediate step (the intermediate-output mode) or
  256. C TOUT, whichever comes first. This is a good way to
  257. C proceed if you want to see the behavior of the solution.
  258. C If you must have solutions at a great many specific
  259. C TOUT points, this code will compute them efficiently.
  260. C
  261. C **** Do you want the solution only at
  262. C TOUT (and NOT at the next intermediate step) ...
  263. C YES -- Set INFO(3) = 0
  264. C NO -- Set INFO(3) = 1 ****
  265. C
  266. C INFO(4) -- To handle solutions at a great many specific
  267. C values TOUT efficiently, this code may integrate past
  268. C TOUT and interpolate to obtain the result at TOUT.
  269. C Sometimes it is not possible to integrate beyond some
  270. C point TSTOP because the equation changes there or it is
  271. C not defined past TSTOP. Then you must tell the code
  272. C not to go past.
  273. C
  274. C **** Can the integration be carried out without any
  275. C restrictions on the independent variable T ...
  276. C YES -- Set INFO(4)=0
  277. C NO -- Set INFO(4)=1
  278. C and define the stopping point TSTOP by
  279. C setting RWORK(1)=TSTOP ****
  280. C
  281. C INFO(5) -- To solve stiff problems it is necessary to use the
  282. C Jacobian matrix of partial derivatives of the system
  283. C of differential equations. If you do not provide a
  284. C subroutine to evaluate it analytically (see the
  285. C description of the item JAC in the call list), it will
  286. C be approximated by numerical differencing in this code.
  287. C Although it is less trouble for you to have the code
  288. C compute partial derivatives by numerical differencing,
  289. C the solution will be more reliable if you provide the
  290. C derivatives via JAC. Sometimes numerical differencing
  291. C is cheaper than evaluating derivatives in JAC and
  292. C sometimes it is not - this depends on your problem.
  293. C
  294. C If your problem is linear, i.e. has the form
  295. C DU/DX = F(X,U) = J(X)*U + G(X) for some matrix J(X)
  296. C and vector G(X), the Jacobian matrix DF/DU = J(X).
  297. C Since you must provide a subroutine to evaluate F(X,U)
  298. C analytically, it is little extra trouble to provide
  299. C subroutine JAC for evaluating J(X) analytically.
  300. C Furthermore, in such cases, numerical differencing is
  301. C much more expensive than analytic evaluation.
  302. C
  303. C **** Do you want the code to evaluate the partial
  304. C derivatives automatically by numerical differences ...
  305. C YES -- Set INFO(5)=0
  306. C NO -- Set INFO(5)=1
  307. C and provide subroutine JAC for evaluating the
  308. C Jacobian matrix ****
  309. C
  310. C INFO(6) -- DEBDF will perform much better if the Jacobian
  311. C matrix is banded and the code is told this. In this
  312. C case, the storage needed will be greatly reduced,
  313. C numerical differencing will be performed more cheaply,
  314. C and a number of important algorithms will execute much
  315. C faster. The differential equation is said to have
  316. C half-bandwidths ML (lower) and MU (upper) if equation I
  317. C involves only unknowns Y(J) with
  318. C I-ML .LE. J .LE. I+MU
  319. C for all I=1,2,...,NEQ. Thus, ML and MU are the widths
  320. C of the lower and upper parts of the band, respectively,
  321. C with the main diagonal being excluded. If you do not
  322. C indicate that the equation has a banded Jacobian,
  323. C the code works with a full matrix of NEQ**2 elements
  324. C (stored in the conventional way). Computations with
  325. C banded matrices cost less time and storage than with
  326. C full matrices if 2*ML+MU .LT. NEQ. If you tell the
  327. C code that the Jacobian matrix has a banded structure and
  328. C you want to provide subroutine JAC to compute the
  329. C partial derivatives, then you must be careful to store
  330. C the elements of the Jacobian matrix in the special form
  331. C indicated in the description of JAC.
  332. C
  333. C **** Do you want to solve the problem using a full
  334. C (dense) Jacobian matrix (and not a special banded
  335. C structure) ...
  336. C YES -- Set INFO(6)=0
  337. C NO -- Set INFO(6)=1
  338. C and provide the lower (ML) and upper (MU)
  339. C bandwidths by setting
  340. C IWORK(1)=ML
  341. C IWORK(2)=MU ****
  342. C
  343. C RTOL, ATOL -- You must assign relative (RTOL) and absolute (ATOL)
  344. C error tolerances to tell the code how accurately you want
  345. C the solution to be computed. They must be defined as
  346. C program variables because the code may change them. You
  347. C have two choices --
  348. C Both RTOL and ATOL are scalars. (INFO(2)=0)
  349. C Both RTOL and ATOL are vectors. (INFO(2)=1)
  350. C In either case all components must be non-negative.
  351. C
  352. C The tolerances are used by the code in a local error test
  353. C at each step which requires roughly that
  354. C ABS(LOCAL ERROR) .LE. RTOL*ABS(Y)+ATOL
  355. C for each vector component.
  356. C (More specifically, a root-mean-square norm is used to
  357. C measure the size of vectors, and the error test uses the
  358. C magnitude of the solution at the beginning of the step.)
  359. C
  360. C The true (global) error is the difference between the true
  361. C solution of the initial value problem and the computed
  362. C approximation. Practically all present day codes,
  363. C including this one, control the local error at each step
  364. C and do not even attempt to control the global error
  365. C directly. Roughly speaking, they produce a solution Y(T)
  366. C which satisfies the differential equations with a
  367. C residual R(T), DY(T)/DT = F(T,Y(T)) + R(T) ,
  368. C and, almost always, R(T) is bounded by the error
  369. C tolerances. Usually, but not always, the true accuracy of
  370. C the computed Y is comparable to the error tolerances. This
  371. C code will usually, but not always, deliver a more accurate
  372. C solution if you reduce the tolerances and integrate again.
  373. C By comparing two such solutions you can get a fairly
  374. C reliable idea of the true error in the solution at the
  375. C bigger tolerances.
  376. C
  377. C Setting ATOL=0. results in a pure relative error test on
  378. C that component. Setting RTOL=0. results in a pure abso-
  379. C lute error test on that component. A mixed test with non-
  380. C zero RTOL and ATOL corresponds roughly to a relative error
  381. C test when the solution component is much bigger than ATOL
  382. C and to an absolute error test when the solution component
  383. C is smaller than the threshold ATOL.
  384. C
  385. C Proper selection of the absolute error control parameters
  386. C ATOL requires you to have some idea of the scale of the
  387. C solution components. To acquire this information may mean
  388. C that you will have to solve the problem more than once. In
  389. C the absence of scale information, you should ask for some
  390. C relative accuracy in all the components (by setting RTOL
  391. C values non-zero) and perhaps impose extremely small
  392. C absolute error tolerances to protect against the danger of
  393. C a solution component becoming zero.
  394. C
  395. C The code will not attempt to compute a solution at an
  396. C accuracy unreasonable for the machine being used. It will
  397. C advise you if you ask for too much accuracy and inform
  398. C you as to the maximum accuracy it believes possible.
  399. C
  400. C RWORK(*) -- Dimension this REAL work array of length LRW in your
  401. C calling program.
  402. C
  403. C RWORK(1) -- If you have set INFO(4)=0, you can ignore this
  404. C optional input parameter. Otherwise you must define a
  405. C stopping point TSTOP by setting RWORK(1) = TSTOP.
  406. C (For some problems it may not be permissible to integrate
  407. C past a point TSTOP because a discontinuity occurs there
  408. C or the solution or its derivative is not defined beyond
  409. C TSTOP.)
  410. C
  411. C LRW -- Set it to the declared length of the RWORK array.
  412. C You must have
  413. C LRW .GE. 250+10*NEQ+NEQ**2
  414. C for the full (dense) Jacobian case (when INFO(6)=0), or
  415. C LRW .GE. 250+10*NEQ+(2*ML+MU+1)*NEQ
  416. C for the banded Jacobian case (when INFO(6)=1).
  417. C
  418. C IWORK(*) -- Dimension this INTEGER work array of length LIW in
  419. C your calling program.
  420. C
  421. C IWORK(1), IWORK(2) -- If you have set INFO(6)=0, you can ignore
  422. C these optional input parameters. Otherwise you must define
  423. C the half-bandwidths ML (lower) and MU (upper) of the
  424. C Jacobian matrix by setting IWORK(1) = ML and
  425. C IWORK(2) = MU. (The code will work with a full matrix
  426. C of NEQ**2 elements unless it is told that the problem has
  427. C a banded Jacobian, in which case the code will work with
  428. C a matrix containing at most (2*ML+MU+1)*NEQ elements.)
  429. C
  430. C LIW -- Set it to the declared length of the IWORK array.
  431. C You must have LIW .GE. 56+NEQ.
  432. C
  433. C RPAR, IPAR -- These are parameter arrays, of REAL and INTEGER
  434. C type, respectively. You can use them for communication
  435. C between your program that calls DEBDF and the F
  436. C subroutine (and the JAC subroutine). They are not used or
  437. C altered by DEBDF. If you do not need RPAR or IPAR, ignore
  438. C these parameters by treating them as dummy arguments. If
  439. C you do choose to use them, dimension them in your calling
  440. C program and in F (and in JAC) as arrays of appropriate
  441. C length.
  442. C
  443. C JAC -- If you have set INFO(5)=0, you can ignore this parameter
  444. C by treating it as a dummy argument. (For some compilers
  445. C you may have to write a dummy subroutine named JAC in
  446. C order to avoid problems associated with missing external
  447. C routine names.) Otherwise, you must provide a subroutine
  448. C of the form
  449. C JAC(X,U,PD,NROWPD,RPAR,IPAR)
  450. C to define the Jacobian matrix of partial derivatives DF/DU
  451. C of the system of differential equations DU/DX = F(X,U).
  452. C For the given values of X and the vector
  453. C U(*)=(U(1),U(2),...,U(NEQ)), the subroutine must evaluate
  454. C the non-zero partial derivatives DF(I)/DU(J) for each
  455. C differential equation I=1,...,NEQ and each solution
  456. C component J=1,...,NEQ , and store these values in the
  457. C matrix PD. The elements of PD are set to zero before each
  458. C call to JAC so only non-zero elements need to be defined.
  459. C
  460. C Subroutine JAC must not alter X, U(*), or NROWPD. You
  461. C must declare the name JAC in an EXTERNAL statement in your
  462. C program that calls DEBDF. NROWPD is the row dimension of
  463. C the PD matrix and is assigned by the code. Therefore you
  464. C must dimension PD in JAC according to
  465. C DIMENSION PD(NROWPD,1)
  466. C You must also dimension U in JAC.
  467. C
  468. C The way you must store the elements into the PD matrix
  469. C depends on the structure of the Jacobian which you
  470. C indicated by INFO(6).
  471. C *** INFO(6)=0 -- Full (Dense) Jacobian ***
  472. C When you evaluate the (non-zero) partial derivative
  473. C of equation I with respect to variable J, you must
  474. C store it in PD according to
  475. C PD(I,J) = * DF(I)/DU(J) *
  476. C *** INFO(6)=1 -- Banded Jacobian with ML Lower and MU
  477. C Upper Diagonal Bands (refer to INFO(6) description of
  478. C ML and MU) ***
  479. C When you evaluate the (non-zero) partial derivative
  480. C of equation I with respect to variable J, you must
  481. C store it in PD according to
  482. C IROW = I - J + ML + MU + 1
  483. C PD(IROW,J) = * DF(I)/DU(J) *
  484. C
  485. C RPAR and IPAR are REAL and INTEGER parameter
  486. C arrays which you can use for communication between your
  487. C calling program and your Jacobian subroutine JAC. They
  488. C are not altered by DEBDF. If you do not need RPAR or
  489. C IPAR, ignore these parameters by treating them as dummy
  490. C arguments. If you do choose to use them, dimension them
  491. C in your calling program and in JAC as arrays of
  492. C appropriate length.
  493. C
  494. C **********************************************************************
  495. C * OUTPUT -- After any return from DDEBDF *
  496. C **********************************************************************
  497. C
  498. C The principal aim of the code is to return a computed solution at
  499. C TOUT, although it is also possible to obtain intermediate results
  500. C along the way. To find out whether the code achieved its goal
  501. C or if the integration process was interrupted before the task was
  502. C completed, you must check the IDID parameter.
  503. C
  504. C
  505. C T -- The solution was successfully advanced to the
  506. C output value of T.
  507. C
  508. C Y(*) -- Contains the computed solution approximation at T.
  509. C You may also be interested in the approximate derivative
  510. C of the solution at T. It is contained in
  511. C RWORK(21),...,RWORK(20+NEQ).
  512. C
  513. C IDID -- Reports what the code did
  514. C
  515. C *** Task Completed ***
  516. C Reported by positive values of IDID
  517. C
  518. C IDID = 1 -- A step was successfully taken in the
  519. C intermediate-output mode. The code has not
  520. C yet reached TOUT.
  521. C
  522. C IDID = 2 -- The integration to TOUT was successfully
  523. C completed (T=TOUT) by stepping exactly to TOUT.
  524. C
  525. C IDID = 3 -- The integration to TOUT was successfully
  526. C completed (T=TOUT) by stepping past TOUT.
  527. C Y(*) is obtained by interpolation.
  528. C
  529. C *** Task Interrupted ***
  530. C Reported by negative values of IDID
  531. C
  532. C IDID = -1 -- A large amount of work has been expended.
  533. C (500 steps attempted)
  534. C
  535. C IDID = -2 -- The error tolerances are too stringent.
  536. C
  537. C IDID = -3 -- The local error test cannot be satisfied
  538. C because you specified a zero component in ATOL
  539. C and the corresponding computed solution
  540. C component is zero. Thus, a pure relative error
  541. C test is impossible for this component.
  542. C
  543. C IDID = -4,-5 -- Not applicable for this code but used
  544. C by other members of DEPAC.
  545. C
  546. C IDID = -6 -- DEBDF had repeated convergence test failures
  547. C on the last attempted step.
  548. C
  549. C IDID = -7 -- DEBDF had repeated error test failures on
  550. C the last attempted step.
  551. C
  552. C IDID = -8,..,-32 -- Not applicable for this code but
  553. C used by other members of DEPAC or possible
  554. C future extensions.
  555. C
  556. C *** Task Terminated ***
  557. C Reported by the value of IDID=-33
  558. C
  559. C IDID = -33 -- The code has encountered trouble from which
  560. C it cannot recover. A message is printed
  561. C explaining the trouble and control is returned
  562. C to the calling program. For example, this
  563. C occurs when invalid input is detected.
  564. C
  565. C RTOL, ATOL -- These quantities remain unchanged except when
  566. C IDID = -2. In this case, the error tolerances have been
  567. C increased by the code to values which are estimated to be
  568. C appropriate for continuing the integration. However, the
  569. C reported solution at T was obtained using the input values
  570. C of RTOL and ATOL.
  571. C
  572. C RWORK, IWORK -- Contain information which is usually of no
  573. C interest to the user but necessary for subsequent calls.
  574. C However, you may find use for
  575. C
  576. C RWORK(11)--which contains the step size H to be
  577. C attempted on the next step.
  578. C
  579. C RWORK(12)--If the tolerances have been increased by the
  580. C code (IDID = -2) , they were multiplied by the
  581. C value in RWORK(12).
  582. C
  583. C RWORK(13)--which contains the current value of the
  584. C independent variable, i.e. the farthest point
  585. C integration has reached. This will be
  586. C different from T only when interpolation has
  587. C been performed (IDID=3).
  588. C
  589. C RWORK(20+I)--which contains the approximate derivative
  590. C of the solution component Y(I). In DEBDF, it
  591. C is never obtained by calling subroutine F to
  592. C evaluate the differential equation using T and
  593. C Y(*), except at the initial point of
  594. C integration.
  595. C
  596. C **********************************************************************
  597. C ** INPUT -- What To Do To Continue The Integration **
  598. C ** (calls after the first) **
  599. C **********************************************************************
  600. C
  601. C This code is organized so that subsequent calls to continue the
  602. C integration involve little (if any) additional effort on your
  603. C part. You must monitor the IDID parameter in order to determine
  604. C what to do next.
  605. C
  606. C Recalling that the principal task of the code is to integrate
  607. C from T to TOUT (the interval mode), usually all you will need
  608. C to do is specify a new TOUT upon reaching the current TOUT.
  609. C
  610. C Do not alter any quantity not specifically permitted below,
  611. C in particular do not alter NEQ, T, Y(*), RWORK(*), IWORK(*) or
  612. C the differential equation in subroutine F. Any such alteration
  613. C constitutes a new problem and must be treated as such, i.e.
  614. C you must start afresh.
  615. C
  616. C You cannot change from vector to scalar error control or vice
  617. C versa (INFO(2)) but you can change the size of the entries of
  618. C RTOL, ATOL. Increasing a tolerance makes the equation easier
  619. C to integrate. Decreasing a tolerance will make the equation
  620. C harder to integrate and should generally be avoided.
  621. C
  622. C You can switch from the intermediate-output mode to the
  623. C interval mode (INFO(3)) or vice versa at any time.
  624. C
  625. C If it has been necessary to prevent the integration from going
  626. C past a point TSTOP (INFO(4), RWORK(1)), keep in mind that the
  627. C code will not integrate to any TOUT beyond the currently
  628. C specified TSTOP. Once TSTOP has been reached you must change
  629. C the value of TSTOP or set INFO(4)=0. You may change INFO(4)
  630. C or TSTOP at any time but you must supply the value of TSTOP in
  631. C RWORK(1) whenever you set INFO(4)=1.
  632. C
  633. C Do not change INFO(5), INFO(6), IWORK(1), or IWORK(2)
  634. C unless you are going to restart the code.
  635. C
  636. C The parameter INFO(1) is used by the code to indicate the
  637. C beginning of a new problem and to indicate whether integration
  638. C is to be continued. You must input the value INFO(1) = 0
  639. C when starting a new problem. You must input the value
  640. C INFO(1) = 1 if you wish to continue after an interrupted task.
  641. C Do not set INFO(1) = 0 on a continuation call unless you
  642. C want the code to restart at the current T.
  643. C
  644. C *** Following a Completed Task ***
  645. C If
  646. C IDID = 1, call the code again to continue the integration
  647. C another step in the direction of TOUT.
  648. C
  649. C IDID = 2 or 3, define a new TOUT and call the code again.
  650. C TOUT must be different from T. You cannot change
  651. C the direction of integration without restarting.
  652. C
  653. C *** Following an Interrupted Task ***
  654. C To show the code that you realize the task was
  655. C interrupted and that you want to continue, you
  656. C must take appropriate action and reset INFO(1) = 1
  657. C If
  658. C IDID = -1, the code has attempted 500 steps.
  659. C If you want to continue, set INFO(1) = 1 and
  660. C call the code again. An additional 500 steps
  661. C will be allowed.
  662. C
  663. C IDID = -2, the error tolerances RTOL, ATOL have been
  664. C increased to values the code estimates appropriate
  665. C for continuing. You may want to change them
  666. C yourself. If you are sure you want to continue
  667. C with relaxed error tolerances, set INFO(1)=1 and
  668. C call the code again.
  669. C
  670. C IDID = -3, a solution component is zero and you set the
  671. C corresponding component of ATOL to zero. If you
  672. C are sure you want to continue, you must first
  673. C alter the error criterion to use positive values
  674. C for those components of ATOL corresponding to zero
  675. C solution components, then set INFO(1)=1 and call
  676. C the code again.
  677. C
  678. C IDID = -4,-5 --- cannot occur with this code but used
  679. C by other members of DEPAC.
  680. C
  681. C IDID = -6, repeated convergence test failures occurred
  682. C on the last attempted step in DEBDF. An inaccu-
  683. C rate Jacobian may be the problem. If you are
  684. C absolutely certain you want to continue, restart
  685. C the integration at the current T by setting
  686. C INFO(1)=0 and call the code again.
  687. C
  688. C IDID = -7, repeated error test failures occurred on the
  689. C last attempted step in DEBDF. A singularity in
  690. C the solution may be present. You should re-
  691. C examine the problem being solved. If you are
  692. C absolutely certain you want to continue, restart
  693. C the integration at the current T by setting
  694. C INFO(1)=0 and call the code again.
  695. C
  696. C IDID = -8,..,-32 --- cannot occur with this code but
  697. C used by other members of DEPAC or possible future
  698. C extensions.
  699. C
  700. C *** Following a Terminated Task ***
  701. C If
  702. C IDID = -33, you cannot continue the solution of this
  703. C problem. An attempt to do so will result in your
  704. C run being terminated.
  705. C
  706. C **********************************************************************
  707. C
  708. C ***** Warning *****
  709. C
  710. C If DEBDF is to be used in an overlay situation, you must save and
  711. C restore certain items used internally by DEBDF (values in the
  712. C common block DEBDF1). This can be accomplished as follows.
  713. C
  714. C To save the necessary values upon return from DEBDF, simply call
  715. C SVCO(RWORK(22+NEQ),IWORK(21+NEQ)).
  716. C
  717. C To restore the necessary values before the next call to DEBDF,
  718. C simply call RSCO(RWORK(22+NEQ),IWORK(21+NEQ)).
  719. C
  720. C***REFERENCES L. F. Shampine and H. A. Watts, DEPAC - design of a user
  721. C oriented package of ODE solvers, Report SAND79-2374,
  722. C Sandia Laboratories, 1979.
  723. C***ROUTINES CALLED LSOD, XERMSG
  724. C***COMMON BLOCKS DEBDF1
  725. C***REVISION HISTORY (YYMMDD)
  726. C 800901 DATE WRITTEN
  727. C 890831 Modified array declarations. (WRB)
  728. C 891024 Changed references from VNORM to HVNRM. (WRB)
  729. C 891024 REVISION DATE from Version 3.2
  730. C 891214 Prologue converted to Version 4.0 format. (BAB)
  731. C 900326 Removed duplicate information from DESCRIPTION section.
  732. C (WRB)
  733. C 900510 Convert XERRWV calls to XERMSG calls, change Prologue
  734. C comments to agree with DDEBDF. (RWC)
  735. C 920501 Reformatted the REFERENCES section. (WRB)
  736. C***END PROLOGUE DEBDF
  737. C
  738. C
  739. LOGICAL INTOUT
  740. CHARACTER*8 XERN1, XERN2
  741. CHARACTER*16 XERN3
  742. C
  743. DIMENSION Y(*),INFO(15),RTOL(*),ATOL(*),RWORK(*),IWORK(*),
  744. 1 RPAR(*),IPAR(*)
  745. C
  746. COMMON /DEBDF1/ TOLD, ROWNS(210),
  747. 1 EL0, H, HMIN, HMXI, HU, TN, UROUND,
  748. 2 IQUIT, INIT, IYH, IEWT, IACOR, ISAVF, IWM, KSTEPS,
  749. 3 IBEGIN, ITOL, IINTEG, ITSTOP, IJAC, IBAND, IOWNS(6),
  750. 4 IER, JSTART, KFLAG, L, METH, MITER, MAXORD, N, NQ, NST, NFE,
  751. 5 NJE, NQU
  752. C
  753. EXTERNAL F, JAC
  754. C
  755. C CHECK FOR AN APPARENT INFINITE LOOP
  756. C
  757. C***FIRST EXECUTABLE STATEMENT DEBDF
  758. IF (INFO(1) .EQ. 0) IWORK(LIW) = 0
  759. C
  760. IF (IWORK(LIW).GE. 5) THEN
  761. IF (T .EQ. RWORK(21+NEQ)) THEN
  762. WRITE (XERN3, '(1PE15.6)') T
  763. CALL XERMSG ('SLATEC', 'DEBDF',
  764. * 'AN APPARENT INFINITE LOOP HAS BEEN DETECTED.$$' //
  765. * 'YOU HAVE MADE REPEATED CALLS AT T = ' // XERN3 //
  766. * ' AND THE INTEGRATION HAS NOT ADVANCED. CHECK THE ' //
  767. * 'WAY YOU HAVE SET PARAMETERS FOR THE CALL TO THE ' //
  768. * 'CODE PARTICULARLY INFO(1).', 13, 2)
  769. RETURN
  770. ENDIF
  771. ENDIF
  772. C
  773. IDID = 0
  774. C
  775. C CHECK VALIDITY OF INFO PARAMETERS
  776. C
  777. IF (INFO(1) .NE. 0 .AND. INFO(1) .NE. 1) THEN
  778. WRITE (XERN1, '(I8)') INFO(1)
  779. CALL XERMSG ('SLATEC', 'DEBDF', 'INFO(1) MUST BE SET TO 0 ' //
  780. * 'FOR THE START OF A NEW PROBLEM, AND MUST BE SET TO 1 ' //
  781. * 'FOLLOWING AN INTERRUPTED TASK. YOU ARE ATTEMPTING TO ' //
  782. * 'CONTINUE THE INTEGRATION ILLEGALLY BY CALLING THE ' //
  783. * 'CODE WITH INFO(1) = ' // XERN1, 3, 1)
  784. IDID = -33
  785. ENDIF
  786. C
  787. IF (INFO(2) .NE. 0 .AND. INFO(2) .NE. 1) THEN
  788. WRITE (XERN1, '(I8)') INFO(2)
  789. CALL XERMSG ('SLATEC', 'DEBDF', 'INFO(2) MUST BE 0 OR 1 ' //
  790. * 'INDICATING SCALAR AND VECTOR ERROR TOLERANCES, ' //
  791. * 'RESPECTIVELY. YOU HAVE CALLED THE CODE WITH INFO(2) = ' //
  792. * XERN1, 4, 1)
  793. IDID = -33
  794. ENDIF
  795. C
  796. IF (INFO(3) .NE. 0 .AND. INFO(3) .NE. 1) THEN
  797. WRITE (XERN1, '(I8)') INFO(3)
  798. CALL XERMSG ('SLATEC', 'DEBDF', 'INFO(3) MUST BE 0 OR 1 ' //
  799. * 'INDICATING THE INTERVAL OR INTERMEDIATE-OUTPUT MODE OF ' //
  800. * 'INTEGRATION, RESPECTIVELY. YOU HAVE CALLED THE CODE ' //
  801. * 'WITH INFO(3) = ' // XERN1, 5, 1)
  802. IDID = -33
  803. ENDIF
  804. C
  805. IF (INFO(4) .NE. 0 .AND. INFO(4) .NE. 1) THEN
  806. WRITE (XERN1, '(I8)') INFO(4)
  807. CALL XERMSG ('SLATEC', 'DEBDF', 'INFO(4) MUST BE 0 OR 1 ' //
  808. * 'INDICATING WHETHER OR NOT THE INTEGRATION INTERVAL IS ' //
  809. * 'TO BE RESTRICTED BY A POINT TSTOP. YOU HAVE CALLED ' //
  810. * 'THE CODE WITH INFO(4) = ' // XERN1, 14, 1)
  811. IDID = -33
  812. ENDIF
  813. C
  814. IF (INFO(5) .NE. 0 .AND. INFO(5) .NE. 1) THEN
  815. WRITE (XERN1, '(I8)') INFO(5)
  816. CALL XERMSG ('SLATEC', 'DEBDF', 'INFO(5) MUST BE 0 OR 1 ' //
  817. * 'INDICATING WHETHER THE CODE IS TOLD TO FORM THE ' //
  818. * 'JACOBIAN MATRIX BY NUMERICAL DIFFERENCING OR YOU ' //
  819. * 'PROVIDE A SUBROUTINE TO EVALUATE IT ANALYTICALLY. ' //
  820. * 'YOU HAVE CALLED THE CODE WITH INFO(5) = ' // XERN1, 15, 1)
  821. IDID = -33
  822. ENDIF
  823. C
  824. IF (INFO(6) .NE. 0 .AND. INFO(6) .NE. 1) THEN
  825. WRITE (XERN1, '(I8)') INFO(6)
  826. CALL XERMSG ('SLATEC', 'DEBDF', 'INFO(6) MUST BE 0 OR 1 ' //
  827. * 'INDICATING WHETHER THE CODE IS TOLD TO TREAT THE ' //
  828. * 'JACOBIAN AS A FULL (DENSE) MATRIX OR AS HAVING A ' //
  829. * 'SPECIAL BANDED STRUCTURE. YOU HAVE CALLED THE CODE ' //
  830. * 'WITH INFO(6) = ' // XERN1, 16, 1)
  831. IDID = -33
  832. ENDIF
  833. C
  834. ILRW = NEQ
  835. IF (INFO(6) .NE. 0) THEN
  836. C
  837. C CHECK BANDWIDTH PARAMETERS
  838. C
  839. ML = IWORK(1)
  840. MU = IWORK(2)
  841. ILRW = 2*ML + MU + 1
  842. C
  843. IF (ML.LT.0 .OR. ML.GE.NEQ .OR. MU.LT.0 .OR. MU.GE.NEQ) THEN
  844. WRITE (XERN1, '(I8)') ML
  845. WRITE (XERN2, '(I8)') MU
  846. CALL XERMSG ('SLATEC', 'DEBDF', 'YOU HAVE SET INFO(6) ' //
  847. * '= 1, TELLING THE CODE THAT THE JACOBIAN MATRIX HAS ' //
  848. * 'A SPECIAL BANDED STRUCTURE. HOWEVER, THE LOWER ' //
  849. * '(UPPER) BANDWIDTHS ML (MU) VIOLATE THE CONSTRAINTS ' //
  850. * 'ML,MU .GE. 0 AND ML,MU .LT. NEQ. YOU HAVE CALLED ' //
  851. * 'THE CODE WITH ML = ' // XERN1 // ' AND MU = ' // XERN2,
  852. * 17, 1)
  853. IDID = -33
  854. ENDIF
  855. ENDIF
  856. C
  857. C CHECK LRW AND LIW FOR SUFFICIENT STORAGE ALLOCATION
  858. C
  859. IF (LRW .LT. 250 + (10 + ILRW)*NEQ) THEN
  860. WRITE (XERN1, '(I8)') LRW
  861. IF (INFO(6) .EQ. 0) THEN
  862. CALL XERMSG ('SLATEC', 'DEBDF', 'LENGTH OF ARRAY RWORK ' //
  863. * 'MUST BE AT LEAST 250 + 10*NEQ + NEQ*NEQ.$$' //
  864. * 'YOU HAVE CALLED THE CODE WITH LRW = ' // XERN1, 1, 1)
  865. ELSE
  866. CALL XERMSG ('SLATEC', 'DEBDF', 'LENGTH OF ARRAY RWORK ' //
  867. * 'MUST BE AT LEAST 250 + 10*NEQ + (2*ML+MU+1)*NEQ.$$' //
  868. * 'YOU HAVE CALLED THE CODE WITH LRW = ' // XERN1, 18, 1)
  869. ENDIF
  870. IDID = -33
  871. ENDIF
  872. C
  873. IF (LIW .LT. 56 + NEQ) THEN
  874. WRITE (XERN1, '(I8)') LIW
  875. CALL XERMSG ('SLATEC', 'DEBDF', 'LENGTH OF ARRAY IWORK ' //
  876. * 'BE AT LEAST 56 + NEQ. YOU HAVE CALLED THE CODE WITH ' //
  877. * 'LIW = ' // XERN1, 2, 1)
  878. IDID = -33
  879. ENDIF
  880. C
  881. C COMPUTE THE INDICES FOR THE ARRAYS TO BE STORED IN THE WORK
  882. C ARRAY AND RESTORE COMMON BLOCK DATA
  883. C
  884. ICOMI = 21 + NEQ
  885. IINOUT = ICOMI + 33
  886. C
  887. IYPOUT = 21
  888. ITSTAR = 21 + NEQ
  889. ICOMR = 22 + NEQ
  890. C
  891. IF (INFO(1) .NE. 0) INTOUT = IWORK(IINOUT) .NE. (-1)
  892. C CALL RSCO(RWORK(ICOMR),IWORK(ICOMI))
  893. C
  894. IYH = ICOMR + 218
  895. IEWT = IYH + 6*NEQ
  896. ISAVF = IEWT + NEQ
  897. IACOR = ISAVF + NEQ
  898. IWM = IACOR + NEQ
  899. IDELSN = IWM + 2 + ILRW*NEQ
  900. C
  901. IBEGIN = INFO(1)
  902. ITOL = INFO(2)
  903. IINTEG = INFO(3)
  904. ITSTOP = INFO(4)
  905. IJAC = INFO(5)
  906. IBAND = INFO(6)
  907. RWORK(ITSTAR) = T
  908. C
  909. CALL LSOD(F,NEQ,T,Y,TOUT,RTOL,ATOL,IDID,RWORK(IYPOUT),
  910. 1 RWORK(IYH),RWORK(IYH),RWORK(IEWT),RWORK(ISAVF),
  911. 2 RWORK(IACOR),RWORK(IWM),IWORK(1),JAC,INTOUT,
  912. 3 RWORK(1),RWORK(12),RWORK(IDELSN),RPAR,IPAR)
  913. C
  914. IWORK(IINOUT) = -1
  915. IF (INTOUT) IWORK(IINOUT) = 1
  916. C
  917. IF (IDID .NE. (-2)) IWORK(LIW) = IWORK(LIW) + 1
  918. IF (T .NE. RWORK(ITSTAR)) IWORK(LIW) = 0
  919. C CALL SVCO(RWORK(ICOMR),IWORK(ICOMI))
  920. RWORK(11) = H
  921. RWORK(13) = TN
  922. INFO(1) = IBEGIN
  923. C
  924. RETURN
  925. END