dsdgmr.f 19 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386
  1. *DECK DSDGMR
  2. SUBROUTINE DSDGMR (N, B, X, NELT, IA, JA, A, ISYM, NSAVE, ITOL,
  3. + TOL, ITMAX, ITER, ERR, IERR, IUNIT, RWORK, LENW, IWORK, LENIW)
  4. C***BEGIN PROLOGUE DSDGMR
  5. C***PURPOSE Diagonally scaled GMRES iterative sparse Ax=b solver.
  6. C This routine uses the generalized minimum residual
  7. C (GMRES) method with diagonal scaling to solve possibly
  8. C non-symmetric linear systems of the form: Ax = b.
  9. C***LIBRARY SLATEC (SLAP)
  10. C***CATEGORY D2A4, D2B4
  11. C***TYPE DOUBLE PRECISION (SSDGMR-S, DSDGMR-D)
  12. C***KEYWORDS GENERALIZED MINIMUM RESIDUAL, ITERATIVE PRECONDITION,
  13. C NON-SYMMETRIC LINEAR SYSTEM, SLAP, SPARSE
  14. C***AUTHOR Brown, Peter, (LLNL), pnbrown@llnl.gov
  15. C Hindmarsh, Alan, (LLNL), alanh@llnl.gov
  16. C Seager, Mark K., (LLNL), seager@llnl.gov
  17. C Lawrence Livermore National Laboratory
  18. C PO Box 808, L-60
  19. C Livermore, CA 94550 (510) 423-3141
  20. C***DESCRIPTION
  21. C
  22. C *Usage:
  23. C INTEGER N, NELT, IA(NELT), JA(NELT), ISYM, NSAVE, ITOL
  24. C INTEGER ITMAX, ITER, IERR, IUNIT, LENW, IWORK(LENIW), LENIW
  25. C DOUBLE PRECISION B(N), X(N), A(NELT), TOL, ERR, RWORK(LENW)
  26. C
  27. C CALL DSDGMR(N, B, X, NELT, IA, JA, A, ISYM, NSAVE,
  28. C $ ITOL, TOL, ITMAX, ITER, ERR, IERR, IUNIT,
  29. C $ RWORK, LENW, IWORK, LENIW)
  30. C
  31. C *Arguments:
  32. C N :IN Integer.
  33. C Order of the Matrix.
  34. C B :IN Double Precision B(N).
  35. C Right-hand side vector.
  36. C X :INOUT Double Precision X(N).
  37. C On input X is your initial guess for solution vector.
  38. C On output X is the final approximate solution.
  39. C NELT :IN Integer.
  40. C Number of Non-Zeros stored in A.
  41. C IA :IN Integer IA(NELT).
  42. C JA :IN Integer JA(NELT).
  43. C A :IN Double Precision A(NELT).
  44. C These arrays should hold the matrix A in either the SLAP
  45. C Triad format or the SLAP Column format. See "Description",
  46. C below. If the SLAP Triad format is chosen it is changed
  47. C internally to the SLAP Column format.
  48. C ISYM :IN Integer.
  49. C Flag to indicate symmetric storage format.
  50. C If ISYM=0, all non-zero entries of the matrix are stored.
  51. C If ISYM=1, the matrix is symmetric, and only the upper
  52. C or lower triangle of the matrix is stored.
  53. C NSAVE :IN Integer.
  54. C Number of direction vectors to save and orthogonalize against.
  55. C Must be greater than 1.
  56. C ITOL :IN Integer.
  57. C Flag to indicate the type of convergence criterion used.
  58. C ITOL=0 Means the iteration stops when the test described
  59. C below on the residual RL is satisfied. This is
  60. C the "Natural Stopping Criteria" for this routine.
  61. C Other values of ITOL cause extra, otherwise
  62. C unnecessary, computation per iteration and are
  63. C therefore much less efficient. See ISDGMR (the
  64. C stop test routine) for more information.
  65. C ITOL=1 Means the iteration stops when the first test
  66. C described below on the residual RL is satisfied,
  67. C and there is either right or no preconditioning
  68. C being used.
  69. C ITOL=2 Implies that the user is using left
  70. C preconditioning, and the second stopping criterion
  71. C below is used.
  72. C ITOL=3 Means the iteration stops when the third test
  73. C described below on Minv*Residual is satisfied, and
  74. C there is either left or no preconditioning begin
  75. C used.
  76. C ITOL=11 is often useful for checking and comparing
  77. C different routines. For this case, the user must
  78. C supply the "exact" solution or a very accurate
  79. C approximation (one with an error much less than
  80. C TOL) through a common block,
  81. C COMMON /DSLBLK/ SOLN( )
  82. C If ITOL=11, iteration stops when the 2-norm of the
  83. C difference between the iterative approximation and
  84. C the user-supplied solution divided by the 2-norm
  85. C of the user-supplied solution is less than TOL.
  86. C Note that this requires the user to set up the
  87. C "COMMON /DSLBLK/ SOLN(LENGTH)" in the calling
  88. C routine. The routine with this declaration should
  89. C be loaded before the stop test so that the correct
  90. C length is used by the loader. This procedure is
  91. C not standard Fortran and may not work correctly on
  92. C your system (although it has worked on every
  93. C system the authors have tried). If ITOL is not 11
  94. C then this common block is indeed standard Fortran.
  95. C TOL :INOUT Double Precision.
  96. C Convergence criterion, as described below. If TOL is set
  97. C to zero on input, then a default value of 500*(the smallest
  98. C positive magnitude, machine epsilon) is used.
  99. C ITMAX :IN Integer.
  100. C Maximum number of iterations. This routine uses the default
  101. C of NRMAX = ITMAX/NSAVE to determine when each restart
  102. C should occur. See the description of NRMAX and MAXL in
  103. C DGMRES for a full and frightfully interesting discussion of
  104. C this topic.
  105. C ITER :OUT Integer.
  106. C Number of iterations required to reach convergence, or
  107. C ITMAX+1 if convergence criterion could not be achieved in
  108. C ITMAX iterations.
  109. C ERR :OUT Double Precision.
  110. C Error estimate of error in final approximate solution, as
  111. C defined by ITOL. Letting norm() denote the Euclidean
  112. C norm, ERR is defined as follows...
  113. C If ITOL=0, then ERR = norm(SB*(B-A*X(L)))/norm(SB*B),
  114. C for right or no preconditioning, and
  115. C ERR = norm(SB*(M-inverse)*(B-A*X(L)))/
  116. C norm(SB*(M-inverse)*B),
  117. C for left preconditioning.
  118. C If ITOL=1, then ERR = norm(SB*(B-A*X(L)))/norm(SB*B),
  119. C since right or no preconditioning
  120. C being used.
  121. C If ITOL=2, then ERR = norm(SB*(M-inverse)*(B-A*X(L)))/
  122. C norm(SB*(M-inverse)*B),
  123. C since left preconditioning is being
  124. C used.
  125. C If ITOL=3, then ERR = Max |(Minv*(B-A*X(L)))(i)/x(i)|
  126. C i=1,n
  127. C If ITOL=11, then ERR = norm(SB*(X(L)-SOLN))/norm(SB*SOLN).
  128. C IERR :OUT Integer.
  129. C Return error flag.
  130. C IERR = 0 => All went well.
  131. C IERR = 1 => Insufficient storage allocated for
  132. C RGWK or IGWK.
  133. C IERR = 2 => Routine DPIGMR failed to reduce the norm
  134. C of the current residual on its last call,
  135. C and so the iteration has stalled. In
  136. C this case, X equals the last computed
  137. C approximation. The user must either
  138. C increase MAXL, or choose a different
  139. C initial guess.
  140. C IERR =-1 => Insufficient length for RGWK array.
  141. C IGWK(6) contains the required minimum
  142. C length of the RGWK array.
  143. C IERR =-2 => Inconsistent ITOL and JPRE values.
  144. C For IERR <= 2, RGWK(1) = RHOL, which is the norm on the
  145. C left-hand-side of the relevant stopping test defined
  146. C below associated with the residual for the current
  147. C approximation X(L).
  148. C IUNIT :IN Integer.
  149. C Unit number on which to write the error at each iteration,
  150. C if this is desired for monitoring convergence. If unit
  151. C number is 0, no writing will occur.
  152. C RWORK :WORK Double Precision RWORK(LENW).
  153. C Double Precision array of size LENW.
  154. C LENW :IN Integer.
  155. C Length of the double precision workspace, RWORK.
  156. C LENW >= 1 + N*(NSAVE+7) + NSAVE*(NSAVE+3).
  157. C For the recommended values of NSAVE (10), RWORK has size at
  158. C least 131 + 17*N.
  159. C IWORK :INOUT Integer IWORK(USER DEFINED >= 30).
  160. C Used to hold pointers into the RWORK array.
  161. C Upon return the following locations of IWORK hold information
  162. C which may be of use to the user:
  163. C IWORK(9) Amount of Integer workspace actually used.
  164. C IWORK(10) Amount of Double Precision workspace actually used.
  165. C LENIW :IN Integer.
  166. C Length of the integer workspace IWORK. LENIW >= 30.
  167. C
  168. C *Description:
  169. C DSDGMR solves a linear system A*X = B rewritten in the form:
  170. C
  171. C (SB*A*(M-inverse)*(SX-inverse))*(SX*M*X) = SB*B,
  172. C
  173. C with right preconditioning, or
  174. C
  175. C (SB*(M-inverse)*A*(SX-inverse))*(SX*X) = SB*(M-inverse)*B,
  176. C
  177. C with left preconditioning, where A is an n-by-n double precision
  178. C matrix, X and B are N-vectors, SB and SX are diagonal scaling
  179. C matrices, and M is the diagonal of A. It uses
  180. C preconditioned Krylov subpace methods based on the
  181. C generalized minimum residual method (GMRES). This routine
  182. C is a driver routine which assumes a SLAP matrix data
  183. C structure and sets up the necessary information to do
  184. C diagonal preconditioning and calls the main GMRES routine
  185. C DGMRES for the solution of the linear system. DGMRES
  186. C optionally performs either the full orthogonalization
  187. C version of the GMRES algorithm or an incomplete variant of
  188. C it. Both versions use restarting of the linear iteration by
  189. C default, although the user can disable this feature.
  190. C
  191. C The GMRES algorithm generates a sequence of approximations
  192. C X(L) to the true solution of the above linear system. The
  193. C convergence criteria for stopping the iteration is based on
  194. C the size of the scaled norm of the residual R(L) = B -
  195. C A*X(L). The actual stopping test is either:
  196. C
  197. C norm(SB*(B-A*X(L))) .le. TOL*norm(SB*B),
  198. C
  199. C for right preconditioning, or
  200. C
  201. C norm(SB*(M-inverse)*(B-A*X(L))) .le.
  202. C TOL*norm(SB*(M-inverse)*B),
  203. C
  204. C for left preconditioning, where norm() denotes the Euclidean
  205. C norm, and TOL is a positive scalar less than one input by
  206. C the user. If TOL equals zero when DSDGMR is called, then a
  207. C default value of 500*(the smallest positive magnitude,
  208. C machine epsilon) is used. If the scaling arrays SB and SX
  209. C are used, then ideally they should be chosen so that the
  210. C vectors SX*X(or SX*M*X) and SB*B have all their components
  211. C approximately equal to one in magnitude. If one wants to
  212. C use the same scaling in X and B, then SB and SX can be the
  213. C same array in the calling program.
  214. C
  215. C The following is a list of the other routines and their
  216. C functions used by GMRES:
  217. C DGMRES Contains the matrix structure independent driver
  218. C routine for GMRES.
  219. C DPIGMR Contains the main iteration loop for GMRES.
  220. C DORTH Orthogonalizes a new vector against older basis vectors.
  221. C DHEQR Computes a QR decomposition of a Hessenberg matrix.
  222. C DHELS Solves a Hessenberg least-squares system, using QR
  223. C factors.
  224. C RLCALC Computes the scaled residual RL.
  225. C XLCALC Computes the solution XL.
  226. C ISDGMR User-replaceable stopping routine.
  227. C
  228. C The Sparse Linear Algebra Package (SLAP) utilizes two matrix
  229. C data structures: 1) the SLAP Triad format or 2) the SLAP
  230. C Column format. The user can hand this routine either of the
  231. C of these data structures and SLAP will figure out which on
  232. C is being used and act accordingly.
  233. C
  234. C =================== S L A P Triad format ===================
  235. C This routine requires that the matrix A be stored in the
  236. C SLAP Triad format. In this format only the non-zeros are
  237. C stored. They may appear in *ANY* order. The user supplies
  238. C three arrays of length NELT, where NELT is the number of
  239. C non-zeros in the matrix: (IA(NELT), JA(NELT), A(NELT)). For
  240. C each non-zero the user puts the row and column index of that
  241. C matrix element in the IA and JA arrays. The value of the
  242. C non-zero matrix element is placed in the corresponding
  243. C location of the A array. This is an extremely easy data
  244. C structure to generate. On the other hand it is not too
  245. C efficient on vector computers for the iterative solution of
  246. C linear systems. Hence, SLAP changes this input data
  247. C structure to the SLAP Column format for the iteration (but
  248. C does not change it back).
  249. C
  250. C Here is an example of the SLAP Triad storage format for a
  251. C 5x5 Matrix. Recall that the entries may appear in any order.
  252. C
  253. C 5x5 Matrix SLAP Triad format for 5x5 matrix on left.
  254. C 1 2 3 4 5 6 7 8 9 10 11
  255. C |11 12 0 0 15| A: 51 12 11 33 15 53 55 22 35 44 21
  256. C |21 22 0 0 0| IA: 5 1 1 3 1 5 5 2 3 4 2
  257. C | 0 0 33 0 35| JA: 1 2 1 3 5 3 5 2 5 4 1
  258. C | 0 0 0 44 0|
  259. C |51 0 53 0 55|
  260. C
  261. C =================== S L A P Column format ==================
  262. C
  263. C This routine requires that the matrix A be stored in the
  264. C SLAP Column format. In this format the non-zeros are stored
  265. C counting down columns (except for the diagonal entry, which
  266. C must appear first in each "column") and are stored in the
  267. C double precision array A. In other words, for each column
  268. C in the matrix put the diagonal entry in A. Then put in the
  269. C other non-zero elements going down the column (except the
  270. C diagonal) in order. The IA array holds the row index for
  271. C each non-zero. The JA array holds the offsets into the IA,
  272. C A arrays for the beginning of each column. That is,
  273. C IA(JA(ICOL)), A(JA(ICOL)) points to the beginning of the
  274. C ICOL-th column in IA and A. IA(JA(ICOL+1)-1),
  275. C A(JA(ICOL+1)-1) points to the end of the ICOL-th column.
  276. C Note that we always have JA(N+1) = NELT+1, where N is the
  277. C number of columns in the matrix and NELT is the number of
  278. C non-zeros in the matrix.
  279. C
  280. C Here is an example of the SLAP Column storage format for a
  281. C 5x5 Matrix (in the A and IA arrays '|' denotes the end of a
  282. C column):
  283. C
  284. C 5x5 Matrix SLAP Column format for 5x5 matrix on left.
  285. C 1 2 3 4 5 6 7 8 9 10 11
  286. C |11 12 0 0 15| A: 11 21 51 | 22 12 | 33 53 | 44 | 55 15 35
  287. C |21 22 0 0 0| IA: 1 2 5 | 2 1 | 3 5 | 4 | 5 1 3
  288. C | 0 0 33 0 35| JA: 1 4 6 8 9 12
  289. C | 0 0 0 44 0|
  290. C |51 0 53 0 55|
  291. C
  292. C *Side Effects:
  293. C The SLAP Triad format (IA, JA, A) is modified internally to be
  294. C the SLAP Column format. See above.
  295. C
  296. C *Cautions:
  297. C This routine will attempt to write to the Fortran logical output
  298. C unit IUNIT, if IUNIT .ne. 0. Thus, the user must make sure that
  299. C this logical unit is attached to a file or terminal before calling
  300. C this routine with a non-zero value for IUNIT. This routine does
  301. C not check for the validity of a non-zero IUNIT unit number.
  302. C
  303. C***REFERENCES 1. Peter N. Brown and A. C. Hindmarsh, Reduced Storage
  304. C Matrix Methods in Stiff ODE Systems, Lawrence Liver-
  305. C more National Laboratory Report UCRL-95088, Rev. 1,
  306. C Livermore, California, June 1987.
  307. C***ROUTINES CALLED DCHKW, DGMRES, DS2Y, DSDI, DSDS, DSMV
  308. C***REVISION HISTORY (YYMMDD)
  309. C 890404 DATE WRITTEN
  310. C 890404 Previous REVISION DATE
  311. C 890915 Made changes requested at July 1989 CML Meeting. (MKS)
  312. C 890922 Numerous changes to prologue to make closer to SLATEC
  313. C standard. (FNF)
  314. C 890929 Numerous changes to reduce SP/DP differences. (FNF)
  315. C 910411 Prologue converted to Version 4.0 format. (BAB)
  316. C 920407 COMMON BLOCK renamed DSLBLK. (WRB)
  317. C 920511 Added complete declaration section. (WRB)
  318. C 920929 Corrected format of references. (FNF)
  319. C***END PROLOGUE DSDGMR
  320. C The following is for optimized compilation on LLNL/LTSS Crays.
  321. CLLL. OPTIMIZE
  322. C .. Parameters ..
  323. INTEGER LOCRB, LOCIB
  324. PARAMETER (LOCRB=1, LOCIB=11)
  325. C .. Scalar Arguments ..
  326. DOUBLE PRECISION ERR, TOL
  327. INTEGER IERR, ISYM, ITER, ITMAX, ITOL, IUNIT, LENIW, LENW, N,
  328. + NELT, NSAVE
  329. C .. Array Arguments ..
  330. DOUBLE PRECISION A(NELT), B(N), RWORK(LENW), X(N)
  331. INTEGER IA(NELT), IWORK(LENIW), JA(NELT)
  332. C .. Local Scalars ..
  333. INTEGER LOCDIN, LOCIGW, LOCIW, LOCRGW, LOCW, MYITOL
  334. C .. External Subroutines ..
  335. EXTERNAL DCHKW, DGMRES, DS2Y, DSDI, DSDS, DSMV
  336. C***FIRST EXECUTABLE STATEMENT DSDGMR
  337. C
  338. IERR = 0
  339. ERR = 0
  340. IF( NSAVE.LE.1 ) THEN
  341. IERR = 3
  342. RETURN
  343. ENDIF
  344. C
  345. C Change the SLAP input matrix IA, JA, A to SLAP-Column format.
  346. CALL DS2Y( N, NELT, IA, JA, A, ISYM )
  347. C
  348. C Set up the workspace. We assume MAXL=KMP=NSAVE.
  349. LOCIGW = LOCIB
  350. LOCIW = LOCIGW + 20
  351. C
  352. LOCDIN = LOCRB
  353. LOCRGW = LOCDIN + N
  354. LOCW = LOCRGW + 1+N*(NSAVE+6)+NSAVE*(NSAVE+3)
  355. C
  356. IWORK(4) = LOCDIN
  357. IWORK(9) = LOCIW
  358. IWORK(10) = LOCW
  359. C
  360. C Check the workspace allocations.
  361. CALL DCHKW( 'DSDGMR', LOCIW, LENIW, LOCW, LENW, IERR, ITER, ERR )
  362. IF( IERR.NE.0 ) RETURN
  363. C
  364. C Compute the inverse of the diagonal of the matrix.
  365. CALL DSDS(N, NELT, IA, JA, A, ISYM, RWORK(LOCDIN))
  366. C
  367. C Perform the Diagonally Scaled Generalized Minimum
  368. C Residual iteration algorithm. The following DGMRES
  369. C defaults are used MAXL = KMP = NSAVE, JSCAL = 0,
  370. C JPRE = -1, NRMAX = ITMAX/NSAVE
  371. IWORK(LOCIGW ) = NSAVE
  372. IWORK(LOCIGW+1) = NSAVE
  373. IWORK(LOCIGW+2) = 0
  374. IWORK(LOCIGW+3) = -1
  375. IWORK(LOCIGW+4) = ITMAX/NSAVE
  376. MYITOL = 0
  377. C
  378. CALL DGMRES( N, B, X, NELT, IA, JA, A, ISYM, DSMV, DSDI,
  379. $ MYITOL, TOL, ITMAX, ITER, ERR, IERR, IUNIT, RWORK, RWORK,
  380. $ RWORK(LOCRGW), LENW-LOCRGW, IWORK(LOCIGW), 20,
  381. $ RWORK, IWORK )
  382. C
  383. IF( ITER.GT.ITMAX ) IERR = 2
  384. RETURN
  385. C------------- LAST LINE OF DSDGMR FOLLOWS ----------------------------
  386. END