LAPACK 3.3.1
Linear Algebra PACKage

sla_gerfsx_extended.f

Go to the documentation of this file.
00001       SUBROUTINE SLA_GERFSX_EXTENDED( PREC_TYPE, TRANS_TYPE, N, NRHS, A,
00002      $                                LDA, AF, LDAF, IPIV, COLEQU, C, B,
00003      $                                LDB, Y, LDY, BERR_OUT, N_NORMS,
00004      $                                ERR_BNDS_NORM, ERR_BNDS_COMP, RES,
00005      $                                AYB, DY, Y_TAIL, RCOND, ITHRESH,
00006      $                                RTHRESH, DZ_UB, IGNORE_CWISE,
00007      $                                INFO )
00008 *
00009 *     -- LAPACK routine (version 3.2.1)                                 --
00010 *     -- Contributed by James Demmel, Deaglan Halligan, Yozo Hida and --
00011 *     -- Jason Riedy of Univ. of California Berkeley.                 --
00012 *     -- April 2009                                                   --
00013 *
00014 *     -- LAPACK is a software package provided by Univ. of Tennessee, --
00015 *     -- Univ. of California Berkeley and NAG Ltd.                    --
00016 *
00017       IMPLICIT NONE
00018 *     ..
00019 *     .. Scalar Arguments ..
00020       INTEGER            INFO, LDA, LDAF, LDB, LDY, N, NRHS, PREC_TYPE,
00021      $                   TRANS_TYPE, N_NORMS, ITHRESH
00022       LOGICAL            COLEQU, IGNORE_CWISE
00023       REAL               RTHRESH, DZ_UB
00024 *     ..
00025 *     .. Array Arguments ..
00026       INTEGER            IPIV( * )
00027       REAL               A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
00028      $                   Y( LDY, * ), RES( * ), DY( * ), Y_TAIL( * )
00029       REAL               C( * ), AYB( * ), RCOND, BERR_OUT( * ),
00030      $                   ERR_BNDS_NORM( NRHS, * ),
00031      $                   ERR_BNDS_COMP( NRHS, * )
00032 *     ..
00033 *
00034 *  Purpose
00035 *  =======
00036 *
00037 *  SLA_GERFSX_EXTENDED improves the computed solution to a system of
00038 *  linear equations by performing extra-precise iterative refinement
00039 *  and provides error bounds and backward error estimates for the solution.
00040 *  This subroutine is called by SGERFSX to perform iterative refinement.
00041 *  In addition to normwise error bound, the code provides maximum
00042 *  componentwise error bound if possible. See comments for ERR_BNDS_NORM
00043 *  and ERR_BNDS_COMP for details of the error bounds. Note that this
00044 *  subroutine is only resonsible for setting the second fields of
00045 *  ERR_BNDS_NORM and ERR_BNDS_COMP.
00046 *
00047 *  Arguments
00048 *  =========
00049 *
00050 *     PREC_TYPE      (input) INTEGER
00051 *     Specifies the intermediate precision to be used in refinement.
00052 *     The value is defined by ILAPREC(P) where P is a CHARACTER and
00053 *     P    = 'S':  Single
00054 *          = 'D':  Double
00055 *          = 'I':  Indigenous
00056 *          = 'X', 'E':  Extra
00057 *
00058 *     TRANS_TYPE     (input) INTEGER
00059 *     Specifies the transposition operation on A.
00060 *     The value is defined by ILATRANS(T) where T is a CHARACTER and
00061 *     T    = 'N':  No transpose
00062 *          = 'T':  Transpose
00063 *          = 'C':  Conjugate transpose
00064 *
00065 *     N              (input) INTEGER
00066 *     The number of linear equations, i.e., the order of the
00067 *     matrix A.  N >= 0.
00068 *
00069 *     NRHS           (input) INTEGER
00070 *     The number of right-hand-sides, i.e., the number of columns of the
00071 *     matrix B.
00072 *
00073 *     A              (input) REAL array, dimension (LDA,N)
00074 *     On entry, the N-by-N matrix A.
00075 *
00076 *     LDA            (input) INTEGER
00077 *     The leading dimension of the array A.  LDA >= max(1,N).
00078 *
00079 *     AF             (input) REAL array, dimension (LDAF,N)
00080 *     The factors L and U from the factorization
00081 *     A = P*L*U as computed by SGETRF.
00082 *
00083 *     LDAF           (input) INTEGER
00084 *     The leading dimension of the array AF.  LDAF >= max(1,N).
00085 *
00086 *     IPIV           (input) INTEGER array, dimension (N)
00087 *     The pivot indices from the factorization A = P*L*U
00088 *     as computed by SGETRF; row i of the matrix was interchanged
00089 *     with row IPIV(i).
00090 *
00091 *     COLEQU         (input) LOGICAL
00092 *     If .TRUE. then column equilibration was done to A before calling
00093 *     this routine. This is needed to compute the solution and error
00094 *     bounds correctly.
00095 *
00096 *     C              (input) REAL array, dimension (N)
00097 *     The column scale factors for A. If COLEQU = .FALSE., C
00098 *     is not accessed. If C is input, each element of C should be a power
00099 *     of the radix to ensure a reliable solution and error estimates.
00100 *     Scaling by powers of the radix does not cause rounding errors unless
00101 *     the result underflows or overflows. Rounding errors during scaling
00102 *     lead to refining with a matrix that is not equivalent to the
00103 *     input matrix, producing error estimates that may not be
00104 *     reliable.
00105 *
00106 *     B              (input) REAL array, dimension (LDB,NRHS)
00107 *     The right-hand-side matrix B.
00108 *
00109 *     LDB            (input) INTEGER
00110 *     The leading dimension of the array B.  LDB >= max(1,N).
00111 *
00112 *     Y              (input/output) REAL array, dimension (LDY,NRHS)
00113 *     On entry, the solution matrix X, as computed by SGETRS.
00114 *     On exit, the improved solution matrix Y.
00115 *
00116 *     LDY            (input) INTEGER
00117 *     The leading dimension of the array Y.  LDY >= max(1,N).
00118 *
00119 *     BERR_OUT       (output) REAL array, dimension (NRHS)
00120 *     On exit, BERR_OUT(j) contains the componentwise relative backward
00121 *     error for right-hand-side j from the formula
00122 *         max(i) ( abs(RES(i)) / ( abs(op(A_s))*abs(Y) + abs(B_s) )(i) )
00123 *     where abs(Z) is the componentwise absolute value of the matrix
00124 *     or vector Z. This is computed by SLA_LIN_BERR.
00125 *
00126 *     N_NORMS        (input) INTEGER
00127 *     Determines which error bounds to return (see ERR_BNDS_NORM
00128 *     and ERR_BNDS_COMP).
00129 *     If N_NORMS >= 1 return normwise error bounds.
00130 *     If N_NORMS >= 2 return componentwise error bounds.
00131 *
00132 *     ERR_BNDS_NORM  (input/output) REAL array, dimension (NRHS, N_ERR_BNDS)
00133 *     For each right-hand side, this array contains information about
00134 *     various error bounds and condition numbers corresponding to the
00135 *     normwise relative error, which is defined as follows:
00136 *
00137 *     Normwise relative error in the ith solution vector:
00138 *             max_j (abs(XTRUE(j,i) - X(j,i)))
00139 *            ------------------------------
00140 *                  max_j abs(X(j,i))
00141 *
00142 *     The array is indexed by the type of error information as described
00143 *     below. There currently are up to three pieces of information
00144 *     returned.
00145 *
00146 *     The first index in ERR_BNDS_NORM(i,:) corresponds to the ith
00147 *     right-hand side.
00148 *
00149 *     The second index in ERR_BNDS_NORM(:,err) contains the following
00150 *     three fields:
00151 *     err = 1 "Trust/don't trust" boolean. Trust the answer if the
00152 *              reciprocal condition number is less than the threshold
00153 *              sqrt(n) * slamch('Epsilon').
00154 *
00155 *     err = 2 "Guaranteed" error bound: The estimated forward error,
00156 *              almost certainly within a factor of 10 of the true error
00157 *              so long as the next entry is greater than the threshold
00158 *              sqrt(n) * slamch('Epsilon'). This error bound should only
00159 *              be trusted if the previous boolean is true.
00160 *
00161 *     err = 3  Reciprocal condition number: Estimated normwise
00162 *              reciprocal condition number.  Compared with the threshold
00163 *              sqrt(n) * slamch('Epsilon') to determine if the error
00164 *              estimate is "guaranteed". These reciprocal condition
00165 *              numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
00166 *              appropriately scaled matrix Z.
00167 *              Let Z = S*A, where S scales each row by a power of the
00168 *              radix so all absolute row sums of Z are approximately 1.
00169 *
00170 *     This subroutine is only responsible for setting the second field
00171 *     above.
00172 *     See Lapack Working Note 165 for further details and extra
00173 *     cautions.
00174 *
00175 *     ERR_BNDS_COMP  (input/output) REAL array, dimension (NRHS, N_ERR_BNDS)
00176 *     For each right-hand side, this array contains information about
00177 *     various error bounds and condition numbers corresponding to the
00178 *     componentwise relative error, which is defined as follows:
00179 *
00180 *     Componentwise relative error in the ith solution vector:
00181 *                    abs(XTRUE(j,i) - X(j,i))
00182 *             max_j ----------------------
00183 *                         abs(X(j,i))
00184 *
00185 *     The array is indexed by the right-hand side i (on which the
00186 *     componentwise relative error depends), and the type of error
00187 *     information as described below. There currently are up to three
00188 *     pieces of information returned for each right-hand side. If
00189 *     componentwise accuracy is not requested (PARAMS(3) = 0.0), then
00190 *     ERR_BNDS_COMP is not accessed.  If N_ERR_BNDS .LT. 3, then at most
00191 *     the first (:,N_ERR_BNDS) entries are returned.
00192 *
00193 *     The first index in ERR_BNDS_COMP(i,:) corresponds to the ith
00194 *     right-hand side.
00195 *
00196 *     The second index in ERR_BNDS_COMP(:,err) contains the following
00197 *     three fields:
00198 *     err = 1 "Trust/don't trust" boolean. Trust the answer if the
00199 *              reciprocal condition number is less than the threshold
00200 *              sqrt(n) * slamch('Epsilon').
00201 *
00202 *     err = 2 "Guaranteed" error bound: The estimated forward error,
00203 *              almost certainly within a factor of 10 of the true error
00204 *              so long as the next entry is greater than the threshold
00205 *              sqrt(n) * slamch('Epsilon'). This error bound should only
00206 *              be trusted if the previous boolean is true.
00207 *
00208 *     err = 3  Reciprocal condition number: Estimated componentwise
00209 *              reciprocal condition number.  Compared with the threshold
00210 *              sqrt(n) * slamch('Epsilon') to determine if the error
00211 *              estimate is "guaranteed". These reciprocal condition
00212 *              numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
00213 *              appropriately scaled matrix Z.
00214 *              Let Z = S*(A*diag(x)), where x is the solution for the
00215 *              current right-hand side and S scales each row of
00216 *              A*diag(x) by a power of the radix so all absolute row
00217 *              sums of Z are approximately 1.
00218 *
00219 *     This subroutine is only responsible for setting the second field
00220 *     above.
00221 *     See Lapack Working Note 165 for further details and extra
00222 *     cautions.
00223 *
00224 *     RES            (input) REAL array, dimension (N)
00225 *     Workspace to hold the intermediate residual.
00226 *
00227 *     AYB            (input) REAL array, dimension (N)
00228 *     Workspace. This can be the same workspace passed for Y_TAIL.
00229 *
00230 *     DY             (input) REAL array, dimension (N)
00231 *     Workspace to hold the intermediate solution.
00232 *
00233 *     Y_TAIL         (input) REAL array, dimension (N)
00234 *     Workspace to hold the trailing bits of the intermediate solution.
00235 *
00236 *     RCOND          (input) REAL
00237 *     Reciprocal scaled condition number.  This is an estimate of the
00238 *     reciprocal Skeel condition number of the matrix A after
00239 *     equilibration (if done).  If this is less than the machine
00240 *     precision (in particular, if it is zero), the matrix is singular
00241 *     to working precision.  Note that the error may still be small even
00242 *     if this number is very small and the matrix appears ill-
00243 *     conditioned.
00244 *
00245 *     ITHRESH        (input) INTEGER
00246 *     The maximum number of residual computations allowed for
00247 *     refinement. The default is 10. For 'aggressive' set to 100 to
00248 *     permit convergence using approximate factorizations or
00249 *     factorizations other than LU. If the factorization uses a
00250 *     technique other than Gaussian elimination, the guarantees in
00251 *     ERR_BNDS_NORM and ERR_BNDS_COMP may no longer be trustworthy.
00252 *
00253 *     RTHRESH        (input) REAL
00254 *     Determines when to stop refinement if the error estimate stops
00255 *     decreasing. Refinement will stop when the next solution no longer
00256 *     satisfies norm(dx_{i+1}) < RTHRESH * norm(dx_i) where norm(Z) is
00257 *     the infinity norm of Z. RTHRESH satisfies 0 < RTHRESH <= 1. The
00258 *     default value is 0.5. For 'aggressive' set to 0.9 to permit
00259 *     convergence on extremely ill-conditioned matrices. See LAWN 165
00260 *     for more details.
00261 *
00262 *     DZ_UB          (input) REAL
00263 *     Determines when to start considering componentwise convergence.
00264 *     Componentwise convergence is only considered after each component
00265 *     of the solution Y is stable, which we definte as the relative
00266 *     change in each component being less than DZ_UB. The default value
00267 *     is 0.25, requiring the first bit to be stable. See LAWN 165 for
00268 *     more details.
00269 *
00270 *     IGNORE_CWISE   (input) LOGICAL
00271 *     If .TRUE. then ignore componentwise convergence. Default value
00272 *     is .FALSE..
00273 *
00274 *     INFO           (output) INTEGER
00275 *       = 0:  Successful exit.
00276 *       < 0:  if INFO = -i, the ith argument to SGETRS had an illegal
00277 *             value
00278 *
00279 *  =====================================================================
00280 *
00281 *     .. Local Scalars ..
00282       CHARACTER          TRANS
00283       INTEGER            CNT, I, J, X_STATE, Z_STATE, Y_PREC_STATE
00284       REAL               YK, DYK, YMIN, NORMY, NORMX, NORMDX, DXRAT,
00285      $                   DZRAT, PREVNORMDX, PREV_DZ_Z, DXRATMAX,
00286      $                   DZRATMAX, DX_X, DZ_Z, FINAL_DX_X, FINAL_DZ_Z,
00287      $                   EPS, HUGEVAL, INCR_THRESH
00288       LOGICAL            INCR_PREC
00289 *     ..
00290 *     .. Parameters ..
00291       INTEGER            UNSTABLE_STATE, WORKING_STATE, CONV_STATE,
00292      $                   NOPROG_STATE, BASE_RESIDUAL, EXTRA_RESIDUAL,
00293      $                   EXTRA_Y
00294       PARAMETER          ( UNSTABLE_STATE = 0, WORKING_STATE = 1,
00295      $                   CONV_STATE = 2, NOPROG_STATE = 3 )
00296       PARAMETER          ( BASE_RESIDUAL = 0, EXTRA_RESIDUAL = 1,
00297      $                   EXTRA_Y = 2 )
00298       INTEGER            FINAL_NRM_ERR_I, FINAL_CMP_ERR_I, BERR_I
00299       INTEGER            RCOND_I, NRM_RCOND_I, NRM_ERR_I, CMP_RCOND_I
00300       INTEGER            CMP_ERR_I, PIV_GROWTH_I
00301       PARAMETER          ( FINAL_NRM_ERR_I = 1, FINAL_CMP_ERR_I = 2,
00302      $                   BERR_I = 3 )
00303       PARAMETER          ( RCOND_I = 4, NRM_RCOND_I = 5, NRM_ERR_I = 6 )
00304       PARAMETER          ( CMP_RCOND_I = 7, CMP_ERR_I = 8,
00305      $                   PIV_GROWTH_I = 9 )
00306       INTEGER            LA_LINRX_ITREF_I, LA_LINRX_ITHRESH_I,
00307      $                   LA_LINRX_CWISE_I
00308       PARAMETER          ( LA_LINRX_ITREF_I = 1,
00309      $                   LA_LINRX_ITHRESH_I = 2 )
00310       PARAMETER          ( LA_LINRX_CWISE_I = 3 )
00311       INTEGER            LA_LINRX_TRUST_I, LA_LINRX_ERR_I,
00312      $                   LA_LINRX_RCOND_I
00313       PARAMETER          ( LA_LINRX_TRUST_I = 1, LA_LINRX_ERR_I = 2 )
00314       PARAMETER          ( LA_LINRX_RCOND_I = 3 )
00315 *     ..
00316 *     .. External Subroutines ..
00317       EXTERNAL           SAXPY, SCOPY, SGETRS, SGEMV, BLAS_SGEMV_X,
00318      $                   BLAS_SGEMV2_X, SLA_GEAMV, SLA_WWADDW, SLAMCH,
00319      $                   CHLA_TRANSTYPE, SLA_LIN_BERR
00320       REAL               SLAMCH
00321       CHARACTER          CHLA_TRANSTYPE
00322 *     ..
00323 *     .. Intrinsic Functions ..
00324       INTRINSIC          ABS, MAX, MIN
00325 *     ..
00326 *     .. Executable Statements ..
00327 *
00328       IF ( INFO.NE.0 ) RETURN
00329       TRANS = CHLA_TRANSTYPE(TRANS_TYPE)
00330       EPS = SLAMCH( 'Epsilon' )
00331       HUGEVAL = SLAMCH( 'Overflow' )
00332 *     Force HUGEVAL to Inf
00333       HUGEVAL = HUGEVAL * HUGEVAL
00334 *     Using HUGEVAL may lead to spurious underflows.
00335       INCR_THRESH = REAL( N ) * EPS
00336 *
00337       DO J = 1, NRHS
00338          Y_PREC_STATE = EXTRA_RESIDUAL
00339          IF ( Y_PREC_STATE .EQ. EXTRA_Y ) THEN
00340             DO I = 1, N
00341                Y_TAIL( I ) = 0.0
00342             END DO
00343          END IF
00344 
00345          DXRAT = 0.0
00346          DXRATMAX = 0.0
00347          DZRAT = 0.0
00348          DZRATMAX = 0.0
00349          FINAL_DX_X = HUGEVAL
00350          FINAL_DZ_Z = HUGEVAL
00351          PREVNORMDX = HUGEVAL
00352          PREV_DZ_Z = HUGEVAL
00353          DZ_Z = HUGEVAL
00354          DX_X = HUGEVAL
00355 
00356          X_STATE = WORKING_STATE
00357          Z_STATE = UNSTABLE_STATE
00358          INCR_PREC = .FALSE.
00359 
00360          DO CNT = 1, ITHRESH
00361 *
00362 *         Compute residual RES = B_s - op(A_s) * Y,
00363 *             op(A) = A, A**T, or A**H depending on TRANS (and type).
00364 *
00365             CALL SCOPY( N, B( 1, J ), 1, RES, 1 )
00366             IF ( Y_PREC_STATE .EQ. BASE_RESIDUAL ) THEN
00367                CALL SGEMV( TRANS, N, N, -1.0, A, LDA, Y( 1, J ), 1,
00368      $              1.0, RES, 1 )
00369             ELSE IF ( Y_PREC_STATE .EQ. EXTRA_RESIDUAL ) THEN
00370                CALL BLAS_SGEMV_X( TRANS_TYPE, N, N, -1.0, A, LDA,
00371      $              Y( 1, J ), 1, 1.0, RES, 1, PREC_TYPE )
00372             ELSE
00373                CALL BLAS_SGEMV2_X( TRANS_TYPE, N, N, -1.0, A, LDA,
00374      $              Y( 1, J ), Y_TAIL, 1, 1.0, RES, 1, PREC_TYPE )
00375             END IF
00376 
00377 !        XXX: RES is no longer needed.
00378             CALL SCOPY( N, RES, 1, DY, 1 )
00379             CALL SGETRS( TRANS, N, 1, AF, LDAF, IPIV, DY, N, INFO )
00380 *
00381 *         Calculate relative changes DX_X, DZ_Z and ratios DXRAT, DZRAT.
00382 *
00383             NORMX = 0.0
00384             NORMY = 0.0
00385             NORMDX = 0.0
00386             DZ_Z = 0.0
00387             YMIN = HUGEVAL
00388 *
00389             DO I = 1, N
00390                YK = ABS( Y( I, J ) )
00391                DYK = ABS( DY( I ) )
00392 
00393                IF ( YK .NE. 0.0 ) THEN
00394                   DZ_Z = MAX( DZ_Z, DYK / YK )
00395                ELSE IF ( DYK .NE. 0.0 ) THEN
00396                   DZ_Z = HUGEVAL
00397                END IF
00398 
00399                YMIN = MIN( YMIN, YK )
00400 
00401                NORMY = MAX( NORMY, YK )
00402 
00403                IF ( COLEQU ) THEN
00404                   NORMX = MAX( NORMX, YK * C( I ) )
00405                   NORMDX = MAX( NORMDX, DYK * C( I ) )
00406                ELSE
00407                   NORMX = NORMY
00408                   NORMDX = MAX( NORMDX, DYK )
00409                END IF
00410             END DO
00411 
00412             IF ( NORMX .NE. 0.0 ) THEN
00413                DX_X = NORMDX / NORMX
00414             ELSE IF ( NORMDX .EQ. 0.0 ) THEN
00415                DX_X = 0.0
00416             ELSE
00417                DX_X = HUGEVAL
00418             END IF
00419 
00420             DXRAT = NORMDX / PREVNORMDX
00421             DZRAT = DZ_Z / PREV_DZ_Z
00422 *
00423 *         Check termination criteria
00424 *
00425             IF (.NOT.IGNORE_CWISE
00426      $           .AND. YMIN*RCOND .LT. INCR_THRESH*NORMY
00427      $           .AND. Y_PREC_STATE .LT. EXTRA_Y)
00428      $           INCR_PREC = .TRUE.
00429 
00430             IF ( X_STATE .EQ. NOPROG_STATE .AND. DXRAT .LE. RTHRESH )
00431      $           X_STATE = WORKING_STATE
00432             IF ( X_STATE .EQ. WORKING_STATE ) THEN
00433                IF ( DX_X .LE. EPS ) THEN
00434                   X_STATE = CONV_STATE
00435                ELSE IF ( DXRAT .GT. RTHRESH ) THEN
00436                   IF ( Y_PREC_STATE .NE. EXTRA_Y ) THEN
00437                      INCR_PREC = .TRUE.
00438                   ELSE
00439                      X_STATE = NOPROG_STATE
00440                   END IF
00441                ELSE
00442                   IF ( DXRAT .GT. DXRATMAX ) DXRATMAX = DXRAT
00443                END IF
00444                IF ( X_STATE .GT. WORKING_STATE ) FINAL_DX_X = DX_X
00445             END IF
00446 
00447             IF ( Z_STATE .EQ. UNSTABLE_STATE .AND. DZ_Z .LE. DZ_UB )
00448      $           Z_STATE = WORKING_STATE
00449             IF ( Z_STATE .EQ. NOPROG_STATE .AND. DZRAT .LE. RTHRESH )
00450      $           Z_STATE = WORKING_STATE
00451             IF ( Z_STATE .EQ. WORKING_STATE ) THEN
00452                IF ( DZ_Z .LE. EPS ) THEN
00453                   Z_STATE = CONV_STATE
00454                ELSE IF ( DZ_Z .GT. DZ_UB ) THEN
00455                   Z_STATE = UNSTABLE_STATE
00456                   DZRATMAX = 0.0
00457                   FINAL_DZ_Z = HUGEVAL
00458                ELSE IF ( DZRAT .GT. RTHRESH ) THEN
00459                   IF ( Y_PREC_STATE .NE. EXTRA_Y ) THEN
00460                      INCR_PREC = .TRUE.
00461                   ELSE
00462                      Z_STATE = NOPROG_STATE
00463                   END IF
00464                ELSE
00465                   IF ( DZRAT .GT. DZRATMAX ) DZRATMAX = DZRAT
00466                END IF
00467                IF ( Z_STATE .GT. WORKING_STATE ) FINAL_DZ_Z = DZ_Z
00468             END IF
00469 *
00470 *           Exit if both normwise and componentwise stopped working,
00471 *           but if componentwise is unstable, let it go at least two
00472 *           iterations.
00473 *
00474             IF ( X_STATE.NE.WORKING_STATE ) THEN
00475                IF ( IGNORE_CWISE) GOTO 666
00476                IF ( Z_STATE.EQ.NOPROG_STATE .OR. Z_STATE.EQ.CONV_STATE )
00477      $              GOTO 666
00478                IF ( Z_STATE.EQ.UNSTABLE_STATE .AND. CNT.GT.1 ) GOTO 666
00479             END IF
00480 
00481             IF ( INCR_PREC ) THEN
00482                INCR_PREC = .FALSE.
00483                Y_PREC_STATE = Y_PREC_STATE + 1
00484                DO I = 1, N
00485                   Y_TAIL( I ) = 0.0
00486                END DO
00487             END IF
00488 
00489             PREVNORMDX = NORMDX
00490             PREV_DZ_Z = DZ_Z
00491 *
00492 *           Update soluton.
00493 *
00494             IF ( Y_PREC_STATE .LT. EXTRA_Y ) THEN
00495                CALL SAXPY( N, 1.0, DY, 1, Y( 1, J ), 1 )
00496             ELSE
00497                CALL SLA_WWADDW( N, Y( 1, J ), Y_TAIL, DY )
00498             END IF
00499 
00500          END DO
00501 *        Target of "IF (Z_STOP .AND. X_STOP)".  Sun's f77 won't EXIT.
00502  666     CONTINUE
00503 *
00504 *     Set final_* when cnt hits ithresh.
00505 *
00506          IF ( X_STATE .EQ. WORKING_STATE ) FINAL_DX_X = DX_X
00507          IF ( Z_STATE .EQ. WORKING_STATE ) FINAL_DZ_Z = DZ_Z
00508 *
00509 *     Compute error bounds
00510 *
00511          IF (N_NORMS .GE. 1) THEN
00512             ERR_BNDS_NORM( J, LA_LINRX_ERR_I ) =
00513      $           FINAL_DX_X / (1 - DXRATMAX)
00514          END IF
00515          IF ( N_NORMS .GE. 2 ) THEN
00516             ERR_BNDS_COMP( J, LA_LINRX_ERR_I ) =
00517      $           FINAL_DZ_Z / (1 - DZRATMAX)
00518          END IF
00519 *
00520 *     Compute componentwise relative backward error from formula
00521 *         max(i) ( abs(R(i)) / ( abs(op(A_s))*abs(Y) + abs(B_s) )(i) )
00522 *     where abs(Z) is the componentwise absolute value of the matrix
00523 *     or vector Z.
00524 *
00525 *         Compute residual RES = B_s - op(A_s) * Y,
00526 *             op(A) = A, A**T, or A**H depending on TRANS (and type).
00527 *
00528          CALL SCOPY( N, B( 1, J ), 1, RES, 1 )
00529          CALL SGEMV( TRANS, N, N, -1.0, A, LDA, Y(1,J), 1, 1.0, RES, 1 )
00530 
00531          DO I = 1, N
00532             AYB( I ) = ABS( B( I, J ) )
00533          END DO
00534 *
00535 *     Compute abs(op(A_s))*abs(Y) + abs(B_s).
00536 *
00537          CALL SLA_GEAMV ( TRANS_TYPE, N, N, 1.0,
00538      $        A, LDA, Y(1, J), 1, 1.0, AYB, 1 )
00539 
00540          CALL SLA_LIN_BERR ( N, N, 1, RES, AYB, BERR_OUT( J ) )
00541 *
00542 *     End of loop for each RHS.
00543 *
00544       END DO
00545 *
00546       RETURN
00547       END
 All Files Functions