dgmres.f 26 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553
  1. *DECK DGMRES
  2. SUBROUTINE DGMRES (N, B, X, NELT, IA, JA, A, ISYM, MATVEC, MSOLVE,
  3. + ITOL, TOL, ITMAX, ITER, ERR, IERR, IUNIT, SB, SX, RGWK, LRGW,
  4. + IGWK, LIGW, RWORK, IWORK)
  5. C***BEGIN PROLOGUE DGMRES
  6. C***PURPOSE Preconditioned GMRES iterative sparse Ax=b solver.
  7. C This routine uses the generalized minimum residual
  8. C (GMRES) method with preconditioning to solve
  9. C non-symmetric linear systems of the form: Ax = b.
  10. C***LIBRARY SLATEC (SLAP)
  11. C***CATEGORY D2A4, D2B4
  12. C***TYPE DOUBLE PRECISION (SGMRES-S, DGMRES-D)
  13. C***KEYWORDS GENERALIZED MINIMUM RESIDUAL, ITERATIVE PRECONDITION,
  14. C NON-SYMMETRIC LINEAR SYSTEM, SLAP, SPARSE
  15. C***AUTHOR Brown, Peter, (LLNL), pnbrown@llnl.gov
  16. C Hindmarsh, Alan, (LLNL), alanh@llnl.gov
  17. C Seager, Mark K., (LLNL), seager@llnl.gov
  18. C Lawrence Livermore National Laboratory
  19. C PO Box 808, L-60
  20. C Livermore, CA 94550 (510) 423-3141
  21. C***DESCRIPTION
  22. C
  23. C *Usage:
  24. C INTEGER N, NELT, IA(NELT), JA(NELT), ISYM, ITOL, ITMAX
  25. C INTEGER ITER, IERR, IUNIT, LRGW, IGWK(LIGW), LIGW
  26. C INTEGER IWORK(USER DEFINED)
  27. C DOUBLE PRECISION B(N), X(N), A(NELT), TOL, ERR, SB(N), SX(N)
  28. C DOUBLE PRECISION RGWK(LRGW), RWORK(USER DEFINED)
  29. C EXTERNAL MATVEC, MSOLVE
  30. C
  31. C CALL DGMRES(N, B, X, NELT, IA, JA, A, ISYM, MATVEC, MSOLVE,
  32. C $ ITOL, TOL, ITMAX, ITER, ERR, IERR, IUNIT, SB, SX,
  33. C $ RGWK, LRGW, IGWK, LIGW, RWORK, IWORK)
  34. C
  35. C *Arguments:
  36. C N :IN Integer.
  37. C Order of the Matrix.
  38. C B :IN Double Precision B(N).
  39. C Right-hand side vector.
  40. C X :INOUT Double Precision X(N).
  41. C On input X is your initial guess for the solution vector.
  42. C On output X is the final approximate solution.
  43. C NELT :IN Integer.
  44. C Number of Non-Zeros stored in A.
  45. C IA :IN Integer IA(NELT).
  46. C JA :IN Integer JA(NELT).
  47. C A :IN Double Precision A(NELT).
  48. C These arrays contain the matrix data structure for A.
  49. C It could take any form. See "Description", below,
  50. C for more details.
  51. C ISYM :IN Integer.
  52. C Flag to indicate symmetric storage format.
  53. C If ISYM=0, all non-zero entries of the matrix are stored.
  54. C If ISYM=1, the matrix is symmetric, and only the upper
  55. C or lower triangle of the matrix is stored.
  56. C MATVEC :EXT External.
  57. C Name of a routine which performs the matrix vector multiply
  58. C Y = A*X given A and X. The name of the MATVEC routine must
  59. C be declared external in the calling program. The calling
  60. C sequence to MATVEC is:
  61. C CALL MATVEC(N, X, Y, NELT, IA, JA, A, ISYM)
  62. C where N is the number of unknowns, Y is the product A*X
  63. C upon return, X is an input vector, and NELT is the number of
  64. C non-zeros in the SLAP IA, JA, A storage for the matrix A.
  65. C ISYM is a flag which, if non-zero, denotes that A is
  66. C symmetric and only the lower or upper triangle is stored.
  67. C MSOLVE :EXT External.
  68. C Name of the routine which solves a linear system Mz = r for
  69. C z given r with the preconditioning matrix M (M is supplied via
  70. C RWORK and IWORK arrays. The name of the MSOLVE routine must
  71. C be declared external in the calling program. The calling
  72. C sequence to MSOLVE is:
  73. C CALL MSOLVE(N, R, Z, NELT, IA, JA, A, ISYM, RWORK, IWORK)
  74. C Where N is the number of unknowns, R is the right-hand side
  75. C vector and Z is the solution upon return. NELT, IA, JA, A and
  76. C ISYM are defined as above. RWORK is a double precision array
  77. C that can be used to pass necessary preconditioning information
  78. C and/or workspace to MSOLVE. IWORK is an integer work array
  79. C for the same purpose as RWORK.
  80. C ITOL :IN Integer.
  81. C Flag to indicate the type of convergence criterion used.
  82. C ITOL=0 Means the iteration stops when the test described
  83. C below on the residual RL is satisfied. This is
  84. C the "Natural Stopping Criteria" for this routine.
  85. C Other values of ITOL cause extra, otherwise
  86. C unnecessary, computation per iteration and are
  87. C therefore much less efficient. See ISDGMR (the
  88. C stop test routine) for more information.
  89. C ITOL=1 Means the iteration stops when the first test
  90. C described below on the residual RL is satisfied,
  91. C and there is either right or no preconditioning
  92. C being used.
  93. C ITOL=2 Implies that the user is using left
  94. C preconditioning, and the second stopping criterion
  95. C below is used.
  96. C ITOL=3 Means the iteration stops when the third test
  97. C described below on Minv*Residual is satisfied, and
  98. C there is either left or no preconditioning being
  99. C used.
  100. C ITOL=11 is often useful for checking and comparing
  101. C different routines. For this case, the user must
  102. C supply the "exact" solution or a very accurate
  103. C approximation (one with an error much less than
  104. C TOL) through a common block,
  105. C COMMON /DSLBLK/ SOLN( )
  106. C If ITOL=11, iteration stops when the 2-norm of the
  107. C difference between the iterative approximation and
  108. C the user-supplied solution divided by the 2-norm
  109. C of the user-supplied solution is less than TOL.
  110. C Note that this requires the user to set up the
  111. C "COMMON /DSLBLK/ SOLN(LENGTH)" in the calling
  112. C routine. The routine with this declaration should
  113. C be loaded before the stop test so that the correct
  114. C length is used by the loader. This procedure is
  115. C not standard Fortran and may not work correctly on
  116. C your system (although it has worked on every
  117. C system the authors have tried). If ITOL is not 11
  118. C then this common block is indeed standard Fortran.
  119. C TOL :INOUT Double Precision.
  120. C Convergence criterion, as described below. If TOL is set
  121. C to zero on input, then a default value of 500*(the smallest
  122. C positive magnitude, machine epsilon) is used.
  123. C ITMAX :DUMMY Integer.
  124. C Maximum number of iterations in most SLAP routines. In
  125. C this routine this does not make sense. The maximum number
  126. C of iterations here is given by ITMAX = MAXL*(NRMAX+1).
  127. C See IGWK for definitions of MAXL and NRMAX.
  128. C ITER :OUT Integer.
  129. C Number of iterations required to reach convergence, or
  130. C ITMAX if convergence criterion could not be achieved in
  131. C ITMAX iterations.
  132. C ERR :OUT Double Precision.
  133. C Error estimate of error in final approximate solution, as
  134. C defined by ITOL. Letting norm() denote the Euclidean
  135. C norm, ERR is defined as follows..
  136. C
  137. C If ITOL=0, then ERR = norm(SB*(B-A*X(L)))/norm(SB*B),
  138. C for right or no preconditioning, and
  139. C ERR = norm(SB*(M-inverse)*(B-A*X(L)))/
  140. C norm(SB*(M-inverse)*B),
  141. C for left preconditioning.
  142. C If ITOL=1, then ERR = norm(SB*(B-A*X(L)))/norm(SB*B),
  143. C since right or no preconditioning
  144. C being used.
  145. C If ITOL=2, then ERR = norm(SB*(M-inverse)*(B-A*X(L)))/
  146. C norm(SB*(M-inverse)*B),
  147. C since left preconditioning is being
  148. C used.
  149. C If ITOL=3, then ERR = Max |(Minv*(B-A*X(L)))(i)/x(i)|
  150. C i=1,n
  151. C If ITOL=11, then ERR = norm(SB*(X(L)-SOLN))/norm(SB*SOLN).
  152. C IERR :OUT Integer.
  153. C Return error flag.
  154. C IERR = 0 => All went well.
  155. C IERR = 1 => Insufficient storage allocated for
  156. C RGWK or IGWK.
  157. C IERR = 2 => Routine DGMRES failed to reduce the norm
  158. C of the current residual on its last call,
  159. C and so the iteration has stalled. In
  160. C this case, X equals the last computed
  161. C approximation. The user must either
  162. C increase MAXL, or choose a different
  163. C initial guess.
  164. C IERR =-1 => Insufficient length for RGWK array.
  165. C IGWK(6) contains the required minimum
  166. C length of the RGWK array.
  167. C IERR =-2 => Illegal value of ITOL, or ITOL and JPRE
  168. C values are inconsistent.
  169. C For IERR <= 2, RGWK(1) = RHOL, which is the norm on the
  170. C left-hand-side of the relevant stopping test defined
  171. C below associated with the residual for the current
  172. C approximation X(L).
  173. C IUNIT :IN Integer.
  174. C Unit number on which to write the error at each iteration,
  175. C if this is desired for monitoring convergence. If unit
  176. C number is 0, no writing will occur.
  177. C SB :IN Double Precision SB(N).
  178. C Array of length N containing scale factors for the right
  179. C hand side vector B. If JSCAL.eq.0 (see below), SB need
  180. C not be supplied.
  181. C SX :IN Double Precision SX(N).
  182. C Array of length N containing scale factors for the solution
  183. C vector X. If JSCAL.eq.0 (see below), SX need not be
  184. C supplied. SB and SX can be the same array in the calling
  185. C program if desired.
  186. C RGWK :INOUT Double Precision RGWK(LRGW).
  187. C Double Precision array used for workspace by DGMRES.
  188. C On return, RGWK(1) = RHOL. See IERR for definition of RHOL.
  189. C LRGW :IN Integer.
  190. C Length of the double precision workspace, RGWK.
  191. C LRGW >= 1 + N*(MAXL+6) + MAXL*(MAXL+3).
  192. C See below for definition of MAXL.
  193. C For the default values, RGWK has size at least 131 + 16*N.
  194. C IGWK :INOUT Integer IGWK(LIGW).
  195. C The following IGWK parameters should be set by the user
  196. C before calling this routine.
  197. C IGWK(1) = MAXL. Maximum dimension of Krylov subspace in
  198. C which X - X0 is to be found (where, X0 is the initial
  199. C guess). The default value of MAXL is 10.
  200. C IGWK(2) = KMP. Maximum number of previous Krylov basis
  201. C vectors to which each new basis vector is made orthogonal.
  202. C The default value of KMP is MAXL.
  203. C IGWK(3) = JSCAL. Flag indicating whether the scaling
  204. C arrays SB and SX are to be used.
  205. C JSCAL = 0 => SB and SX are not used and the algorithm
  206. C will perform as if all SB(I) = 1 and SX(I) = 1.
  207. C JSCAL = 1 => Only SX is used, and the algorithm
  208. C performs as if all SB(I) = 1.
  209. C JSCAL = 2 => Only SB is used, and the algorithm
  210. C performs as if all SX(I) = 1.
  211. C JSCAL = 3 => Both SB and SX are used.
  212. C IGWK(4) = JPRE. Flag indicating whether preconditioning
  213. C is being used.
  214. C JPRE = 0 => There is no preconditioning.
  215. C JPRE > 0 => There is preconditioning on the right
  216. C only, and the solver will call routine MSOLVE.
  217. C JPRE < 0 => There is preconditioning on the left
  218. C only, and the solver will call routine MSOLVE.
  219. C IGWK(5) = NRMAX. Maximum number of restarts of the
  220. C Krylov iteration. The default value of NRMAX = 10.
  221. C if IWORK(5) = -1, then no restarts are performed (in
  222. C this case, NRMAX is set to zero internally).
  223. C The following IWORK parameters are diagnostic information
  224. C made available to the user after this routine completes.
  225. C IGWK(6) = MLWK. Required minimum length of RGWK array.
  226. C IGWK(7) = NMS. The total number of calls to MSOLVE.
  227. C LIGW :IN Integer.
  228. C Length of the integer workspace, IGWK. LIGW >= 20.
  229. C RWORK :WORK Double Precision RWORK(USER DEFINED).
  230. C Double Precision array that can be used for workspace in
  231. C MSOLVE.
  232. C IWORK :WORK Integer IWORK(USER DEFINED).
  233. C Integer array that can be used for workspace in MSOLVE.
  234. C
  235. C *Description:
  236. C DGMRES solves a linear system A*X = B rewritten in the form:
  237. C
  238. C (SB*A*(M-inverse)*(SX-inverse))*(SX*M*X) = SB*B,
  239. C
  240. C with right preconditioning, or
  241. C
  242. C (SB*(M-inverse)*A*(SX-inverse))*(SX*X) = SB*(M-inverse)*B,
  243. C
  244. C with left preconditioning, where A is an N-by-N double precision
  245. C matrix, X and B are N-vectors, SB and SX are diagonal scaling
  246. C matrices, and M is a preconditioning matrix. It uses
  247. C preconditioned Krylov subpace methods based on the
  248. C generalized minimum residual method (GMRES). This routine
  249. C optionally performs either the full orthogonalization
  250. C version of the GMRES algorithm or an incomplete variant of
  251. C it. Both versions use restarting of the linear iteration by
  252. C default, although the user can disable this feature.
  253. C
  254. C The GMRES algorithm generates a sequence of approximations
  255. C X(L) to the true solution of the above linear system. The
  256. C convergence criteria for stopping the iteration is based on
  257. C the size of the scaled norm of the residual R(L) = B -
  258. C A*X(L). The actual stopping test is either:
  259. C
  260. C norm(SB*(B-A*X(L))) .le. TOL*norm(SB*B),
  261. C
  262. C for right preconditioning, or
  263. C
  264. C norm(SB*(M-inverse)*(B-A*X(L))) .le.
  265. C TOL*norm(SB*(M-inverse)*B),
  266. C
  267. C for left preconditioning, where norm() denotes the Euclidean
  268. C norm, and TOL is a positive scalar less than one input by
  269. C the user. If TOL equals zero when DGMRES is called, then a
  270. C default value of 500*(the smallest positive magnitude,
  271. C machine epsilon) is used. If the scaling arrays SB and SX
  272. C are used, then ideally they should be chosen so that the
  273. C vectors SX*X(or SX*M*X) and SB*B have all their components
  274. C approximately equal to one in magnitude. If one wants to
  275. C use the same scaling in X and B, then SB and SX can be the
  276. C same array in the calling program.
  277. C
  278. C The following is a list of the other routines and their
  279. C functions used by DGMRES:
  280. C DPIGMR Contains the main iteration loop for GMRES.
  281. C DORTH Orthogonalizes a new vector against older basis vectors.
  282. C DHEQR Computes a QR decomposition of a Hessenberg matrix.
  283. C DHELS Solves a Hessenberg least-squares system, using QR
  284. C factors.
  285. C DRLCAL Computes the scaled residual RL.
  286. C DXLCAL Computes the solution XL.
  287. C ISDGMR User-replaceable stopping routine.
  288. C
  289. C This routine does not care what matrix data structure is
  290. C used for A and M. It simply calls the MATVEC and MSOLVE
  291. C routines, with the arguments as described above. The user
  292. C could write any type of structure and the appropriate MATVEC
  293. C and MSOLVE routines. It is assumed that A is stored in the
  294. C IA, JA, A arrays in some fashion and that M (or INV(M)) is
  295. C stored in IWORK and RWORK in some fashion. The SLAP
  296. C routines DSDCG and DSICCG are examples of this procedure.
  297. C
  298. C Two examples of matrix data structures are the: 1) SLAP
  299. C Triad format and 2) SLAP Column format.
  300. C
  301. C =================== S L A P Triad format ===================
  302. C This routine requires that the matrix A be stored in the
  303. C SLAP Triad format. In this format only the non-zeros are
  304. C stored. They may appear in *ANY* order. The user supplies
  305. C three arrays of length NELT, where NELT is the number of
  306. C non-zeros in the matrix: (IA(NELT), JA(NELT), A(NELT)). For
  307. C each non-zero the user puts the row and column index of that
  308. C matrix element in the IA and JA arrays. The value of the
  309. C non-zero matrix element is placed in the corresponding
  310. C location of the A array. This is an extremely easy data
  311. C structure to generate. On the other hand it is not too
  312. C efficient on vector computers for the iterative solution of
  313. C linear systems. Hence, SLAP changes this input data
  314. C structure to the SLAP Column format for the iteration (but
  315. C does not change it back).
  316. C
  317. C Here is an example of the SLAP Triad storage format for a
  318. C 5x5 Matrix. Recall that the entries may appear in any order.
  319. C
  320. C 5x5 Matrix SLAP Triad format for 5x5 matrix on left.
  321. C 1 2 3 4 5 6 7 8 9 10 11
  322. C |11 12 0 0 15| A: 51 12 11 33 15 53 55 22 35 44 21
  323. C |21 22 0 0 0| IA: 5 1 1 3 1 5 5 2 3 4 2
  324. C | 0 0 33 0 35| JA: 1 2 1 3 5 3 5 2 5 4 1
  325. C | 0 0 0 44 0|
  326. C |51 0 53 0 55|
  327. C
  328. C =================== S L A P Column format ==================
  329. C
  330. C This routine requires that the matrix A be stored in the
  331. C SLAP Column format. In this format the non-zeros are stored
  332. C counting down columns (except for the diagonal entry, which
  333. C must appear first in each "column") and are stored in the
  334. C double precision array A. In other words, for each column
  335. C in the matrix put the diagonal entry in A. Then put in the
  336. C other non-zero elements going down the column (except the
  337. C diagonal) in order. The IA array holds the row index for
  338. C each non-zero. The JA array holds the offsets into the IA,
  339. C A arrays for the beginning of each column. That is,
  340. C IA(JA(ICOL)), A(JA(ICOL)) points to the beginning of the
  341. C ICOL-th column in IA and A. IA(JA(ICOL+1)-1),
  342. C A(JA(ICOL+1)-1) points to the end of the ICOL-th column.
  343. C Note that we always have JA(N+1) = NELT+1, where N is the
  344. C number of columns in the matrix and NELT is the number of
  345. C non-zeros in the matrix.
  346. C
  347. C Here is an example of the SLAP Column storage format for a
  348. C 5x5 Matrix (in the A and IA arrays '|' denotes the end of a
  349. C column):
  350. C
  351. C 5x5 Matrix SLAP Column format for 5x5 matrix on left.
  352. C 1 2 3 4 5 6 7 8 9 10 11
  353. C |11 12 0 0 15| A: 11 21 51 | 22 12 | 33 53 | 44 | 55 15 35
  354. C |21 22 0 0 0| IA: 1 2 5 | 2 1 | 3 5 | 4 | 5 1 3
  355. C | 0 0 33 0 35| JA: 1 4 6 8 9 12
  356. C | 0 0 0 44 0|
  357. C |51 0 53 0 55|
  358. C
  359. C *Cautions:
  360. C This routine will attempt to write to the Fortran logical output
  361. C unit IUNIT, if IUNIT .ne. 0. Thus, the user must make sure that
  362. C this logical unit is attached to a file or terminal before calling
  363. C this routine with a non-zero value for IUNIT. This routine does
  364. C not check for the validity of a non-zero IUNIT unit number.
  365. C
  366. C***REFERENCES 1. Peter N. Brown and A. C. Hindmarsh, Reduced Storage
  367. C Matrix Methods in Stiff ODE Systems, Lawrence Liver-
  368. C more National Laboratory Report UCRL-95088, Rev. 1,
  369. C Livermore, California, June 1987.
  370. C 2. Mark K. Seager, A SLAP for the Masses, in
  371. C G. F. Carey, Ed., Parallel Supercomputing: Methods,
  372. C Algorithms and Applications, Wiley, 1989, pp.135-155.
  373. C***ROUTINES CALLED D1MACH, DCOPY, DNRM2, DPIGMR
  374. C***REVISION HISTORY (YYMMDD)
  375. C 890404 DATE WRITTEN
  376. C 890404 Previous REVISION DATE
  377. C 890915 Made changes requested at July 1989 CML Meeting. (MKS)
  378. C 890922 Numerous changes to prologue to make closer to SLATEC
  379. C standard. (FNF)
  380. C 890929 Numerous changes to reduce SP/DP differences. (FNF)
  381. C 891004 Added new reference.
  382. C 910411 Prologue converted to Version 4.0 format. (BAB)
  383. C 910506 Corrected errors in C***ROUTINES CALLED list. (FNF)
  384. C 920407 COMMON BLOCK renamed DSLBLK. (WRB)
  385. C 920511 Added complete declaration section. (WRB)
  386. C 920929 Corrected format of references. (FNF)
  387. C 921019 Changed 500.0 to 500 to reduce SP/DP differences. (FNF)
  388. C 921026 Added check for valid value of ITOL. (FNF)
  389. C***END PROLOGUE DGMRES
  390. C The following is for optimized compilation on LLNL/LTSS Crays.
  391. CLLL. OPTIMIZE
  392. C .. Scalar Arguments ..
  393. DOUBLE PRECISION ERR, TOL
  394. INTEGER IERR, ISYM, ITER, ITMAX, ITOL, IUNIT, LIGW, LRGW, N, NELT
  395. C .. Array Arguments ..
  396. DOUBLE PRECISION A(NELT), B(N), RGWK(LRGW), RWORK(*), SB(N),
  397. + SX(N), X(N)
  398. INTEGER IA(NELT), IGWK(LIGW), IWORK(*), JA(NELT)
  399. C .. Subroutine Arguments ..
  400. EXTERNAL MATVEC, MSOLVE
  401. C .. Local Scalars ..
  402. DOUBLE PRECISION BNRM, RHOL, SUM
  403. INTEGER I, IFLAG, JPRE, JSCAL, KMP, LDL, LGMR, LHES, LQ, LR, LV,
  404. + LW, LXL, LZ, LZM1, MAXL, MAXLP1, NMS, NMSL, NRMAX, NRSTS
  405. C .. External Functions ..
  406. DOUBLE PRECISION D1MACH, DNRM2
  407. EXTERNAL D1MACH, DNRM2
  408. C .. External Subroutines ..
  409. EXTERNAL DCOPY, DPIGMR
  410. C .. Intrinsic Functions ..
  411. INTRINSIC SQRT
  412. C***FIRST EXECUTABLE STATEMENT DGMRES
  413. IERR = 0
  414. C ------------------------------------------------------------------
  415. C Load method parameters with user values or defaults.
  416. C ------------------------------------------------------------------
  417. MAXL = IGWK(1)
  418. IF (MAXL .EQ. 0) MAXL = 10
  419. IF (MAXL .GT. N) MAXL = N
  420. KMP = IGWK(2)
  421. IF (KMP .EQ. 0) KMP = MAXL
  422. IF (KMP .GT. MAXL) KMP = MAXL
  423. JSCAL = IGWK(3)
  424. JPRE = IGWK(4)
  425. C Check for valid value of ITOL.
  426. IF( (ITOL.LT.0) .OR. ((ITOL.GT.3).AND.(ITOL.NE.11)) ) GOTO 650
  427. C Check for consistent values of ITOL and JPRE.
  428. IF( ITOL.EQ.1 .AND. JPRE.LT.0 ) GOTO 650
  429. IF( ITOL.EQ.2 .AND. JPRE.GE.0 ) GOTO 650
  430. NRMAX = IGWK(5)
  431. IF( NRMAX.EQ.0 ) NRMAX = 10
  432. C If NRMAX .eq. -1, then set NRMAX = 0 to turn off restarting.
  433. IF( NRMAX.EQ.-1 ) NRMAX = 0
  434. C If input value of TOL is zero, set it to its default value.
  435. IF( TOL.EQ.0.0D0 ) TOL = 500*D1MACH(3)
  436. C
  437. C Initialize counters.
  438. ITER = 0
  439. NMS = 0
  440. NRSTS = 0
  441. C ------------------------------------------------------------------
  442. C Form work array segment pointers.
  443. C ------------------------------------------------------------------
  444. MAXLP1 = MAXL + 1
  445. LV = 1
  446. LR = LV + N*MAXLP1
  447. LHES = LR + N + 1
  448. LQ = LHES + MAXL*MAXLP1
  449. LDL = LQ + 2*MAXL
  450. LW = LDL + N
  451. LXL = LW + N
  452. LZ = LXL + N
  453. C
  454. C Load IGWK(6) with required minimum length of the RGWK array.
  455. IGWK(6) = LZ + N - 1
  456. IF( LZ+N-1.GT.LRGW ) GOTO 640
  457. C ------------------------------------------------------------------
  458. C Calculate scaled-preconditioned norm of RHS vector b.
  459. C ------------------------------------------------------------------
  460. IF (JPRE .LT. 0) THEN
  461. CALL MSOLVE(N, B, RGWK(LR), NELT, IA, JA, A, ISYM,
  462. $ RWORK, IWORK)
  463. NMS = NMS + 1
  464. ELSE
  465. CALL DCOPY(N, B, 1, RGWK(LR), 1)
  466. ENDIF
  467. IF( JSCAL.EQ.2 .OR. JSCAL.EQ.3 ) THEN
  468. SUM = 0
  469. DO 10 I = 1,N
  470. SUM = SUM + (RGWK(LR-1+I)*SB(I))**2
  471. 10 CONTINUE
  472. BNRM = SQRT(SUM)
  473. ELSE
  474. BNRM = DNRM2(N,RGWK(LR),1)
  475. ENDIF
  476. C ------------------------------------------------------------------
  477. C Calculate initial residual.
  478. C ------------------------------------------------------------------
  479. CALL MATVEC(N, X, RGWK(LR), NELT, IA, JA, A, ISYM)
  480. DO 50 I = 1,N
  481. RGWK(LR-1+I) = B(I) - RGWK(LR-1+I)
  482. 50 CONTINUE
  483. C ------------------------------------------------------------------
  484. C If performing restarting, then load the residual into the
  485. C correct location in the RGWK array.
  486. C ------------------------------------------------------------------
  487. 100 CONTINUE
  488. IF( NRSTS.GT.NRMAX ) GOTO 610
  489. IF( NRSTS.GT.0 ) THEN
  490. C Copy the current residual to a different location in the RGWK
  491. C array.
  492. CALL DCOPY(N, RGWK(LDL), 1, RGWK(LR), 1)
  493. ENDIF
  494. C ------------------------------------------------------------------
  495. C Use the DPIGMR algorithm to solve the linear system A*Z = R.
  496. C ------------------------------------------------------------------
  497. CALL DPIGMR(N, RGWK(LR), SB, SX, JSCAL, MAXL, MAXLP1, KMP,
  498. $ NRSTS, JPRE, MATVEC, MSOLVE, NMSL, RGWK(LZ), RGWK(LV),
  499. $ RGWK(LHES), RGWK(LQ), LGMR, RWORK, IWORK, RGWK(LW),
  500. $ RGWK(LDL), RHOL, NRMAX, B, BNRM, X, RGWK(LXL), ITOL,
  501. $ TOL, NELT, IA, JA, A, ISYM, IUNIT, IFLAG, ERR)
  502. ITER = ITER + LGMR
  503. NMS = NMS + NMSL
  504. C
  505. C Increment X by the current approximate solution Z of A*Z = R.
  506. C
  507. LZM1 = LZ - 1
  508. DO 110 I = 1,N
  509. X(I) = X(I) + RGWK(LZM1+I)
  510. 110 CONTINUE
  511. IF( IFLAG.EQ.0 ) GOTO 600
  512. IF( IFLAG.EQ.1 ) THEN
  513. NRSTS = NRSTS + 1
  514. GOTO 100
  515. ENDIF
  516. IF( IFLAG.EQ.2 ) GOTO 620
  517. C ------------------------------------------------------------------
  518. C All returns are made through this section.
  519. C ------------------------------------------------------------------
  520. C The iteration has converged.
  521. C
  522. 600 CONTINUE
  523. IGWK(7) = NMS
  524. RGWK(1) = RHOL
  525. IERR = 0
  526. RETURN
  527. C
  528. C Max number((NRMAX+1)*MAXL) of linear iterations performed.
  529. 610 CONTINUE
  530. IGWK(7) = NMS
  531. RGWK(1) = RHOL
  532. IERR = 1
  533. RETURN
  534. C
  535. C GMRES failed to reduce last residual in MAXL iterations.
  536. C The iteration has stalled.
  537. 620 CONTINUE
  538. IGWK(7) = NMS
  539. RGWK(1) = RHOL
  540. IERR = 2
  541. RETURN
  542. C Error return. Insufficient length for RGWK array.
  543. 640 CONTINUE
  544. ERR = TOL
  545. IERR = -1
  546. RETURN
  547. C Error return. Inconsistent ITOL and JPRE values.
  548. 650 CONTINUE
  549. ERR = TOL
  550. IERR = -2
  551. RETURN
  552. C------------- LAST LINE OF DGMRES FOLLOWS ----------------------------
  553. END