dlpdoc.f 25 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460
  1. *DECK DLPDOC
  2. SUBROUTINE DLPDOC
  3. C***BEGIN PROLOGUE DLPDOC
  4. C***PURPOSE Sparse Linear Algebra Package Version 2.0.2 Documentation.
  5. C Routines to solve large sparse symmetric and nonsymmetric
  6. C positive definite linear systems, Ax = b, using precondi-
  7. C tioned iterative methods.
  8. C***LIBRARY SLATEC (SLAP)
  9. C***CATEGORY D2A4, D2B4, Z
  10. C***TYPE DOUBLE PRECISION (SLPDOC-S, DLPDOC-D)
  11. C***KEYWORDS BICONJUGATE GRADIENT SQUARED, DOCUMENTATION,
  12. C GENERALIZED MINIMUM RESIDUAL, ITERATIVE IMPROVEMENT,
  13. C NORMAL EQUATIONS, ORTHOMIN,
  14. C PRECONDITIONED CONJUGATE GRADIENT, SLAP,
  15. C SPARSE ITERATIVE METHODS
  16. C***AUTHOR Seager, Mark. K., (LLNL)
  17. C User Systems Division
  18. C Lawrence Livermore National Laboratory
  19. C PO BOX 808, L-60
  20. C Livermore, CA 94550
  21. C (FTS) 543-3141, (510) 423-3141
  22. C seager@llnl.gov
  23. C***DESCRIPTION
  24. C The
  25. C Sparse Linear Algebra Package
  26. C Double Precision Routines
  27. C
  28. C @@@@@@@ @ @@@ @@@@@@@@
  29. C @ @ @ @ @ @ @
  30. C @ @ @ @ @ @
  31. C @@@@@@@ @ @ @ @@@@@@@@
  32. C @ @ @@@@@@@@@ @
  33. C @ @ @ @ @ @
  34. C @@@@@@@ @@@@@@@@@ @ @ @
  35. C
  36. C @ @ @@@@@@@ @@@@@
  37. C @ @ @ @ @ @@
  38. C @ @ @@@@@@@ @ @@ @ @ @ @
  39. C @ @ @ @ @@ @ @@@@@@ @ @ @
  40. C @ @ @@@@@@@@@ @ @ @ @ @
  41. C @ @ @ @ @ @@@ @@ @
  42. C @@@ @@@@@@@ @ @@@@@@@@@ @@@ @@@@@
  43. C
  44. C
  45. C =================================================================
  46. C ========================== Introduction =========================
  47. C =================================================================
  48. C This package was originally derived from a set of iterative
  49. C routines written by Anne Greenbaum, as announced in "Routines
  50. C for Solving Large Sparse Linear Systems", Tentacle, Lawrence
  51. C Livermore National Laboratory, Livermore Computing Center
  52. C (January 1986), pp 15-21.
  53. C
  54. C This document contains the specifications for the SLAP Version
  55. C 2.0 package, a Fortran 77 package for the solution of large
  56. C sparse linear systems, Ax = b, via preconditioned iterative
  57. C methods. Included in this package are "core" routines to do
  58. C Iterative Refinement (Jacobi's method), Conjugate Gradient,
  59. C Conjugate Gradient on the normal equations, AA'y = b, (where x =
  60. C A'y and A' denotes the transpose of A), BiConjugate Gradient,
  61. C BiConjugate Gradient Squared, Orthomin and Generalized Minimum
  62. C Residual Iteration. These "core" routines do not require a
  63. C "fixed" data structure for storing the matrix A and the
  64. C preconditioning matrix M. The user is free to choose any
  65. C structure that facilitates efficient solution of the problem at
  66. C hand. The drawback to this approach is that the user must also
  67. C supply at least two routines (MATVEC and MSOLVE, say). MATVEC
  68. C must calculate, y = Ax, given x and the user's data structure for
  69. C A. MSOLVE must solve, r = Mz, for z (*NOT* r) given r and the
  70. C user's data structure for M (or its inverse). The user should
  71. C choose M so that inv(M)*A is approximately the identity and the
  72. C solution step r = Mz is "easy" to solve. For some of the "core"
  73. C routines (Orthomin, BiConjugate Gradient and Conjugate Gradient
  74. C on the normal equations) the user must also supply a matrix
  75. C transpose times vector routine (MTTVEC, say) and (possibly,
  76. C depending on the "core" method) a routine that solves the
  77. C transpose of the preconditioning step (MTSOLV, say).
  78. C Specifically, MTTVEC is a routine which calculates y = A'x, given
  79. C x and the user's data structure for A (A' is the transpose of A).
  80. C MTSOLV is a routine which solves the system r = M'z for z given r
  81. C and the user's data structure for M.
  82. C
  83. C This process of writing the matrix vector operations can be time
  84. C consuming and error prone. To alleviate these problems we have
  85. C written drivers for the "core" methods that assume the user
  86. C supplies one of two specific data structures (SLAP Triad and SLAP
  87. C Column format), see below. Utilizing these data structures we
  88. C have augmented each "core" method with two preconditioners:
  89. C Diagonal Scaling and Incomplete Factorization. Diagonal scaling
  90. C is easy to implement, vectorizes very well and for problems that
  91. C are not too ill-conditioned reduces the number of iterations
  92. C enough to warrant its use. On the other hand, an Incomplete
  93. C factorization (Incomplete Cholesky for symmetric systems and
  94. C Incomplete LU for nonsymmetric systems) may take much longer to
  95. C calculate, but it reduces the iteration count (for most problems)
  96. C significantly. Our implementations of IC and ILU vectorize for
  97. C machines with hardware gather scatter, but the vector lengths can
  98. C be quite short if the number of non-zeros in a column is not
  99. C large.
  100. C
  101. C =================================================================
  102. C ==================== Supplied Data Structures ===================
  103. C =================================================================
  104. C The following describes the data structures supplied with the
  105. C package: SLAP Triad and Column formats.
  106. C
  107. C ====================== S L A P Triad format =====================
  108. C
  109. C In the SLAP Triad format only the non-zeros are stored. They may
  110. C appear in *ANY* order. The user supplies three arrays of length
  111. C NELT, where NELT is the number of non-zeros in the matrix:
  112. C (IA(NELT), JA(NELT), A(NELT)). If the matrix is symmetric then
  113. C one need only store the lower triangle (including the diagonal)
  114. C and NELT would be the corresponding number of non-zeros stored.
  115. C For each non-zero the user puts the row and column index of that
  116. C matrix element in the IA and JA arrays. The value of the
  117. C non-zero matrix element is placed in the corresponding location
  118. C of the A array. This is an extremely easy data structure to
  119. C generate. On the other hand, it is not very efficient on vector
  120. C computers for the iterative solution of linear systems. Hence,
  121. C SLAP changes this input data structure to the SLAP Column format
  122. C for the iteration (but does not change it back).
  123. C
  124. C Here is an example of the SLAP Triad storage format for a
  125. C nonsymmetric 5x5 Matrix. NELT=11. Recall that the entries may
  126. C appear in any order.
  127. C
  128. C 5x5 Matrix SLAP Triad format for 5x5 matrix on left.
  129. C 1 2 3 4 5 6 7 8 9 10 11
  130. C |11 12 0 0 15| A: 51 12 11 33 15 53 55 22 35 44 21
  131. C |21 22 0 0 0| IA: 5 1 1 3 1 5 5 2 3 4 2
  132. C | 0 0 33 0 35| JA: 1 2 1 3 5 3 5 2 5 4 1
  133. C | 0 0 0 44 0|
  134. C |51 0 53 0 55|
  135. C
  136. C ====================== S L A P Column format ====================
  137. C
  138. C In the SLAP Column format the non-zeros are stored counting down
  139. C columns (except for the diagonal entry, which must appear first
  140. C in each "column") and are stored in the double precision array A.
  141. C In other words, for each column in the matrix first put the
  142. C diagonal entry in A. Then put in the other non-zero elements
  143. C going down the column (except the diagonal) in order. The IA
  144. C array holds the row index for each non-zero. The JA array holds
  145. C the offsets into the IA, A arrays for the beginning of each
  146. C column. That is, IA(JA(ICOL)), A(JA(ICOL)) are the first elements
  147. C of the ICOL-th column in IA and A, and IA(JA(ICOL+1)-1),
  148. C A(JA(ICOL+1)-1) are the last elements of the ICOL-th column. Note
  149. C that we always have JA(N+1) = NELT+1, where N is the number of
  150. C columns in the matrix and NELT is the number of non-zeros in the
  151. C matrix. If the matrix is symmetric one need only store the lower
  152. C triangle (including the diagonal) and NELT would be the corre-
  153. C sponding number of non-zeros stored.
  154. C
  155. C Here is an example of the SLAP Column storage format for a
  156. C nonsymmetric 5x5 Matrix (in the A and IA arrays '|' denotes the
  157. C end of a column):
  158. C
  159. C 5x5 Matrix SLAP Column format for 5x5 matrix on left.
  160. C 1 2 3 4 5 6 7 8 9 10 11
  161. C |11 12 0 0 15| A: 11 21 51 | 22 12 | 33 53 | 44 | 55 15 35
  162. C |21 22 0 0 0| IA: 1 2 5 | 2 1 | 3 5 | 4 | 5 1 3
  163. C | 0 0 33 0 35| JA: 1 4 6 8 9 12
  164. C | 0 0 0 44 0|
  165. C |51 0 53 0 55|
  166. C
  167. C =================================================================
  168. C ====================== Which Method To Use ======================
  169. C =================================================================
  170. C
  171. C BACKGROUND
  172. C In solving a large sparse linear system Ax = b using an iterative
  173. C method, it is not necessary to actually store the matrix A.
  174. C Rather, what is needed is a procedure for multiplying the matrix
  175. C A times a given vector y to obtain the matrix-vector product, Ay.
  176. C SLAP has been written to take advantage of this fact. The higher
  177. C level routines in the package require storage only of the non-zero
  178. C elements of A (and their positions), and even this can be
  179. C avoided, if the user writes his own subroutine for multiplying
  180. C the matrix times a vector and calls the lower-level iterative
  181. C routines in the package.
  182. C
  183. C If the matrix A is ill-conditioned, then most iterative methods
  184. C will be slow to converge (if they converge at all!). To improve
  185. C the convergence rate, one may use a "matrix splitting," or,
  186. C "preconditioning matrix," say, M. It is then necessary to solve,
  187. C at each iteration, a linear system with coefficient matrix M. A
  188. C good preconditioner M should have two properties: (1) M should
  189. C "approximate" A, in the sense that the matrix inv(M)*A (or some
  190. C variant thereof) is better conditioned than the original matrix
  191. C A; and (2) linear systems with coefficient matrix M should be
  192. C much easier to solve than the original system with coefficient
  193. C matrix A. Preconditioning routines in the SLAP package are
  194. C separate from the iterative routines, so that any of the
  195. C preconditioners provided in the package, or one that the user
  196. C codes himself, can be used with any of the iterative routines.
  197. C
  198. C CHOICE OF PRECONDITIONER
  199. C If you willing to live with either the SLAP Triad or Column
  200. C matrix data structure you can then choose one of two types of
  201. C preconditioners to use: diagonal scaling or incomplete
  202. C factorization. To choose between these two methods requires
  203. C knowing something about the computer you're going to run these
  204. C codes on and how well incomplete factorization approximates the
  205. C inverse of your matrix.
  206. C
  207. C Let us suppose you have a scalar machine. Then, unless the
  208. C incomplete factorization is very, very poor this is *GENERALLY*
  209. C the method to choose. It will reduce the number of iterations
  210. C significantly and is not all that expensive to compute. So if
  211. C you have just one linear system to solve and "just want to get
  212. C the job done" then try incomplete factorization first. If you
  213. C are thinking of integrating some SLAP iterative method into your
  214. C favorite "production code" then try incomplete factorization
  215. C first, but also check to see that diagonal scaling is indeed
  216. C slower for a large sample of test problems.
  217. C
  218. C Let us now suppose you have a vector computer with hardware
  219. C gather/scatter support (Cray X-MP, Y-MP, SCS-40 or Cyber 205, ETA
  220. C 10, ETA Piper, Convex C-1, etc.). Then it is much harder to
  221. C choose between the two methods. The versions of incomplete
  222. C factorization in SLAP do in fact vectorize, but have short vector
  223. C lengths and the factorization step is relatively more expensive.
  224. C Hence, for most problems (i.e., unless your problem is ill
  225. C conditioned, sic!) diagonal scaling is faster, with its very
  226. C fast set up time and vectorized (with long vectors)
  227. C preconditioning step (even though it may take more iterations).
  228. C If you have several systems (or right hand sides) to solve that
  229. C can utilize the same preconditioner then the cost of the
  230. C incomplete factorization can be amortized over these several
  231. C solutions. This situation gives more advantage to the incomplete
  232. C factorization methods. If you have a vector machine without
  233. C hardware gather/scatter (Cray 1, Cray 2 & Cray 3) then the
  234. C advantages for incomplete factorization are even less.
  235. C
  236. C If you're trying to shoehorn SLAP into your favorite "production
  237. C code" and can not easily generate either the SLAP Triad or Column
  238. C format then you are left to your own devices in terms of
  239. C preconditioning. Also, you may find that the preconditioners
  240. C supplied with SLAP are not sufficient for your problem. In this
  241. C situation we would recommend that you talk with a numerical
  242. C analyst versed in iterative methods about writing other
  243. C preconditioning subroutines (e.g., polynomial preconditioning,
  244. C shifted incomplete factorization, SOR or SSOR iteration). You
  245. C can always "roll your own" by using the "core" iterative methods
  246. C and supplying your own MSOLVE and MATVEC (and possibly MTSOLV and
  247. C MTTVEC) routines.
  248. C
  249. C SYMMETRIC SYSTEMS
  250. C If your matrix is symmetric then you would want to use one of the
  251. C symmetric system solvers. If your system is also positive
  252. C definite, (Ax,x) (Ax dot product with x) is positive for all
  253. C non-zero vectors x, then use Conjugate Gradient (DCG, DSDCG,
  254. C DSICSG). If you're not sure it's SPD (symmetric and Positive
  255. C Definite) then try DCG anyway and if it works, fine. If you're
  256. C sure your matrix is not positive definite then you may want to
  257. C try the iterative refinement methods (DIR) or the GMRES code
  258. C (DGMRES) if DIR converges too slowly.
  259. C
  260. C NONSYMMETRIC SYSTEMS
  261. C This is currently an area of active research in numerical
  262. C analysis and there are new strategies being developed.
  263. C Consequently take the following advice with a grain of salt. If
  264. C you matrix is positive definite, (Ax,x) (Ax dot product with x
  265. C is positive for all non-zero vectors x), then you can use any of
  266. C the methods for nonsymmetric systems (Orthomin, GMRES,
  267. C BiConjugate Gradient, BiConjugate Gradient Squared and Conjugate
  268. C Gradient applied to the normal equations). If your system is not
  269. C too ill conditioned then try BiConjugate Gradient Squared (BCGS)
  270. C or GMRES (DGMRES). Both of these methods converge very quickly
  271. C and do not require A' or M' (' denotes transpose) information.
  272. C DGMRES does require some additional storage, though. If the
  273. C system is very ill conditioned or nearly positive indefinite
  274. C ((Ax,x) is positive, but may be very small), then GMRES should
  275. C be the first choice, but try the other methods if you have to
  276. C fine tune the solution process for a "production code". If you
  277. C have a great preconditioner for the normal equations (i.e., M is
  278. C an approximation to the inverse of AA' rather than just A) then
  279. C this is not a bad route to travel. Old wisdom would say that the
  280. C normal equations are a disaster (since it squares the condition
  281. C number of the system and DCG convergence is linked to this number
  282. C of infamy), but some preconditioners (like incomplete
  283. C factorization) can reduce the condition number back below that of
  284. C the original system.
  285. C
  286. C =================================================================
  287. C ======================= Naming Conventions ======================
  288. C =================================================================
  289. C SLAP iterative methods, matrix vector and preconditioner
  290. C calculation routines follow a naming convention which, when
  291. C understood, allows one to determine the iterative method and data
  292. C structure(s) used. The subroutine naming convention takes the
  293. C following form:
  294. C P[S][M]DESC
  295. C where
  296. C P stands for the precision (or data type) of the routine and
  297. C is required in all names,
  298. C S denotes whether or not the routine requires the SLAP Triad
  299. C or Column format (it does if the second letter of the name
  300. C is S and does not otherwise),
  301. C M stands for the type of preconditioner used (only appears
  302. C in drivers for "core" routines), and
  303. C DESC is some number of letters describing the method or purpose
  304. C of the routine. The following is a list of the "DESC"
  305. C fields for iterative methods and their meaning:
  306. C BCG,BC: BiConjugate Gradient
  307. C CG: Conjugate Gradient
  308. C CGN,CN: Conjugate Gradient on the Normal equations
  309. C CGS,CS: biConjugate Gradient Squared
  310. C GMRES,GMR,GM: Generalized Minimum RESidual
  311. C IR,R: Iterative Refinement
  312. C JAC: JACobi's method
  313. C GS: Gauss-Seidel
  314. C OMN,OM: OrthoMiN
  315. C
  316. C In the double precision version of SLAP, all routine names start
  317. C with a D. The brackets around the S and M designate that these
  318. C fields are optional.
  319. C
  320. C Here are some examples of the routines:
  321. C 1) DBCG: Double precision BiConjugate Gradient "core" routine.
  322. C One can deduce that this is a "core" routine, because the S and
  323. C M fields are missing and BiConjugate Gradient is an iterative
  324. C method.
  325. C 2) DSDBCG: Double precision, SLAP data structure BCG with Diagonal
  326. C scaling.
  327. C 3) DSLUBC: Double precision, SLAP data structure BCG with incom-
  328. C plete LU factorization as the preconditioning.
  329. C 4) DCG: Double precision Conjugate Gradient "core" routine.
  330. C 5) DSDCG: Double precision, SLAP data structure Conjugate Gradient
  331. C with Diagonal scaling.
  332. C 6) DSICCG: Double precision, SLAP data structure Conjugate Gra-
  333. C dient with Incomplete Cholesky factorization preconditioning.
  334. C
  335. C
  336. C =================================================================
  337. C ===================== USER CALLABLE ROUTINES ====================
  338. C =================================================================
  339. C The following is a list of the "user callable" SLAP routines and
  340. C their one line descriptions. The headers denote the file names
  341. C where the routines can be found, as distributed for UNIX systems.
  342. C
  343. C Note: Each core routine, DXXX, has a corresponding stop routine,
  344. C ISDXXX. If the stop routine does not have the specific stop
  345. C test the user requires (e.g., weighted infinity norm), then
  346. C the user should modify the source for ISDXXX accordingly.
  347. C
  348. C ============================= dir.f =============================
  349. C DIR: Preconditioned Iterative Refinement Sparse Ax = b Solver.
  350. C DSJAC: Jacobi's Method Iterative Sparse Ax = b Solver.
  351. C DSGS: Gauss-Seidel Method Iterative Sparse Ax = b Solver.
  352. C DSILUR: Incomplete LU Iterative Refinement Sparse Ax = b Solver.
  353. C
  354. C ============================= dcg.f =============================
  355. C DCG: Preconditioned Conjugate Gradient Sparse Ax=b Solver.
  356. C DSDCG: Diagonally Scaled Conjugate Gradient Sparse Ax=b Solver.
  357. C DSICCG: Incomplete Cholesky Conjugate Gradient Sparse Ax=b Solver.
  358. C
  359. C ============================= dcgn.f ============================
  360. C DCGN: Preconditioned CG Sparse Ax=b Solver for Normal Equations.
  361. C DSDCGN: Diagonally Scaled CG Sparse Ax=b Solver for Normal Eqn's.
  362. C DSLUCN: Incomplete LU CG Sparse Ax=b Solver for Normal Equations.
  363. C
  364. C ============================= dbcg.f ============================
  365. C DBCG: Preconditioned BiConjugate Gradient Sparse Ax = b Solver.
  366. C DSDBCG: Diagonally Scaled BiConjugate Gradient Sparse Ax=b Solver.
  367. C DSLUBC: Incomplete LU BiConjugate Gradient Sparse Ax=b Solver.
  368. C
  369. C ============================= dcgs.f ============================
  370. C DCGS: Preconditioned BiConjugate Gradient Squared Ax=b Solver.
  371. C DSDCGS: Diagonally Scaled CGS Sparse Ax=b Solver.
  372. C DSLUCS: Incomplete LU BiConjugate Gradient Squared Ax=b Solver.
  373. C
  374. C ============================= domn.f ============================
  375. C DOMN: Preconditioned Orthomin Sparse Iterative Ax=b Solver.
  376. C DSDOMN: Diagonally Scaled Orthomin Sparse Iterative Ax=b Solver.
  377. C DSLUOM: Incomplete LU Orthomin Sparse Iterative Ax=b Solver.
  378. C
  379. C ============================ dgmres.f ===========================
  380. C DGMRES: Preconditioned GMRES Iterative Sparse Ax=b Solver.
  381. C DSDGMR: Diagonally Scaled GMRES Iterative Sparse Ax=b Solver.
  382. C DSLUGM: Incomplete LU GMRES Iterative Sparse Ax=b Solver.
  383. C
  384. C ============================ dmset.f ============================
  385. C The following routines are used to set up preconditioners.
  386. C
  387. C DSDS: Diagonal Scaling Preconditioner SLAP Set Up.
  388. C DSDSCL: Diagonally Scales/Unscales a SLAP Column Matrix.
  389. C DSD2S: Diagonal Scaling Preconditioner SLAP Normal Eqns Set Up.
  390. C DS2LT: Lower Triangle Preconditioner SLAP Set Up.
  391. C DSICS: Incomplete Cholesky Decomp. Preconditioner SLAP Set Up.
  392. C DSILUS: Incomplete LU Decomposition Preconditioner SLAP Set Up.
  393. C
  394. C ============================ dmvops.f ===========================
  395. C Most of the incomplete factorization (LL' and LDU) solvers
  396. C in this file require an intermediate routine to translate
  397. C from the SLAP MSOLVE(N, R, Z, NELT, IA, JA, A, ISYM, RWORK,
  398. C IWORK) calling convention to the calling sequence required
  399. C by the solve routine. This generally is accomplished by
  400. C fishing out pointers to the preconditioner (stored in RWORK)
  401. C from the IWORK array and then making a call to the routine
  402. C that actually does the backsolve.
  403. C
  404. C DSMV: SLAP Column Format Sparse Matrix Vector Product.
  405. C DSMTV: SLAP Column Format Sparse Matrix (transpose) Vector Prod.
  406. C DSDI: Diagonal Matrix Vector Multiply.
  407. C DSLI: SLAP MSOLVE for Lower Triangle Matrix (set up for DSLI2).
  408. C DSLI2: Lower Triangle Matrix Backsolve.
  409. C DSLLTI: SLAP MSOLVE for LDL' (IC) Fact. (set up for DLLTI2).
  410. C DLLTI2: Backsolve routine for LDL' Factorization.
  411. C DSLUI: SLAP MSOLVE for LDU Factorization (set up for DSLUI2).
  412. C DSLUI2: SLAP Backsolve for LDU Factorization.
  413. C DSLUTI: SLAP MTSOLV for LDU Factorization (set up for DSLUI4).
  414. C DSLUI4: SLAP Backsolve for LDU Factorization.
  415. C DSMMTI: SLAP MSOLVE for LDU Fact of Normal Eq (set up for DSMMI2).
  416. C DSMMI2: SLAP Backsolve for LDU Factorization of Normal Equations.
  417. C
  418. C =========================== dlaputil.f ==========================
  419. C The following utility routines are useful additions to SLAP.
  420. C
  421. C DBHIN: Read Sparse Linear System in the Boeing/Harwell Format.
  422. C DCHKW: SLAP WORK/IWORK Array Bounds Checker.
  423. C DCPPLT: Printer Plot of SLAP Column Format Matrix.
  424. C DS2Y: SLAP Triad to SLAP Column Format Converter.
  425. C QS2I1D: Quick Sort Integer array, moving integer and DP arrays.
  426. C (Used by DS2Y.)
  427. C DTIN: Read in SLAP Triad Format Linear System.
  428. C DTOUT: Write out SLAP Triad Format Linear System.
  429. C
  430. C
  431. C***REFERENCES 1. Mark K. Seager, A SLAP for the Masses, in
  432. C G. F. Carey, Ed., Parallel Supercomputing: Methods,
  433. C Algorithms and Applications, Wiley, 1989, pp.135-155.
  434. C***ROUTINES CALLED (NONE)
  435. C***REVISION HISTORY (YYMMDD)
  436. C 890404 DATE WRITTEN
  437. C 890404 Previous REVISION DATE
  438. C 890915 Made changes requested at July 1989 CML Meeting. (MKS)
  439. C 890921 Removed TeX from comments. (FNF)
  440. C 890922 Numerous changes to prologue to make closer to SLATEC
  441. C standard. (FNF)
  442. C 890929 Numerous changes to reduce SP/DP differences. (FNF)
  443. C -----( This produced Version 2.0.1. )-----
  444. C 891003 Rearranged list of user callable routines to agree with
  445. C order in source deck. (FNF)
  446. C 891004 Updated reference.
  447. C 910411 Prologue converted to Version 4.0 format. (BAB)
  448. C -----( This produced Version 2.0.2. )-----
  449. C 910506 Minor improvements to prologue. (FNF)
  450. C 920511 Added complete declaration section. (WRB)
  451. C 920929 Corrected format of reference. (FNF)
  452. C 921019 Improved one-line descriptions, reordering some. (FNF)
  453. C***END PROLOGUE DLPDOC
  454. C***FIRST EXECUTABLE STATEMENT DLPDOC
  455. C
  456. C This is a *DUMMY* subroutine and should never be called.
  457. C
  458. RETURN
  459. C------------- LAST LINE OF DLPDOC FOLLOWS -----------------------------
  460. END