* SUBROUTINE BICGREVCOM( 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 REAL RESID INTEGER NDX1, NDX2 REAL SCLR1, SCLR2 INTEGER IJOB * .. * .. Array Arguments .. REAL X( * ), B( * ), WORK( LDW,* ) * * .. * * Purpose * ======= * * BiCG solves the linear system Ax = b using the * BiConjugate Gradient iterative method with preconditioning. * * Arguments * ========= * * N (input) INTEGER. * On entry, the dimension of the matrix. * Unchanged on exit. * * B (input) REAL array, dimension N. * On entry, right hand side vector B. * Unchanged on exit. * * X (input/output) REAL array, dimension N. * On input, the initial guess; on exit, the iterated solution. * * WORK (workspace) REAL array, dimension (LDW,6). * Workspace for residual, direction vector, etc. * Note that Z and Q, and ZTLD and QTLD share workspace. * * LDW (input) INTEGER * The leading dimension of the array WORK. LDW >= 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) REAL * 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. * * > 0: Convergence to tolerance not achieved. This will be * set to the number of iterations performed. * * < 0: Illegal input parameter, or breakdown occurred * during iteration. * * Illegal parameter: * * -1: matrix dimension N < 0 * -2: LDW < N * -3: Maximum number of iterations ITER <= 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 < BREAKTOL: RHO and RTLD have become * orthogonal. * * BREAKTOL is set in function 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) REAL . * 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: SAXPY, SCOPY, SDOT, SNRM2, * ============================================================== * * .. Parameters .. REAL ZERO, ONE PARAMETER ( ZERO = 0.0E+0, ONE = 1.0E+0 ) * .. * .. Local Scalars .. INTEGER R, RTLD, Z, ZTLD, P, PTLD, Q, QTLD, MAXIT, $ NEED1, NEED2 REAL TOL, ALPHA, BETA, RHO, RHO1, BNRM2, RHOTOL, $ GETBREAK, SDOT, SNRM2 * * indicates where to resume from. Only valid when IJOB = 2! INTEGER RLBL * * saving all. SAVE * * .. * .. External Routines .. EXTERNAL SAXPY, SCOPY, SDOT, SNRM2 * .. * .. 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 * * init. ***************** 1 CONTINUE ***************** * INFO = 0 MAXIT = ITER TOL = RESID * * Alias workspace columns. * R = 1 RTLD = 2 Z = 3 ZTLD = 4 P = 5 PTLD = 6 Q = 3 QTLD = 4 * * 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 = ((Z - 1) * LDW) + 1 ELSEIF( NDX1.EQ.4 ) THEN NEED1 = ((ZTLD - 1) * LDW) + 1 ELSEIF( NDX1.EQ.5 ) THEN NEED1 = ((P - 1) * LDW) + 1 ELSEIF( NDX1.EQ.6 ) THEN NEED1 = ((PTLD - 1) * LDW) + 1 ELSEIF( NDX1.EQ.7 ) THEN NEED1 = ((Q - 1) * LDW) + 1 ELSEIF( NDX1.EQ.8 ) THEN NEED1 = ((QTLD - 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 = ((Z - 1) * LDW) + 1 ELSEIF( NDX2.EQ.4 ) THEN NEED2 = ((ZTLD - 1) * LDW) + 1 ELSEIF( NDX2.EQ.5 ) THEN NEED2 = ((P - 1) * LDW) + 1 ELSEIF( NDX2.EQ.6 ) THEN NEED2 = ((PTLD - 1) * LDW) + 1 ELSEIF( NDX2.EQ.7 ) THEN NEED2 = ((Q - 1) * LDW) + 1 ELSEIF( NDX2.EQ.8 ) THEN NEED2 = ((QTLD - 1) * LDW) + 1 ELSE * report error INFO = -5 GO TO 20 ENDIF ELSE NEED2 = NDX2 ENDIF * * Set breakdown parameters. * RHOTOL = GETBREAK() * * Set initial residual. * CALL SCOPY( N, B, 1, WORK(1,R), 1 ) IF ( SNRM2( N, X, 1 ).NE.ZERO ) THEN *********CALL MATVEC( -ONE, X, ZERO, WORK(1,R) ) * using WORK[RTLD] as temp *********CALL SCOPY( N, X, 1, WORK(1,RTLD), 1 ) SCLR1 = -ONE SCLR2 = ZERO NDX1 = ((RTLD - 1) * LDW) + 1 NDX2 = ((R - 1) * LDW) + 1 RLBL = 2 IJOB = 5 RETURN ***************** 2 CONTINUE ***************** * IF ( SNRM2( N, WORK(1,R), 1 ).LE.TOL ) GO TO 30 * ENDIF CALL SCOPY( N, WORK(1,R), 1, WORK(1,RTLD), 1 ) BNRM2 = SNRM2( N, B, 1 ) IF ( BNRM2.EQ.ZERO ) BNRM2 = ONE * ITER = 0 * 10 CONTINUE * * Perform BiConjugate Gradient iteration. * ITER = ITER + 1 * * Compute direction vectors PK and PTLD. * *********CALL PSOLVE( WORK(1,Z), WORK(1,R) ) * NDX1 = ((Z - 1) * LDW) + 1 NDX2 = ((R - 1) * LDW) + 1 RLBL = 3 IJOB = 3 RETURN ***************** 3 CONTINUE ***************** * *********CALL PSOLVETRANS( WORK(1,ZTLD), WORK(1,RTLD) ) * NDX1 = ((ZTLD - 1) * LDW) + 1 NDX2 = ((RTLD - 1) * LDW) + 1 RLBL = 4 IJOB = 4 RETURN ***************** 4 CONTINUE ***************** * * RHO = SDOT( N, WORK(1,Z), 1, WORK(1,RTLD), 1 ) IF ( ABS( RHO ).LT.RHOTOL ) GO TO 25 * IF ( ITER.GT.1 ) THEN BETA = RHO / RHO1 CALL SAXPY( N, BETA, WORK(1,P), 1, WORK(1,Z), 1 ) CALL SAXPY( N, BETA, WORK(1,PTLD), 1, WORK(1,ZTLD), 1 ) CALL SCOPY( N, WORK(1,Z), 1, WORK(1,P), 1 ) CALL SCOPY( N, WORK(1,ZTLD), 1, WORK(1,PTLD), 1 ) ELSE CALL SCOPY( N, WORK(1,Z), 1, WORK(1,P), 1 ) CALL SCOPY( N, WORK(1,ZTLD), 1, WORK(1,PTLD), 1 ) ENDIF * *********CALL MATVEC( ONE, WORK(1,P), ZERO, WORK(1,Q) ) * SCLR1 = ONE SCLR2 = ZERO NDX1 = ((P - 1) * LDW) + 1 NDX2 = ((Q - 1) * LDW) + 1 RLBL = 5 IJOB = 1 RETURN ***************** 5 CONTINUE ***************** * *********CALL MATVECTRANS( ONE, WORK(1,PTLD), ZERO, WORK(1,QTLD) ) * SCLR1 = ONE SCLR2 = ZERO NDX1 = ((PTLD - 1) * LDW) + 1 NDX2 = ((QTLD - 1) * LDW) + 1 RLBL = 6 IJOB = 2 RETURN ***************** 6 CONTINUE ***************** ALPHA = RHO / SDOT( N, WORK(1,PTLD), 1, WORK(1,Q), 1 ) * * Compute current solution vector x. * CALL SAXPY( N, ALPHA, WORK(1,P), 1, X, 1 ) * * Compute residual vector rk, find norm, * then check for tolerance. * CALL SAXPY( N, -ALPHA, WORK(1,Q), 1, WORK(1,R), 1 ) * *********RESID = SNRM2( N, WORK(1,R), 1 ) / BNRM2 *********IF ( RESID.LE.TOL ) GO TO 30 * NDX1 = NEED1 NDX2 = NEED2 * Prepare for resumption & return RLBL = 7 IJOB = 6 RETURN * ***************** 7 CONTINUE ***************** IF( INFO.EQ.1 ) GO TO 30 * IF ( ITER.EQ.MAXIT ) THEN INFO = 1 GO TO 20 ENDIF * CALL SAXPY( N, -ALPHA, WORK(1,QTLD), 1, WORK(1,RTLD), 1 ) RHO1 = RHO * GO TO 10 * 20 CONTINUE * * Iteration fails. * RLBL = -1 IJOB = -1 RETURN * 25 CONTINUE * * Set breakdown flag. * INFO = -10 RLBL = -1 IJOB = -1 RETURN * 30 CONTINUE * * Iteration successful; return. * INFO = 0 RLBL = -1 IJOB = -1 RETURN * * End of BICGREVCOM * END