File: BiCGSTABREVCOM.f.src

package info (click to toggle)
python-scipy 0.18.1-2
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 75,464 kB
  • ctags: 79,406
  • sloc: python: 143,495; cpp: 89,357; fortran: 81,650; ansic: 79,778; makefile: 364; sh: 265
file content (428 lines) | stat: -rw-r--r-- 11,737 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
* -*- fortran -*-
      SUBROUTINE <_c>BICGSTABREVCOM(N, B, X, WORK, LDW, ITER, RESID, 
     $                    INFO,NDX1, NDX2, SCLR1, SCLR2, IJOB)
*
*  -- Iterative template routine --
*     Univ. of Tennessee and Oak Ridge National Laboratory
*     October 1, 1993
*     Details of this algorithm are described in "Templates for the 
*     Solution of Linear Systems: Building Blocks for Iterative 
*     Methods", Barrett, Berry, Chan, Demmel, Donato, Dongarra, 
*     Eijkhout, Pozo, Romine, and van der Vorst, SIAM Publications,
*     1993. (ftp netlib2.cs.utk.edu; cd linalg; get templates.ps).
*
*     .. Scalar Arguments ..
      INTEGER            N, LDW, ITER, INFO 
      <rt=real,double precision,real,double precision>   RESID
      INTEGER            NDX1, NDX2
      <_t>   SCLR1, SCLR2
      INTEGER            IJOB
*     ..
*     .. Array Arguments ..
      <_t>   X( * ), B( * ), WORK( LDW,* )
*     ..
*
*  Purpose
*  =======
*
*  BICGSTAB solves the linear system A*x = b using the
*  BiConjugate Gradient Stabilized iterative method with 
*  preconditioning.
*
*  Arguments
*  =========
*
*  N       (input) INTEGER. 
*          On entry, the dimension of the matrix.
*          Unchanged on exit.
* 
*  B       (input) DOUBLE PRECISION array, dimension N.
*          On entry, right hand side vector B.
*          Unchanged on exit.
*
*  X       (input/output) DOUBLE PRECISION array, dimension N.
*          On input, the initial guess. This is commonly set to 
*          the zero vector. 
*          On exit, if INFO = 0, the iterated approximate solution.
*
*  WORK    (workspace) DOUBLE PRECISION array, dimension (LDW,7)
*          Workspace for residual, direction vector, etc.
*          Note that vectors R and S shared the same workspace.
*
*  LDW     (input) INTEGER
*          The leading dimension of the array WORK. LDW .gt. = max(1,N).
*
*  ITER    (input/output) INTEGER
*          On input, the maximum iterations to be performed.
*          On output, actual number of iterations performed.
*
*  RESID   (input/output) DOUBLE PRECISION
*          On input, the allowable convergence measure for
*          norm( b - A*x ) / norm( b ).
*          On output, the final value of this measure.
*
*  INFO    (output) INTEGER
*
*          =  0: Successful exit. Iterated approximate solution returned.
*
*          .gt.   0: Convergence to tolerance not achieved. This will be 
*                set to the number of iterations performed.
*
*           .ls.   0: Illegal input parameter, or breakdown occurred
*                during iteration.
*
*                Illegal parameter:
*
*                   -1: matrix dimension N  .ls.  0
*                   -2: LDW  .ls.  N
*                   -3: Maximum number of iterations ITER  .ls. = 0.
*                   -5: Erroneous NDX1/NDX2 in INIT call.
*                   -6: Erroneous RLBL.
*
*                BREAKDOWN: If parameters RHO or OMEGA become smaller
*                   than some tolerance, the program will terminate.
*                   Here we check against tolerance BREAKTOL.
*
*                  -10: RHO  .ls.  BREAKTOL: RHO and RTLD have become 
*                                       orthogonal.
*                  -11: OMEGA  .ls.  BREAKTOL: S and T have become 
*                                         orthogonal relative to T'*T.
*
*                  BREAKTOL is set in func GETBREAK.
*
*  NDX1    (input/output) INTEGER. 
*  NDX2    On entry in INIT call contain indices required by interface
*          level for stopping test.
*          All other times, used as output, to indicate indices into
*          WORK[] for the MATVEC, PSOLVE done by the interface level.
*
*  SCLR1   (output) DOUBLE PRECISION.
*  SCLR2   Used to pass the scalars used in MATVEC. Scalars are reqd because
*          original routines use dgemv.
*
*  IJOB    (input/output) INTEGER. 
*          Used to communicate job code between the two levels.
*
*  BLAS CALLS: DAXPY, DCOPY, DDOT, DNRM2, DSCAL
*  ==============================================================
*
*     .. Parameters ..
      <rt>             ZERO, ONE
      PARAMETER        ( ZERO = 0.0D+0, ONE = 1.0D+0 )
*     ..
*     .. Local Scalars ..
      INTEGER            R, RTLD, P, PHAT, V, S, SHAT, T, MAXIT,
     $                   NEED1, NEED2
      <rt>   TOL, BNRM2, 
     $     RHOTOL, OMEGATOL, <sdsd=s,d,s,d>GETBREAK, 
     $     <rc=ws,d,wsc,dz>NRM2
      <_t>   ALPHA, BETA, RHO, RHO1, OMEGA, TMPVAL,
     $     <xdot=wsdot,ddot,wcdotc,wzdotc>
*     indicates where to resume from. Only valid when IJOB = 2!
      INTEGER RLBL
*
*     saving all.
      SAVE
*     ..
*     .. External Funcs ..
      EXTERNAL           <sdsd>GETBREAK, <_c>AXPY, <_c>COPY, 
     $      <xdot>, <rc>NRM2, <_c>SCAL
*     ..
*     .. Intrinsic Funcs ..
      INTRINSIC          ABS, MAX
*     ..
*     .. Executable Statements ..
*
*     Entry point, so test IJOB
      IF (IJOB .eq. 1) THEN
         GOTO 1
      ELSEIF (IJOB .eq. 2) THEN
*        here we do resumption handling
         IF (RLBL .eq. 2) GOTO 2
         IF (RLBL .eq. 3) GOTO 3
         IF (RLBL .eq. 4) GOTO 4
         IF (RLBL .eq. 5) GOTO 5
         IF (RLBL .eq. 6) GOTO 6
         IF (RLBL .eq. 7) GOTO 7
*        if neither of these, then error
         INFO = -6
         GOTO 20
      ENDIF
*
*
*****************
 1    CONTINUE
*****************
*
      INFO = 0
      MAXIT = ITER
      TOL   = RESID
*
*     Alias workspace columns.
*
      R    = 1
      RTLD = 2
      P    = 3
      V    = 4
      T    = 5
      PHAT = 6
      SHAT = 7
      S    = 1
*
*     Check if caller will need indexing info.
*
      IF( NDX1.NE.-1 ) THEN
         IF( NDX1.EQ.1 ) THEN
            NEED1 = ((R - 1) * LDW) + 1
         ELSEIF( NDX1.EQ.2 ) THEN
            NEED1 = ((RTLD - 1) * LDW) + 1
         ELSEIF( NDX1.EQ.3 ) THEN
            NEED1 = ((P - 1) * LDW) + 1
         ELSEIF( NDX1.EQ.4 ) THEN
            NEED1 = ((V - 1) * LDW) + 1
         ELSEIF( NDX1.EQ.5 ) THEN
            NEED1 = ((T - 1) * LDW) + 1
         ELSEIF( NDX1.EQ.6 ) THEN
            NEED1 = ((PHAT - 1) * LDW) + 1
         ELSEIF( NDX1.EQ.7 ) THEN
            NEED1 = ((SHAT - 1) * LDW) + 1
         ELSEIF( NDX1.EQ.8 ) THEN
            NEED1 = ((S - 1) * LDW) + 1
         ELSE
*           report error
            INFO = -5
            GO TO 20
         ENDIF
      ELSE
         NEED1 = NDX1
      ENDIF
*
      IF( NDX2.NE.-1 ) THEN
         IF( NDX2.EQ.1 ) THEN
            NEED2 = ((R - 1) * LDW) + 1
         ELSEIF( NDX2.EQ.2 ) THEN
            NEED2 = ((RTLD - 1) * LDW) + 1
         ELSEIF( NDX2.EQ.3 ) THEN
            NEED2 = ((P - 1) * LDW) + 1
         ELSEIF( NDX2.EQ.4 ) THEN
            NEED2 = ((V - 1) * LDW) + 1
         ELSEIF( NDX2.EQ.5 ) THEN
            NEED2 = ((T - 1) * LDW) + 1
         ELSEIF( NDX2.EQ.6 ) THEN
            NEED2 = ((PHAT - 1) * LDW) + 1
         ELSEIF( NDX2.EQ.7 ) THEN
            NEED2 = ((SHAT - 1) * LDW) + 1
         ELSEIF( NDX2.EQ.8 ) THEN
            NEED2 = ((S - 1) * LDW) + 1
         ELSE
*           report error
            INFO = -5
            GO TO 20
         ENDIF
      ELSE
         NEED2 = NDX2
      ENDIF
*
*     Set parameter tolerances.
*
      RHOTOL = <sdsd>GETBREAK()
      OMEGATOL = <sdsd>GETBREAK()
*
*     Set initial residual.
*
      CALL <_c>COPY( N, B, 1, WORK(1,R), 1 )
      IF ( <rc>NRM2( N, X, 1 ).NE.ZERO ) THEN
*********CALL <_c>MATVEC( -ONE, X, ONE, WORK(1,R) )
*        Note: using RTLD[] as temp. storage.
*********CALL <_c>COPY(N, X, 1, WORK(1,RTLD), 1)
         SCLR1 = -ONE
         SCLR2 = ONE
         NDX1 = -1
         NDX2 = ((R - 1) * LDW) + 1
*
*        Prepare for resumption & return
         RLBL = 2
         IJOB = 3
         RETURN
      ENDIF
*
*****************
 2    CONTINUE
*****************
*
      IF ( <rc>NRM2( N, WORK(1,R), 1 ).LE.TOL ) GO TO 30

      CALL <_c>COPY( N, WORK(1,R), 1, WORK(1,RTLD), 1 )
*
      BNRM2 = <rc>NRM2( N, B, 1 )
      IF ( BNRM2 .EQ. ZERO ) BNRM2 = ONE
*
      ITER = 0
*
   10 CONTINUE
*
*     Perform BiConjugate Gradient Stabilized iteration.
*
      ITER = ITER + 1
*     
      RHO = <xdot>( N, WORK(1,RTLD), 1, WORK(1,R), 1 )
      IF ( ABS( RHO ).LT.RHOTOL ) GO TO 25
*     
*        Compute vector P.
*
      IF ( ITER.GT.1 ) THEN
         BETA = ( RHO / RHO1 ) * ( ALPHA / OMEGA )
         CALL <_c>AXPY( N, -OMEGA, WORK(1,V), 1, WORK(1,P), 1 )
         CALL <_c>SCAL( N, BETA, WORK(1,P), 1 )
         TMPVAL = ONE
         CALL <_c>AXPY( N, TMPVAL, WORK(1,R), 1, WORK(1,P), 1 )
      ELSE
         CALL <_c>COPY( N, WORK(1,R), 1, WORK(1,P), 1 )
      ENDIF
*
*        Compute direction adjusting vector PHAT and scalar ALPHA.
*
*********CALL PSOLVE( WORK(1,PHAT), WORK(1,P) )
*
      NDX1 = ((PHAT - 1) * LDW) + 1
      NDX2 = ((P    - 1) * LDW) + 1
*     Prepare for return & return
      RLBL = 3
      IJOB = 2
      RETURN
*
*****************
 3    CONTINUE
*****************
*
*********CALL MATVEC( ONE, WORK(1,PHAT), ZERO, WORK(1,V) )
*
      NDX1 = ((PHAT - 1) * LDW) + 1
      NDX2 = ((V    - 1) * LDW) + 1
*        Prepare for return & return
      SCLR1 = ONE
      SCLR2 = ZERO
      RLBL = 4
      IJOB = 1
      RETURN
*
*****************
 4    CONTINUE
*****************
*
      ALPHA = RHO / <xdot>( N, WORK(1,RTLD), 1, WORK(1,V), 1 )
*
*        Early check for tolerance.
*
      CALL <_c>AXPY( N, -ALPHA, WORK(1,V), 1, WORK(1,R), 1 )
      CALL <_c>COPY( N, WORK(1,R), 1, WORK(1,S), 1 )
      IF ( <rc>NRM2( N, WORK(1,S), 1 ).LE.TOL ) THEN
         CALL <_c>AXPY( N, ALPHA, WORK(1,PHAT), 1, X, 1 )
         RESID = <rc>NRM2( N, WORK(1,S), 1 ) / BNRM2
         GO TO 30
      ELSE
*
*           Compute stabilizer vector SHAT and scalar OMEGA.
*
************CALL PSOLVE( WORK(1,SHAT), WORK(1,S) )
*
         NDX1 = ((SHAT - 1) * LDW) + 1
         NDX2 = ((S    - 1) * LDW) + 1
*     Prepare for return & return
         RLBL = 5
         IJOB = 2
         RETURN
      ENDIF
*
*****************
 5    CONTINUE
*****************
*
************CALL MATVEC( ONE, WORK(1,SHAT), ZERO, WORK(1,T) )
*
      NDX1 = ((SHAT - 1) * LDW) + 1
      NDX2 = ((T    - 1) * LDW) + 1
*           Prepare for return & return
      SCLR1 = ONE
      SCLR2 = ZERO
      RLBL = 6
      IJOB = 1
      RETURN
*
*****************
 6    CONTINUE
*****************
*
      OMEGA = <xdot>( N, WORK(1,T), 1, WORK(1,S), 1 ) / 
     $     <xdot>( N, WORK(1,T), 1, WORK(1,T), 1 )
*
*           Compute new solution approximation vector X.
*
      CALL <_c>AXPY( N, ALPHA, WORK(1,PHAT), 1, X, 1 )
      CALL <_c>AXPY( N, OMEGA, WORK(1,SHAT), 1, X, 1 )
*     
*     Compute residual R, check for tolerance.
*
      CALL <_c>AXPY( N, -OMEGA, WORK(1,T), 1, WORK(1,R), 1 )
*
************RESID = DNRM2( N, WORK(1,R), 1 ) / BNRM2 
************IF ( RESID.LE.TOL  ) GO TO 30
*
      NDX1 = NEED1
      NDX2 = NEED2
*     Prepare for resumption & return
      RLBL = 7
      IJOB = 4
      RETURN
*
*****************
 7    CONTINUE
*****************
      IF( INFO.EQ.1 ) GO TO 30
*     
      IF ( ITER.EQ.MAXIT ) THEN
         INFO = 1
         GO TO 20
      ENDIF
*
      IF ( ABS( OMEGA ).LT.OMEGATOL ) THEN
         GO TO 25
      ELSE
         RHO1 = RHO
         GO TO 10
      ENDIF
*
   20 CONTINUE
*
*     Iteration fails.
*
      RLBL = -1
      IJOB = -1
      RETURN
*
   25 CONTINUE
*
*     Set breakdown flag.
*
      IF ( ABS( RHO ).LT.RHOTOL ) THEN
         INFO = -10
      ELSE IF ( ABS( OMEGA ).LT.OMEGATOL ) THEN
         INFO = -11
      ENDIF
      RLBL = -1
      IJOB = -1
      RETURN
*
   30 CONTINUE
*
*     Iteration successful; return.
*
      INFO = 0
      RLBL = -1
      IJOB = -1
      RETURN
*
*     End of BICGSTABREVCOM
*
      END
*     END SUBROUTINE <_c>BICGSTABREVCOM