2
SUBROUTINE DSLUGM (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 DSLUGM
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 DOUBLE 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
20
C Livermore, CA 94550 (510) 423-3141
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 DOUBLE PRECISION B(N), X(N), A(NELT), TOL, ERR, RWORK(LENW)
28
C CALL DSLUGM(N, B, X, NELT, IA, JA, A, ISYM, NSAVE,
29
C $ ITOL, TOL, ITMAX, ITER, ERR, IERR, IUNIT,
30
C $ RWORK, LENW, IWORK, LENIW)
34
C Order of the Matrix.
35
C B :IN Double Precision B(N).
36
C Right-hand side vector.
37
C X :INOUT Double Precision X(N).
38
C On input X is your initial guess for solution vector.
39
C On output X is the final approximate solution.
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 Double Precision 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.
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.
55
C Number of direction vectors to save and orthogonalize against.
56
C Must be greater than 1.
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 ISDGMR (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
70
C ITOL=2 Implies that the user is using left
71
C preconditioning, and the second stopping criterion
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
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 /DSLBLK/ 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 /DSLBLK/ 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 Double Precision.
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.
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 DGMRES for a full and frightfully interesting discussion of
107
C Number of iterations required to reach convergence, or
108
C ITMAX+1 if convergence criterion could not be achieved in
110
C ERR :OUT Double Precision.
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
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
126
C If ITOL=3, then ERR = Max |(Minv*(B-A*X(L)))(i)/x(i)|
128
C If ITOL=11, then ERR = norm(SB*(X(L)-SOLN))/norm(SB*SOLN).
131
C IERR = 0 => All went well.
132
C IERR = 1 => Insufficient storage allocated for
134
C IERR = 2 => Routine DPIGMR 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
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).
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 Double Precision RWORK(LENW).
154
C Double Precision array of size LENW.
156
C Length of the double precision 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
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 Double Precision workspace actually used.
171
C Length of the integer workspace, IWORK.
172
C LENIW >= NL+NU+4*N+32.
175
C DSLUGM solves a linear system A*X = B rewritten in the form:
177
C (SB*A*(M-inverse)*(SX-inverse))*(SX*M*X) = SB*B,
179
C with right preconditioning, or
181
C (SB*(M-inverse)*A*(SX-inverse))*(SX*X) = SB*(M-inverse)*B,
183
C with left preconditioning, where A is an n-by-n double precision
184
C matrix, 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 DGMRES for the solution of the linear system. DGMRES
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.
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:
203
C norm(SB*(B-A*X(L))) .le. TOL*norm(SB*B),
205
C for right preconditioning, or
207
C norm(SB*(M-inverse)*(B-A*X(L))) .le.
208
C TOL*norm(SB*(M-inverse)*B),
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 DSLUGM 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.
221
C The following is a list of the other routines and their
222
C functions used by GMRES:
223
C DGMRES Contains the matrix structure independent driver
225
C DPIGMR Contains the main iteration loop for GMRES.
226
C DORTH Orthogonalizes a new vector against older basis vectors.
227
C DHEQR Computes a QR decomposition of a Hessenberg matrix.
228
C DHELS Solves a Hessenberg least-squares system, using QR
230
C RLCALC Computes the scaled residual RL.
231
C XLCALC Computes the solution XL.
232
C ISDGMR User-replaceable stopping routine.
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.
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).
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.
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
267
C =================== S L A P Column format ==================
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 double precision array A. In other words, for each column
274
C in the matrix put the diagonal entry in A. Then put in the
275
C other non-zero elements going down the column (except the
276
C diagonal) in order. The IA array holds the row index for
277
C each non-zero. The JA array holds the offsets into the IA,
278
C A arrays for the beginning of each column. That is,
279
C IA(JA(ICOL)), A(JA(ICOL)) points to the beginning of the
280
C ICOL-th column in IA and A. IA(JA(ICOL+1)-1),
281
C A(JA(ICOL+1)-1) points to the end of the ICOL-th column.
282
C Note that we always have JA(N+1) = NELT+1, where N is the
283
C number of columns in the matrix and NELT is the number of
284
C non-zeros in the matrix.
286
C Here is an example of the SLAP Column storage format for a
287
C 5x5 Matrix (in the A and IA arrays '|' denotes the end of a
290
C 5x5 Matrix SLAP Column format for 5x5 matrix on left.
291
C 1 2 3 4 5 6 7 8 9 10 11
292
C |11 12 0 0 15| A: 11 21 51 | 22 12 | 33 53 | 44 | 55 15 35
293
C |21 22 0 0 0| IA: 1 2 5 | 2 1 | 3 5 | 4 | 5 1 3
294
C | 0 0 33 0 35| JA: 1 4 6 8 9 12
299
C The SLAP Triad format (IA, JA, A) is modified internally to be
300
C the SLAP Column format. See above.
303
C This routine will attempt to write to the Fortran logical output
304
C unit IUNIT, if IUNIT .ne. 0. Thus, the user must make sure that
305
C this logical unit is attached to a file or terminal before calling
306
C this routine with a non-zero value for IUNIT. This routine does
307
C not check for the validity of a non-zero IUNIT unit number.
309
C***REFERENCES 1. Peter N. Brown and A. C. Hindmarsh, Reduced Storage
310
C Matrix Methods in Stiff ODE Systems, Lawrence Liver-
311
C more National Laboratory Report UCRL-95088, Rev. 1,
312
C Livermore, California, June 1987.
313
C***ROUTINES CALLED DCHKW, DGMRES, DS2Y, DSILUS, DSLUI, DSMV
314
C***REVISION HISTORY (YYMMDD)
315
C 890404 DATE WRITTEN
316
C 890404 Previous REVISION DATE
317
C 890915 Made changes requested at July 1989 CML Meeting. (MKS)
318
C 890922 Numerous changes to prologue to make closer to SLATEC
320
C 890929 Numerous changes to reduce SP/DP differences. (FNF)
321
C 910411 Prologue converted to Version 4.0 format. (BAB)
322
C 920407 COMMON BLOCK renamed DSLBLK. (WRB)
323
C 920511 Added complete declaration section. (WRB)
324
C 920929 Corrected format of references. (FNF)
325
C 921019 Corrected NEL to NL. (FNF)
326
C***END PROLOGUE DSLUGM
327
C The following is for optimized compilation on LLNL/LTSS Crays.
331
PARAMETER (LOCRB=1, LOCIB=11)
332
C .. Scalar Arguments ..
333
DOUBLE PRECISION ERR, TOL
334
INTEGER IERR, ISYM, ITER, ITMAX, ITOL, IUNIT, LENIW, LENW, N,
336
C .. Array Arguments ..
337
DOUBLE PRECISION A(NELT), B(N), RWORK(LENW), X(N)
338
INTEGER IA(NELT), IWORK(LENIW), JA(NELT)
339
C .. Local Scalars ..
340
INTEGER ICOL, J, JBGN, JEND, LOCDIN, LOCIGW, LOCIL, LOCIU, LOCIW,
341
+ LOCJL, LOCJU, LOCL, LOCNC, LOCNR, LOCRGW, LOCU, LOCW,
343
C .. External Subroutines ..
344
EXTERNAL DCHKW, DGMRES, DS2Y, DSILUS, DSLUI, DSMV
345
C***FIRST EXECUTABLE STATEMENT DSLUGM
349
IF( NSAVE.LE.1 ) THEN
354
C Change the SLAP input matrix IA, JA, A to SLAP-Column format.
355
CALL DS2Y( N, NELT, IA, JA, A, ISYM )
357
C Count number of Non-Zero elements preconditioner ILU matrix.
358
C Then set up the work arrays. We assume MAXL=KMP=NSAVE.
362
C Don't count diagonal.
365
IF( JBGN.LE.JEND ) THEN
368
IF( IA(J).GT.ICOL ) THEN
370
IF( ISYM.NE.0 ) NU = NU + 1
391
LOCW = LOCRGW + 1+N*(NSAVE+6)+NSAVE*(NSAVE+3)
393
C Check the workspace allocations.
394
CALL DCHKW( 'DSLUGM', LOCIW, LENIW, LOCW, LENW, IERR, ITER, ERR )
395
IF( IERR.NE.0 ) RETURN
407
C Compute the Incomplete LU decomposition.
408
CALL DSILUS( N, NELT, IA, JA, A, ISYM, NL, IWORK(LOCIL),
409
$ IWORK(LOCJL), RWORK(LOCL), RWORK(LOCDIN), NU, IWORK(LOCIU),
410
$ IWORK(LOCJU), RWORK(LOCU), IWORK(LOCNR), IWORK(LOCNC) )
412
C Perform the Incomplete LU Preconditioned Generalized Minimum
413
C Residual iteration algorithm. The following DGMRES
414
C defaults are used MAXL = KMP = NSAVE, JSCAL = 0,
415
C JPRE = -1, NRMAX = ITMAX/NSAVE
416
IWORK(LOCIGW ) = NSAVE
417
IWORK(LOCIGW+1) = NSAVE
420
IWORK(LOCIGW+4) = ITMAX/NSAVE
423
CALL DGMRES( N, B, X, NELT, IA, JA, A, ISYM, DSMV, DSLUI,
424
$ MYITOL, TOL, ITMAX, ITER, ERR, IERR, IUNIT, RWORK, RWORK,
425
$ RWORK(LOCRGW), LENW-LOCRGW, IWORK(LOCIGW), 20,
428
IF( ITER.GT.ITMAX ) IERR = 2
430
C------------- LAST LINE OF DSLUGM FOLLOWS ----------------------------