LAPACK 3.3.0

dgbsvx.f

Go to the documentation of this file.
00001       SUBROUTINE DGBSVX( FACT, TRANS, N, KL, KU, NRHS, AB, LDAB, AFB,
00002      $                   LDAFB, IPIV, EQUED, R, C, B, LDB, X, LDX,
00003      $                   RCOND, FERR, BERR, WORK, IWORK, INFO )
00004 *
00005 *  -- LAPACK driver routine (version 3.2) --
00006 *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
00007 *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
00008 *     November 2006
00009 *
00010 *     .. Scalar Arguments ..
00011       CHARACTER          EQUED, FACT, TRANS
00012       INTEGER            INFO, KL, KU, LDAB, LDAFB, LDB, LDX, N, NRHS
00013       DOUBLE PRECISION   RCOND
00014 *     ..
00015 *     .. Array Arguments ..
00016       INTEGER            IPIV( * ), IWORK( * )
00017       DOUBLE PRECISION   AB( LDAB, * ), AFB( LDAFB, * ), B( LDB, * ),
00018      $                   BERR( * ), C( * ), FERR( * ), R( * ),
00019      $                   WORK( * ), X( LDX, * )
00020 *     ..
00021 *
00022 *  Purpose
00023 *  =======
00024 *
00025 *  DGBSVX uses the LU factorization to compute the solution to a real
00026 *  system of linear equations A * X = B, A**T * X = B, or A**H * X = B,
00027 *  where A is a band matrix of order N with KL subdiagonals and KU
00028 *  superdiagonals, and X and B are N-by-NRHS matrices.
00029 *
00030 *  Error bounds on the solution and a condition estimate are also
00031 *  provided.
00032 *
00033 *  Description
00034 *  ===========
00035 *
00036 *  The following steps are performed by this subroutine:
00037 *
00038 *  1. If FACT = 'E', real scaling factors are computed to equilibrate
00039 *     the system:
00040 *        TRANS = 'N':  diag(R)*A*diag(C)     *inv(diag(C))*X = diag(R)*B
00041 *        TRANS = 'T': (diag(R)*A*diag(C))**T *inv(diag(R))*X = diag(C)*B
00042 *        TRANS = 'C': (diag(R)*A*diag(C))**H *inv(diag(R))*X = diag(C)*B
00043 *     Whether or not the system will be equilibrated depends on the
00044 *     scaling of the matrix A, but if equilibration is used, A is
00045 *     overwritten by diag(R)*A*diag(C) and B by diag(R)*B (if TRANS='N')
00046 *     or diag(C)*B (if TRANS = 'T' or 'C').
00047 *
00048 *  2. If FACT = 'N' or 'E', the LU decomposition is used to factor the
00049 *     matrix A (after equilibration if FACT = 'E') as
00050 *        A = L * U,
00051 *     where L is a product of permutation and unit lower triangular
00052 *     matrices with KL subdiagonals, and U is upper triangular with
00053 *     KL+KU superdiagonals.
00054 *
00055 *  3. If some U(i,i)=0, so that U is exactly singular, then the routine
00056 *     returns with INFO = i. Otherwise, the factored form of A is used
00057 *     to estimate the condition number of the matrix A.  If the
00058 *     reciprocal of the condition number is less than machine precision,
00059 *     INFO = N+1 is returned as a warning, but the routine still goes on
00060 *     to solve for X and compute error bounds as described below.
00061 *
00062 *  4. The system of equations is solved for X using the factored form
00063 *     of A.
00064 *
00065 *  5. Iterative refinement is applied to improve the computed solution
00066 *     matrix and calculate error bounds and backward error estimates
00067 *     for it.
00068 *
00069 *  6. If equilibration was used, the matrix X is premultiplied by
00070 *     diag(C) (if TRANS = 'N') or diag(R) (if TRANS = 'T' or 'C') so
00071 *     that it solves the original system before equilibration.
00072 *
00073 *  Arguments
00074 *  =========
00075 *
00076 *  FACT    (input) CHARACTER*1
00077 *          Specifies whether or not the factored form of the matrix A is
00078 *          supplied on entry, and if not, whether the matrix A should be
00079 *          equilibrated before it is factored.
00080 *          = 'F':  On entry, AFB and IPIV contain the factored form of
00081 *                  A.  If EQUED is not 'N', the matrix A has been
00082 *                  equilibrated with scaling factors given by R and C.
00083 *                  AB, AFB, and IPIV are not modified.
00084 *          = 'N':  The matrix A will be copied to AFB and factored.
00085 *          = 'E':  The matrix A will be equilibrated if necessary, then
00086 *                  copied to AFB and factored.
00087 *
00088 *  TRANS   (input) CHARACTER*1
00089 *          Specifies the form of the system of equations.
00090 *          = 'N':  A * X = B     (No transpose)
00091 *          = 'T':  A**T * X = B  (Transpose)
00092 *          = 'C':  A**H * X = B  (Transpose)
00093 *
00094 *  N       (input) INTEGER
00095 *          The number of linear equations, i.e., the order of the
00096 *          matrix A.  N >= 0.
00097 *
00098 *  KL      (input) INTEGER
00099 *          The number of subdiagonals within the band of A.  KL >= 0.
00100 *
00101 *  KU      (input) INTEGER
00102 *          The number of superdiagonals within the band of A.  KU >= 0.
00103 *
00104 *  NRHS    (input) INTEGER
00105 *          The number of right hand sides, i.e., the number of columns
00106 *          of the matrices B and X.  NRHS >= 0.
00107 *
00108 *  AB      (input/output) DOUBLE PRECISION array, dimension (LDAB,N)
00109 *          On entry, the matrix A in band storage, in rows 1 to KL+KU+1.
00110 *          The j-th column of A is stored in the j-th column of the
00111 *          array AB as follows:
00112 *          AB(KU+1+i-j,j) = A(i,j) for max(1,j-KU)<=i<=min(N,j+kl)
00113 *
00114 *          If FACT = 'F' and EQUED is not 'N', then A must have been
00115 *          equilibrated by the scaling factors in R and/or C.  AB is not
00116 *          modified if FACT = 'F' or 'N', or if FACT = 'E' and
00117 *          EQUED = 'N' on exit.
00118 *
00119 *          On exit, if EQUED .ne. 'N', A is scaled as follows:
00120 *          EQUED = 'R':  A := diag(R) * A
00121 *          EQUED = 'C':  A := A * diag(C)
00122 *          EQUED = 'B':  A := diag(R) * A * diag(C).
00123 *
00124 *  LDAB    (input) INTEGER
00125 *          The leading dimension of the array AB.  LDAB >= KL+KU+1.
00126 *
00127 *  AFB     (input or output) DOUBLE PRECISION array, dimension (LDAFB,N)
00128 *          If FACT = 'F', then AFB is an input argument and on entry
00129 *          contains details of the LU factorization of the band matrix
00130 *          A, as computed by DGBTRF.  U is stored as an upper triangular
00131 *          band matrix with KL+KU superdiagonals in rows 1 to KL+KU+1,
00132 *          and the multipliers used during the factorization are stored
00133 *          in rows KL+KU+2 to 2*KL+KU+1.  If EQUED .ne. 'N', then AFB is
00134 *          the factored form of the equilibrated matrix A.
00135 *
00136 *          If FACT = 'N', then AFB is an output argument and on exit
00137 *          returns details of the LU factorization of A.
00138 *
00139 *          If FACT = 'E', then AFB is an output argument and on exit
00140 *          returns details of the LU factorization of the equilibrated
00141 *          matrix A (see the description of AB for the form of the
00142 *          equilibrated matrix).
00143 *
00144 *  LDAFB   (input) INTEGER
00145 *          The leading dimension of the array AFB.  LDAFB >= 2*KL+KU+1.
00146 *
00147 *  IPIV    (input or output) INTEGER array, dimension (N)
00148 *          If FACT = 'F', then IPIV is an input argument and on entry
00149 *          contains the pivot indices from the factorization A = L*U
00150 *          as computed by DGBTRF; row i of the matrix was interchanged
00151 *          with row IPIV(i).
00152 *
00153 *          If FACT = 'N', then IPIV is an output argument and on exit
00154 *          contains the pivot indices from the factorization A = L*U
00155 *          of the original matrix A.
00156 *
00157 *          If FACT = 'E', then IPIV is an output argument and on exit
00158 *          contains the pivot indices from the factorization A = L*U
00159 *          of the equilibrated matrix A.
00160 *
00161 *  EQUED   (input or output) CHARACTER*1
00162 *          Specifies the form of equilibration that was done.
00163 *          = 'N':  No equilibration (always true if FACT = 'N').
00164 *          = 'R':  Row equilibration, i.e., A has been premultiplied by
00165 *                  diag(R).
00166 *          = 'C':  Column equilibration, i.e., A has been postmultiplied
00167 *                  by diag(C).
00168 *          = 'B':  Both row and column equilibration, i.e., A has been
00169 *                  replaced by diag(R) * A * diag(C).
00170 *          EQUED is an input argument if FACT = 'F'; otherwise, it is an
00171 *          output argument.
00172 *
00173 *  R       (input or output) DOUBLE PRECISION array, dimension (N)
00174 *          The row scale factors for A.  If EQUED = 'R' or 'B', A is
00175 *          multiplied on the left by diag(R); if EQUED = 'N' or 'C', R
00176 *          is not accessed.  R is an input argument if FACT = 'F';
00177 *          otherwise, R is an output argument.  If FACT = 'F' and
00178 *          EQUED = 'R' or 'B', each element of R must be positive.
00179 *
00180 *  C       (input or output) DOUBLE PRECISION array, dimension (N)
00181 *          The column scale factors for A.  If EQUED = 'C' or 'B', A is
00182 *          multiplied on the right by diag(C); if EQUED = 'N' or 'R', C
00183 *          is not accessed.  C is an input argument if FACT = 'F';
00184 *          otherwise, C is an output argument.  If FACT = 'F' and
00185 *          EQUED = 'C' or 'B', each element of C must be positive.
00186 *
00187 *  B       (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS)
00188 *          On entry, the right hand side matrix B.
00189 *          On exit,
00190 *          if EQUED = 'N', B is not modified;
00191 *          if TRANS = 'N' and EQUED = 'R' or 'B', B is overwritten by
00192 *          diag(R)*B;
00193 *          if TRANS = 'T' or 'C' and EQUED = 'C' or 'B', B is
00194 *          overwritten by diag(C)*B.
00195 *
00196 *  LDB     (input) INTEGER
00197 *          The leading dimension of the array B.  LDB >= max(1,N).
00198 *
00199 *  X       (output) DOUBLE PRECISION array, dimension (LDX,NRHS)
00200 *          If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X
00201 *          to the original system of equations.  Note that A and B are
00202 *          modified on exit if EQUED .ne. 'N', and the solution to the
00203 *          equilibrated system is inv(diag(C))*X if TRANS = 'N' and
00204 *          EQUED = 'C' or 'B', or inv(diag(R))*X if TRANS = 'T' or 'C'
00205 *          and EQUED = 'R' or 'B'.
00206 *
00207 *  LDX     (input) INTEGER
00208 *          The leading dimension of the array X.  LDX >= max(1,N).
00209 *
00210 *  RCOND   (output) DOUBLE PRECISION
00211 *          The estimate of the reciprocal condition number of the matrix
00212 *          A after equilibration (if done).  If RCOND is less than the
00213 *          machine precision (in particular, if RCOND = 0), the matrix
00214 *          is singular to working precision.  This condition is
00215 *          indicated by a return code of INFO > 0.
00216 *
00217 *  FERR    (output) DOUBLE PRECISION array, dimension (NRHS)
00218 *          The estimated forward error bound for each solution vector
00219 *          X(j) (the j-th column of the solution matrix X).
00220 *          If XTRUE is the true solution corresponding to X(j), FERR(j)
00221 *          is an estimated upper bound for the magnitude of the largest
00222 *          element in (X(j) - XTRUE) divided by the magnitude of the
00223 *          largest element in X(j).  The estimate is as reliable as
00224 *          the estimate for RCOND, and is almost always a slight
00225 *          overestimate of the true error.
00226 *
00227 *  BERR    (output) DOUBLE PRECISION array, dimension (NRHS)
00228 *          The componentwise relative backward error of each solution
00229 *          vector X(j) (i.e., the smallest relative change in
00230 *          any element of A or B that makes X(j) an exact solution).
00231 *
00232 *  WORK    (workspace/output) DOUBLE PRECISION array, dimension (3*N)
00233 *          On exit, WORK(1) contains the reciprocal pivot growth
00234 *          factor norm(A)/norm(U). The "max absolute element" norm is
00235 *          used. If WORK(1) is much less than 1, then the stability
00236 *          of the LU factorization of the (equilibrated) matrix A
00237 *          could be poor. This also means that the solution X, condition
00238 *          estimator RCOND, and forward error bound FERR could be
00239 *          unreliable. If factorization fails with 0<INFO<=N, then
00240 *          WORK(1) contains the reciprocal pivot growth factor for the
00241 *          leading INFO columns of A.
00242 *
00243 *  IWORK   (workspace) INTEGER array, dimension (N)
00244 *
00245 *  INFO    (output) INTEGER
00246 *          = 0:  successful exit
00247 *          < 0:  if INFO = -i, the i-th argument had an illegal value
00248 *          > 0:  if INFO = i, and i is
00249 *                <= N:  U(i,i) is exactly zero.  The factorization
00250 *                       has been completed, but the factor U is exactly
00251 *                       singular, so the solution and error bounds
00252 *                       could not be computed. RCOND = 0 is returned.
00253 *                = N+1: U is nonsingular, but RCOND is less than machine
00254 *                       precision, meaning that the matrix is singular
00255 *                       to working precision.  Nevertheless, the
00256 *                       solution and error bounds are computed because
00257 *                       there are a number of situations where the
00258 *                       computed solution can be more accurate than the
00259 *                       value of RCOND would suggest.
00260 *
00261 *  =====================================================================
00262 *
00263 *     .. Parameters ..
00264       DOUBLE PRECISION   ZERO, ONE
00265       PARAMETER          ( ZERO = 0.0D+0, ONE = 1.0D+0 )
00266 *     ..
00267 *     .. Local Scalars ..
00268       LOGICAL            COLEQU, EQUIL, NOFACT, NOTRAN, ROWEQU
00269       CHARACTER          NORM
00270       INTEGER            I, INFEQU, J, J1, J2
00271       DOUBLE PRECISION   AMAX, ANORM, BIGNUM, COLCND, RCMAX, RCMIN,
00272      $                   ROWCND, RPVGRW, SMLNUM
00273 *     ..
00274 *     .. External Functions ..
00275       LOGICAL            LSAME
00276       DOUBLE PRECISION   DLAMCH, DLANGB, DLANTB
00277       EXTERNAL           LSAME, DLAMCH, DLANGB, DLANTB
00278 *     ..
00279 *     .. External Subroutines ..
00280       EXTERNAL           DCOPY, DGBCON, DGBEQU, DGBRFS, DGBTRF, DGBTRS,
00281      $                   DLACPY, DLAQGB, XERBLA
00282 *     ..
00283 *     .. Intrinsic Functions ..
00284       INTRINSIC          ABS, MAX, MIN
00285 *     ..
00286 *     .. Executable Statements ..
00287 *
00288       INFO = 0
00289       NOFACT = LSAME( FACT, 'N' )
00290       EQUIL = LSAME( FACT, 'E' )
00291       NOTRAN = LSAME( TRANS, 'N' )
00292       IF( NOFACT .OR. EQUIL ) THEN
00293          EQUED = 'N'
00294          ROWEQU = .FALSE.
00295          COLEQU = .FALSE.
00296       ELSE
00297          ROWEQU = LSAME( EQUED, 'R' ) .OR. LSAME( EQUED, 'B' )
00298          COLEQU = LSAME( EQUED, 'C' ) .OR. LSAME( EQUED, 'B' )
00299          SMLNUM = DLAMCH( 'Safe minimum' )
00300          BIGNUM = ONE / SMLNUM
00301       END IF
00302 *
00303 *     Test the input parameters.
00304 *
00305       IF( .NOT.NOFACT .AND. .NOT.EQUIL .AND. .NOT.LSAME( FACT, 'F' ) )
00306      $     THEN
00307          INFO = -1
00308       ELSE IF( .NOT.NOTRAN .AND. .NOT.LSAME( TRANS, 'T' ) .AND. .NOT.
00309      $         LSAME( TRANS, 'C' ) ) THEN
00310          INFO = -2
00311       ELSE IF( N.LT.0 ) THEN
00312          INFO = -3
00313       ELSE IF( KL.LT.0 ) THEN
00314          INFO = -4
00315       ELSE IF( KU.LT.0 ) THEN
00316          INFO = -5
00317       ELSE IF( NRHS.LT.0 ) THEN
00318          INFO = -6
00319       ELSE IF( LDAB.LT.KL+KU+1 ) THEN
00320          INFO = -8
00321       ELSE IF( LDAFB.LT.2*KL+KU+1 ) THEN
00322          INFO = -10
00323       ELSE IF( LSAME( FACT, 'F' ) .AND. .NOT.
00324      $         ( ROWEQU .OR. COLEQU .OR. LSAME( EQUED, 'N' ) ) ) THEN
00325          INFO = -12
00326       ELSE
00327          IF( ROWEQU ) THEN
00328             RCMIN = BIGNUM
00329             RCMAX = ZERO
00330             DO 10 J = 1, N
00331                RCMIN = MIN( RCMIN, R( J ) )
00332                RCMAX = MAX( RCMAX, R( J ) )
00333    10       CONTINUE
00334             IF( RCMIN.LE.ZERO ) THEN
00335                INFO = -13
00336             ELSE IF( N.GT.0 ) THEN
00337                ROWCND = MAX( RCMIN, SMLNUM ) / MIN( RCMAX, BIGNUM )
00338             ELSE
00339                ROWCND = ONE
00340             END IF
00341          END IF
00342          IF( COLEQU .AND. INFO.EQ.0 ) THEN
00343             RCMIN = BIGNUM
00344             RCMAX = ZERO
00345             DO 20 J = 1, N
00346                RCMIN = MIN( RCMIN, C( J ) )
00347                RCMAX = MAX( RCMAX, C( J ) )
00348    20       CONTINUE
00349             IF( RCMIN.LE.ZERO ) THEN
00350                INFO = -14
00351             ELSE IF( N.GT.0 ) THEN
00352                COLCND = MAX( RCMIN, SMLNUM ) / MIN( RCMAX, BIGNUM )
00353             ELSE
00354                COLCND = ONE
00355             END IF
00356          END IF
00357          IF( INFO.EQ.0 ) THEN
00358             IF( LDB.LT.MAX( 1, N ) ) THEN
00359                INFO = -16
00360             ELSE IF( LDX.LT.MAX( 1, N ) ) THEN
00361                INFO = -18
00362             END IF
00363          END IF
00364       END IF
00365 *
00366       IF( INFO.NE.0 ) THEN
00367          CALL XERBLA( 'DGBSVX', -INFO )
00368          RETURN
00369       END IF
00370 *
00371       IF( EQUIL ) THEN
00372 *
00373 *        Compute row and column scalings to equilibrate the matrix A.
00374 *
00375          CALL DGBEQU( N, N, KL, KU, AB, LDAB, R, C, ROWCND, COLCND,
00376      $                AMAX, INFEQU )
00377          IF( INFEQU.EQ.0 ) THEN
00378 *
00379 *           Equilibrate the matrix.
00380 *
00381             CALL DLAQGB( N, N, KL, KU, AB, LDAB, R, C, ROWCND, COLCND,
00382      $                   AMAX, EQUED )
00383             ROWEQU = LSAME( EQUED, 'R' ) .OR. LSAME( EQUED, 'B' )
00384             COLEQU = LSAME( EQUED, 'C' ) .OR. LSAME( EQUED, 'B' )
00385          END IF
00386       END IF
00387 *
00388 *     Scale the right hand side.
00389 *
00390       IF( NOTRAN ) THEN
00391          IF( ROWEQU ) THEN
00392             DO 40 J = 1, NRHS
00393                DO 30 I = 1, N
00394                   B( I, J ) = R( I )*B( I, J )
00395    30          CONTINUE
00396    40       CONTINUE
00397          END IF
00398       ELSE IF( COLEQU ) THEN
00399          DO 60 J = 1, NRHS
00400             DO 50 I = 1, N
00401                B( I, J ) = C( I )*B( I, J )
00402    50       CONTINUE
00403    60    CONTINUE
00404       END IF
00405 *
00406       IF( NOFACT .OR. EQUIL ) THEN
00407 *
00408 *        Compute the LU factorization of the band matrix A.
00409 *
00410          DO 70 J = 1, N
00411             J1 = MAX( J-KU, 1 )
00412             J2 = MIN( J+KL, N )
00413             CALL DCOPY( J2-J1+1, AB( KU+1-J+J1, J ), 1,
00414      $                  AFB( KL+KU+1-J+J1, J ), 1 )
00415    70    CONTINUE
00416 *
00417          CALL DGBTRF( N, N, KL, KU, AFB, LDAFB, IPIV, INFO )
00418 *
00419 *        Return if INFO is non-zero.
00420 *
00421          IF( INFO.GT.0 ) THEN
00422 *
00423 *           Compute the reciprocal pivot growth factor of the
00424 *           leading rank-deficient INFO columns of A.
00425 *
00426             ANORM = ZERO
00427             DO 90 J = 1, INFO
00428                DO 80 I = MAX( KU+2-J, 1 ), MIN( N+KU+1-J, KL+KU+1 )
00429                   ANORM = MAX( ANORM, ABS( AB( I, J ) ) )
00430    80          CONTINUE
00431    90       CONTINUE
00432             RPVGRW = DLANTB( 'M', 'U', 'N', INFO, MIN( INFO-1, KL+KU ),
00433      $                       AFB( MAX( 1, KL+KU+2-INFO ), 1 ), LDAFB,
00434      $                       WORK )
00435             IF( RPVGRW.EQ.ZERO ) THEN
00436                RPVGRW = ONE
00437             ELSE
00438                RPVGRW = ANORM / RPVGRW
00439             END IF
00440             WORK( 1 ) = RPVGRW
00441             RCOND = ZERO
00442             RETURN
00443          END IF
00444       END IF
00445 *
00446 *     Compute the norm of the matrix A and the
00447 *     reciprocal pivot growth factor RPVGRW.
00448 *
00449       IF( NOTRAN ) THEN
00450          NORM = '1'
00451       ELSE
00452          NORM = 'I'
00453       END IF
00454       ANORM = DLANGB( NORM, N, KL, KU, AB, LDAB, WORK )
00455       RPVGRW = DLANTB( 'M', 'U', 'N', N, KL+KU, AFB, LDAFB, WORK )
00456       IF( RPVGRW.EQ.ZERO ) THEN
00457          RPVGRW = ONE
00458       ELSE
00459          RPVGRW = DLANGB( 'M', N, KL, KU, AB, LDAB, WORK ) / RPVGRW
00460       END IF
00461 *
00462 *     Compute the reciprocal of the condition number of A.
00463 *
00464       CALL DGBCON( NORM, N, KL, KU, AFB, LDAFB, IPIV, ANORM, RCOND,
00465      $             WORK, IWORK, INFO )
00466 *
00467 *     Compute the solution matrix X.
00468 *
00469       CALL DLACPY( 'Full', N, NRHS, B, LDB, X, LDX )
00470       CALL DGBTRS( TRANS, N, KL, KU, NRHS, AFB, LDAFB, IPIV, X, LDX,
00471      $             INFO )
00472 *
00473 *     Use iterative refinement to improve the computed solution and
00474 *     compute error bounds and backward error estimates for it.
00475 *
00476       CALL DGBRFS( TRANS, N, KL, KU, NRHS, AB, LDAB, AFB, LDAFB, IPIV,
00477      $             B, LDB, X, LDX, FERR, BERR, WORK, IWORK, INFO )
00478 *
00479 *     Transform the solution matrix X to a solution of the original
00480 *     system.
00481 *
00482       IF( NOTRAN ) THEN
00483          IF( COLEQU ) THEN
00484             DO 110 J = 1, NRHS
00485                DO 100 I = 1, N
00486                   X( I, J ) = C( I )*X( I, J )
00487   100          CONTINUE
00488   110       CONTINUE
00489             DO 120 J = 1, NRHS
00490                FERR( J ) = FERR( J ) / COLCND
00491   120       CONTINUE
00492          END IF
00493       ELSE IF( ROWEQU ) THEN
00494          DO 140 J = 1, NRHS
00495             DO 130 I = 1, N
00496                X( I, J ) = R( I )*X( I, J )
00497   130       CONTINUE
00498   140    CONTINUE
00499          DO 150 J = 1, NRHS
00500             FERR( J ) = FERR( J ) / ROWCND
00501   150    CONTINUE
00502       END IF
00503 *
00504 *     Set INFO = N+1 if the matrix is singular to working precision.
00505 *
00506       IF( RCOND.LT.DLAMCH( 'Epsilon' ) )
00507      $   INFO = N + 1
00508 *
00509       WORK( 1 ) = RPVGRW
00510       RETURN
00511 *
00512 *     End of DGBSVX
00513 *
00514       END
 All Files Functions