sslugm.f 20 KB

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