Annotation of rpl/lapack/lapack/dgbsvx.f, revision 1.2

1.1       bertrand    1:       SUBROUTINE DGBSVX( FACT, TRANS, N, KL, KU, NRHS, AB, LDAB, AFB,
                      2:      $                   LDAFB, IPIV, EQUED, R, C, B, LDB, X, LDX,
                      3:      $                   RCOND, FERR, BERR, WORK, IWORK, INFO )
                      4: *
                      5: *  -- LAPACK driver routine (version 3.2) --
                      6: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
                      7: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
                      8: *     November 2006
                      9: *
                     10: *     .. Scalar Arguments ..
                     11:       CHARACTER          EQUED, FACT, TRANS
                     12:       INTEGER            INFO, KL, KU, LDAB, LDAFB, LDB, LDX, N, NRHS
                     13:       DOUBLE PRECISION   RCOND
                     14: *     ..
                     15: *     .. Array Arguments ..
                     16:       INTEGER            IPIV( * ), IWORK( * )
                     17:       DOUBLE PRECISION   AB( LDAB, * ), AFB( LDAFB, * ), B( LDB, * ),
                     18:      $                   BERR( * ), C( * ), FERR( * ), R( * ),
                     19:      $                   WORK( * ), X( LDX, * )
                     20: *     ..
                     21: *
                     22: *  Purpose
                     23: *  =======
                     24: *
                     25: *  DGBSVX uses the LU factorization to compute the solution to a real
                     26: *  system of linear equations A * X = B, A**T * X = B, or A**H * X = B,
                     27: *  where A is a band matrix of order N with KL subdiagonals and KU
                     28: *  superdiagonals, and X and B are N-by-NRHS matrices.
                     29: *
                     30: *  Error bounds on the solution and a condition estimate are also
                     31: *  provided.
                     32: *
                     33: *  Description
                     34: *  ===========
                     35: *
                     36: *  The following steps are performed by this subroutine:
                     37: *
                     38: *  1. If FACT = 'E', real scaling factors are computed to equilibrate
                     39: *     the system:
                     40: *        TRANS = 'N':  diag(R)*A*diag(C)     *inv(diag(C))*X = diag(R)*B
                     41: *        TRANS = 'T': (diag(R)*A*diag(C))**T *inv(diag(R))*X = diag(C)*B
                     42: *        TRANS = 'C': (diag(R)*A*diag(C))**H *inv(diag(R))*X = diag(C)*B
                     43: *     Whether or not the system will be equilibrated depends on the
                     44: *     scaling of the matrix A, but if equilibration is used, A is
                     45: *     overwritten by diag(R)*A*diag(C) and B by diag(R)*B (if TRANS='N')
                     46: *     or diag(C)*B (if TRANS = 'T' or 'C').
                     47: *
                     48: *  2. If FACT = 'N' or 'E', the LU decomposition is used to factor the
                     49: *     matrix A (after equilibration if FACT = 'E') as
                     50: *        A = L * U,
                     51: *     where L is a product of permutation and unit lower triangular
                     52: *     matrices with KL subdiagonals, and U is upper triangular with
                     53: *     KL+KU superdiagonals.
                     54: *
                     55: *  3. If some U(i,i)=0, so that U is exactly singular, then the routine
                     56: *     returns with INFO = i. Otherwise, the factored form of A is used
                     57: *     to estimate the condition number of the matrix A.  If the
                     58: *     reciprocal of the condition number is less than machine precision,
                     59: *     INFO = N+1 is returned as a warning, but the routine still goes on
                     60: *     to solve for X and compute error bounds as described below.
                     61: *
                     62: *  4. The system of equations is solved for X using the factored form
                     63: *     of A.
                     64: *
                     65: *  5. Iterative refinement is applied to improve the computed solution
                     66: *     matrix and calculate error bounds and backward error estimates
                     67: *     for it.
                     68: *
                     69: *  6. If equilibration was used, the matrix X is premultiplied by
                     70: *     diag(C) (if TRANS = 'N') or diag(R) (if TRANS = 'T' or 'C') so
                     71: *     that it solves the original system before equilibration.
                     72: *
                     73: *  Arguments
                     74: *  =========
                     75: *
                     76: *  FACT    (input) CHARACTER*1
                     77: *          Specifies whether or not the factored form of the matrix A is
                     78: *          supplied on entry, and if not, whether the matrix A should be
                     79: *          equilibrated before it is factored.
                     80: *          = 'F':  On entry, AFB and IPIV contain the factored form of
                     81: *                  A.  If EQUED is not 'N', the matrix A has been
                     82: *                  equilibrated with scaling factors given by R and C.
                     83: *                  AB, AFB, and IPIV are not modified.
                     84: *          = 'N':  The matrix A will be copied to AFB and factored.
                     85: *          = 'E':  The matrix A will be equilibrated if necessary, then
                     86: *                  copied to AFB and factored.
                     87: *
                     88: *  TRANS   (input) CHARACTER*1
                     89: *          Specifies the form of the system of equations.
                     90: *          = 'N':  A * X = B     (No transpose)
                     91: *          = 'T':  A**T * X = B  (Transpose)
                     92: *          = 'C':  A**H * X = B  (Transpose)
                     93: *
                     94: *  N       (input) INTEGER
                     95: *          The number of linear equations, i.e., the order of the
                     96: *          matrix A.  N >= 0.
                     97: *
                     98: *  KL      (input) INTEGER
                     99: *          The number of subdiagonals within the band of A.  KL >= 0.
                    100: *
                    101: *  KU      (input) INTEGER
                    102: *          The number of superdiagonals within the band of A.  KU >= 0.
                    103: *
                    104: *  NRHS    (input) INTEGER
                    105: *          The number of right hand sides, i.e., the number of columns
                    106: *          of the matrices B and X.  NRHS >= 0.
                    107: *
                    108: *  AB      (input/output) DOUBLE PRECISION array, dimension (LDAB,N)
                    109: *          On entry, the matrix A in band storage, in rows 1 to KL+KU+1.
                    110: *          The j-th column of A is stored in the j-th column of the
                    111: *          array AB as follows:
                    112: *          AB(KU+1+i-j,j) = A(i,j) for max(1,j-KU)<=i<=min(N,j+kl)
                    113: *
                    114: *          If FACT = 'F' and EQUED is not 'N', then A must have been
                    115: *          equilibrated by the scaling factors in R and/or C.  AB is not
                    116: *          modified if FACT = 'F' or 'N', or if FACT = 'E' and
                    117: *          EQUED = 'N' on exit.
                    118: *
                    119: *          On exit, if EQUED .ne. 'N', A is scaled as follows:
                    120: *          EQUED = 'R':  A := diag(R) * A
                    121: *          EQUED = 'C':  A := A * diag(C)
                    122: *          EQUED = 'B':  A := diag(R) * A * diag(C).
                    123: *
                    124: *  LDAB    (input) INTEGER
                    125: *          The leading dimension of the array AB.  LDAB >= KL+KU+1.
                    126: *
                    127: *  AFB     (input or output) DOUBLE PRECISION array, dimension (LDAFB,N)
                    128: *          If FACT = 'F', then AFB is an input argument and on entry
                    129: *          contains details of the LU factorization of the band matrix
                    130: *          A, as computed by DGBTRF.  U is stored as an upper triangular
                    131: *          band matrix with KL+KU superdiagonals in rows 1 to KL+KU+1,
                    132: *          and the multipliers used during the factorization are stored
                    133: *          in rows KL+KU+2 to 2*KL+KU+1.  If EQUED .ne. 'N', then AFB is
                    134: *          the factored form of the equilibrated matrix A.
                    135: *
                    136: *          If FACT = 'N', then AFB is an output argument and on exit
                    137: *          returns details of the LU factorization of A.
                    138: *
                    139: *          If FACT = 'E', then AFB is an output argument and on exit
                    140: *          returns details of the LU factorization of the equilibrated
                    141: *          matrix A (see the description of AB for the form of the
                    142: *          equilibrated matrix).
                    143: *
                    144: *  LDAFB   (input) INTEGER
                    145: *          The leading dimension of the array AFB.  LDAFB >= 2*KL+KU+1.
                    146: *
                    147: *  IPIV    (input or output) INTEGER array, dimension (N)
                    148: *          If FACT = 'F', then IPIV is an input argument and on entry
                    149: *          contains the pivot indices from the factorization A = L*U
                    150: *          as computed by DGBTRF; row i of the matrix was interchanged
                    151: *          with row IPIV(i).
                    152: *
                    153: *          If FACT = 'N', then IPIV is an output argument and on exit
                    154: *          contains the pivot indices from the factorization A = L*U
                    155: *          of the original matrix A.
                    156: *
                    157: *          If FACT = 'E', then IPIV is an output argument and on exit
                    158: *          contains the pivot indices from the factorization A = L*U
                    159: *          of the equilibrated matrix A.
                    160: *
                    161: *  EQUED   (input or output) CHARACTER*1
                    162: *          Specifies the form of equilibration that was done.
                    163: *          = 'N':  No equilibration (always true if FACT = 'N').
                    164: *          = 'R':  Row equilibration, i.e., A has been premultiplied by
                    165: *                  diag(R).
                    166: *          = 'C':  Column equilibration, i.e., A has been postmultiplied
                    167: *                  by diag(C).
                    168: *          = 'B':  Both row and column equilibration, i.e., A has been
                    169: *                  replaced by diag(R) * A * diag(C).
                    170: *          EQUED is an input argument if FACT = 'F'; otherwise, it is an
                    171: *          output argument.
                    172: *
                    173: *  R       (input or output) DOUBLE PRECISION array, dimension (N)
                    174: *          The row scale factors for A.  If EQUED = 'R' or 'B', A is
                    175: *          multiplied on the left by diag(R); if EQUED = 'N' or 'C', R
                    176: *          is not accessed.  R is an input argument if FACT = 'F';
                    177: *          otherwise, R is an output argument.  If FACT = 'F' and
                    178: *          EQUED = 'R' or 'B', each element of R must be positive.
                    179: *
                    180: *  C       (input or output) DOUBLE PRECISION array, dimension (N)
                    181: *          The column scale factors for A.  If EQUED = 'C' or 'B', A is
                    182: *          multiplied on the right by diag(C); if EQUED = 'N' or 'R', C
                    183: *          is not accessed.  C is an input argument if FACT = 'F';
                    184: *          otherwise, C is an output argument.  If FACT = 'F' and
                    185: *          EQUED = 'C' or 'B', each element of C must be positive.
                    186: *
                    187: *  B       (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS)
                    188: *          On entry, the right hand side matrix B.
                    189: *          On exit,
                    190: *          if EQUED = 'N', B is not modified;
                    191: *          if TRANS = 'N' and EQUED = 'R' or 'B', B is overwritten by
                    192: *          diag(R)*B;
                    193: *          if TRANS = 'T' or 'C' and EQUED = 'C' or 'B', B is
                    194: *          overwritten by diag(C)*B.
                    195: *
                    196: *  LDB     (input) INTEGER
                    197: *          The leading dimension of the array B.  LDB >= max(1,N).
                    198: *
                    199: *  X       (output) DOUBLE PRECISION array, dimension (LDX,NRHS)
                    200: *          If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X
                    201: *          to the original system of equations.  Note that A and B are
                    202: *          modified on exit if EQUED .ne. 'N', and the solution to the
                    203: *          equilibrated system is inv(diag(C))*X if TRANS = 'N' and
                    204: *          EQUED = 'C' or 'B', or inv(diag(R))*X if TRANS = 'T' or 'C'
                    205: *          and EQUED = 'R' or 'B'.
                    206: *
                    207: *  LDX     (input) INTEGER
                    208: *          The leading dimension of the array X.  LDX >= max(1,N).
                    209: *
                    210: *  RCOND   (output) DOUBLE PRECISION
                    211: *          The estimate of the reciprocal condition number of the matrix
                    212: *          A after equilibration (if done).  If RCOND is less than the
                    213: *          machine precision (in particular, if RCOND = 0), the matrix
                    214: *          is singular to working precision.  This condition is
                    215: *          indicated by a return code of INFO > 0.
                    216: *
                    217: *  FERR    (output) DOUBLE PRECISION array, dimension (NRHS)
                    218: *          The estimated forward error bound for each solution vector
                    219: *          X(j) (the j-th column of the solution matrix X).
                    220: *          If XTRUE is the true solution corresponding to X(j), FERR(j)
                    221: *          is an estimated upper bound for the magnitude of the largest
                    222: *          element in (X(j) - XTRUE) divided by the magnitude of the
                    223: *          largest element in X(j).  The estimate is as reliable as
                    224: *          the estimate for RCOND, and is almost always a slight
                    225: *          overestimate of the true error.
                    226: *
                    227: *  BERR    (output) DOUBLE PRECISION array, dimension (NRHS)
                    228: *          The componentwise relative backward error of each solution
                    229: *          vector X(j) (i.e., the smallest relative change in
                    230: *          any element of A or B that makes X(j) an exact solution).
                    231: *
                    232: *  WORK    (workspace/output) DOUBLE PRECISION array, dimension (3*N)
                    233: *          On exit, WORK(1) contains the reciprocal pivot growth
                    234: *          factor norm(A)/norm(U). The "max absolute element" norm is
                    235: *          used. If WORK(1) is much less than 1, then the stability
                    236: *          of the LU factorization of the (equilibrated) matrix A
                    237: *          could be poor. This also means that the solution X, condition
                    238: *          estimator RCOND, and forward error bound FERR could be
                    239: *          unreliable. If factorization fails with 0<INFO<=N, then
                    240: *          WORK(1) contains the reciprocal pivot growth factor for the
                    241: *          leading INFO columns of A.
                    242: *
                    243: *  IWORK   (workspace) INTEGER array, dimension (N)
                    244: *
                    245: *  INFO    (output) INTEGER
                    246: *          = 0:  successful exit
                    247: *          < 0:  if INFO = -i, the i-th argument had an illegal value
                    248: *          > 0:  if INFO = i, and i is
                    249: *                <= N:  U(i,i) is exactly zero.  The factorization
                    250: *                       has been completed, but the factor U is exactly
                    251: *                       singular, so the solution and error bounds
                    252: *                       could not be computed. RCOND = 0 is returned.
                    253: *                = N+1: U is nonsingular, but RCOND is less than machine
                    254: *                       precision, meaning that the matrix is singular
                    255: *                       to working precision.  Nevertheless, the
                    256: *                       solution and error bounds are computed because
                    257: *                       there are a number of situations where the
                    258: *                       computed solution can be more accurate than the
                    259: *                       value of RCOND would suggest.
                    260: *
                    261: *  =====================================================================
                    262: *
                    263: *     .. Parameters ..
                    264:       DOUBLE PRECISION   ZERO, ONE
                    265:       PARAMETER          ( ZERO = 0.0D+0, ONE = 1.0D+0 )
                    266: *     ..
                    267: *     .. Local Scalars ..
                    268:       LOGICAL            COLEQU, EQUIL, NOFACT, NOTRAN, ROWEQU
                    269:       CHARACTER          NORM
                    270:       INTEGER            I, INFEQU, J, J1, J2
                    271:       DOUBLE PRECISION   AMAX, ANORM, BIGNUM, COLCND, RCMAX, RCMIN,
                    272:      $                   ROWCND, RPVGRW, SMLNUM
                    273: *     ..
                    274: *     .. External Functions ..
                    275:       LOGICAL            LSAME
                    276:       DOUBLE PRECISION   DLAMCH, DLANGB, DLANTB
                    277:       EXTERNAL           LSAME, DLAMCH, DLANGB, DLANTB
                    278: *     ..
                    279: *     .. External Subroutines ..
                    280:       EXTERNAL           DCOPY, DGBCON, DGBEQU, DGBRFS, DGBTRF, DGBTRS,
                    281:      $                   DLACPY, DLAQGB, XERBLA
                    282: *     ..
                    283: *     .. Intrinsic Functions ..
                    284:       INTRINSIC          ABS, MAX, MIN
                    285: *     ..
                    286: *     .. Executable Statements ..
                    287: *
                    288:       INFO = 0
                    289:       NOFACT = LSAME( FACT, 'N' )
                    290:       EQUIL = LSAME( FACT, 'E' )
                    291:       NOTRAN = LSAME( TRANS, 'N' )
                    292:       IF( NOFACT .OR. EQUIL ) THEN
                    293:          EQUED = 'N'
                    294:          ROWEQU = .FALSE.
                    295:          COLEQU = .FALSE.
                    296:       ELSE
                    297:          ROWEQU = LSAME( EQUED, 'R' ) .OR. LSAME( EQUED, 'B' )
                    298:          COLEQU = LSAME( EQUED, 'C' ) .OR. LSAME( EQUED, 'B' )
                    299:          SMLNUM = DLAMCH( 'Safe minimum' )
                    300:          BIGNUM = ONE / SMLNUM
                    301:       END IF
                    302: *
                    303: *     Test the input parameters.
                    304: *
                    305:       IF( .NOT.NOFACT .AND. .NOT.EQUIL .AND. .NOT.LSAME( FACT, 'F' ) )
                    306:      $     THEN
                    307:          INFO = -1
                    308:       ELSE IF( .NOT.NOTRAN .AND. .NOT.LSAME( TRANS, 'T' ) .AND. .NOT.
                    309:      $         LSAME( TRANS, 'C' ) ) THEN
                    310:          INFO = -2
                    311:       ELSE IF( N.LT.0 ) THEN
                    312:          INFO = -3
                    313:       ELSE IF( KL.LT.0 ) THEN
                    314:          INFO = -4
                    315:       ELSE IF( KU.LT.0 ) THEN
                    316:          INFO = -5
                    317:       ELSE IF( NRHS.LT.0 ) THEN
                    318:          INFO = -6
                    319:       ELSE IF( LDAB.LT.KL+KU+1 ) THEN
                    320:          INFO = -8
                    321:       ELSE IF( LDAFB.LT.2*KL+KU+1 ) THEN
                    322:          INFO = -10
                    323:       ELSE IF( LSAME( FACT, 'F' ) .AND. .NOT.
                    324:      $         ( ROWEQU .OR. COLEQU .OR. LSAME( EQUED, 'N' ) ) ) THEN
                    325:          INFO = -12
                    326:       ELSE
                    327:          IF( ROWEQU ) THEN
                    328:             RCMIN = BIGNUM
                    329:             RCMAX = ZERO
                    330:             DO 10 J = 1, N
                    331:                RCMIN = MIN( RCMIN, R( J ) )
                    332:                RCMAX = MAX( RCMAX, R( J ) )
                    333:    10       CONTINUE
                    334:             IF( RCMIN.LE.ZERO ) THEN
                    335:                INFO = -13
                    336:             ELSE IF( N.GT.0 ) THEN
                    337:                ROWCND = MAX( RCMIN, SMLNUM ) / MIN( RCMAX, BIGNUM )
                    338:             ELSE
                    339:                ROWCND = ONE
                    340:             END IF
                    341:          END IF
                    342:          IF( COLEQU .AND. INFO.EQ.0 ) THEN
                    343:             RCMIN = BIGNUM
                    344:             RCMAX = ZERO
                    345:             DO 20 J = 1, N
                    346:                RCMIN = MIN( RCMIN, C( J ) )
                    347:                RCMAX = MAX( RCMAX, C( J ) )
                    348:    20       CONTINUE
                    349:             IF( RCMIN.LE.ZERO ) THEN
                    350:                INFO = -14
                    351:             ELSE IF( N.GT.0 ) THEN
                    352:                COLCND = MAX( RCMIN, SMLNUM ) / MIN( RCMAX, BIGNUM )
                    353:             ELSE
                    354:                COLCND = ONE
                    355:             END IF
                    356:          END IF
                    357:          IF( INFO.EQ.0 ) THEN
                    358:             IF( LDB.LT.MAX( 1, N ) ) THEN
                    359:                INFO = -16
                    360:             ELSE IF( LDX.LT.MAX( 1, N ) ) THEN
                    361:                INFO = -18
                    362:             END IF
                    363:          END IF
                    364:       END IF
                    365: *
                    366:       IF( INFO.NE.0 ) THEN
                    367:          CALL XERBLA( 'DGBSVX', -INFO )
                    368:          RETURN
                    369:       END IF
                    370: *
                    371:       IF( EQUIL ) THEN
                    372: *
                    373: *        Compute row and column scalings to equilibrate the matrix A.
                    374: *
                    375:          CALL DGBEQU( N, N, KL, KU, AB, LDAB, R, C, ROWCND, COLCND,
                    376:      $                AMAX, INFEQU )
                    377:          IF( INFEQU.EQ.0 ) THEN
                    378: *
                    379: *           Equilibrate the matrix.
                    380: *
                    381:             CALL DLAQGB( N, N, KL, KU, AB, LDAB, R, C, ROWCND, COLCND,
                    382:      $                   AMAX, EQUED )
                    383:             ROWEQU = LSAME( EQUED, 'R' ) .OR. LSAME( EQUED, 'B' )
                    384:             COLEQU = LSAME( EQUED, 'C' ) .OR. LSAME( EQUED, 'B' )
                    385:          END IF
                    386:       END IF
                    387: *
                    388: *     Scale the right hand side.
                    389: *
                    390:       IF( NOTRAN ) THEN
                    391:          IF( ROWEQU ) THEN
                    392:             DO 40 J = 1, NRHS
                    393:                DO 30 I = 1, N
                    394:                   B( I, J ) = R( I )*B( I, J )
                    395:    30          CONTINUE
                    396:    40       CONTINUE
                    397:          END IF
                    398:       ELSE IF( COLEQU ) THEN
                    399:          DO 60 J = 1, NRHS
                    400:             DO 50 I = 1, N
                    401:                B( I, J ) = C( I )*B( I, J )
                    402:    50       CONTINUE
                    403:    60    CONTINUE
                    404:       END IF
                    405: *
                    406:       IF( NOFACT .OR. EQUIL ) THEN
                    407: *
                    408: *        Compute the LU factorization of the band matrix A.
                    409: *
                    410:          DO 70 J = 1, N
                    411:             J1 = MAX( J-KU, 1 )
                    412:             J2 = MIN( J+KL, N )
                    413:             CALL DCOPY( J2-J1+1, AB( KU+1-J+J1, J ), 1,
                    414:      $                  AFB( KL+KU+1-J+J1, J ), 1 )
                    415:    70    CONTINUE
                    416: *
                    417:          CALL DGBTRF( N, N, KL, KU, AFB, LDAFB, IPIV, INFO )
                    418: *
                    419: *        Return if INFO is non-zero.
                    420: *
                    421:          IF( INFO.GT.0 ) THEN
                    422: *
                    423: *           Compute the reciprocal pivot growth factor of the
                    424: *           leading rank-deficient INFO columns of A.
                    425: *
                    426:             ANORM = ZERO
                    427:             DO 90 J = 1, INFO
                    428:                DO 80 I = MAX( KU+2-J, 1 ), MIN( N+KU+1-J, KL+KU+1 )
                    429:                   ANORM = MAX( ANORM, ABS( AB( I, J ) ) )
                    430:    80          CONTINUE
                    431:    90       CONTINUE
                    432:             RPVGRW = DLANTB( 'M', 'U', 'N', INFO, MIN( INFO-1, KL+KU ),
                    433:      $                       AFB( MAX( 1, KL+KU+2-INFO ), 1 ), LDAFB,
                    434:      $                       WORK )
                    435:             IF( RPVGRW.EQ.ZERO ) THEN
                    436:                RPVGRW = ONE
                    437:             ELSE
                    438:                RPVGRW = ANORM / RPVGRW
                    439:             END IF
                    440:             WORK( 1 ) = RPVGRW
                    441:             RCOND = ZERO
                    442:             RETURN
                    443:          END IF
                    444:       END IF
                    445: *
                    446: *     Compute the norm of the matrix A and the
                    447: *     reciprocal pivot growth factor RPVGRW.
                    448: *
                    449:       IF( NOTRAN ) THEN
                    450:          NORM = '1'
                    451:       ELSE
                    452:          NORM = 'I'
                    453:       END IF
                    454:       ANORM = DLANGB( NORM, N, KL, KU, AB, LDAB, WORK )
                    455:       RPVGRW = DLANTB( 'M', 'U', 'N', N, KL+KU, AFB, LDAFB, WORK )
                    456:       IF( RPVGRW.EQ.ZERO ) THEN
                    457:          RPVGRW = ONE
                    458:       ELSE
                    459:          RPVGRW = DLANGB( 'M', N, KL, KU, AB, LDAB, WORK ) / RPVGRW
                    460:       END IF
                    461: *
                    462: *     Compute the reciprocal of the condition number of A.
                    463: *
                    464:       CALL DGBCON( NORM, N, KL, KU, AFB, LDAFB, IPIV, ANORM, RCOND,
                    465:      $             WORK, IWORK, INFO )
                    466: *
                    467: *     Compute the solution matrix X.
                    468: *
                    469:       CALL DLACPY( 'Full', N, NRHS, B, LDB, X, LDX )
                    470:       CALL DGBTRS( TRANS, N, KL, KU, NRHS, AFB, LDAFB, IPIV, X, LDX,
                    471:      $             INFO )
                    472: *
                    473: *     Use iterative refinement to improve the computed solution and
                    474: *     compute error bounds and backward error estimates for it.
                    475: *
                    476:       CALL DGBRFS( TRANS, N, KL, KU, NRHS, AB, LDAB, AFB, LDAFB, IPIV,
                    477:      $             B, LDB, X, LDX, FERR, BERR, WORK, IWORK, INFO )
                    478: *
                    479: *     Transform the solution matrix X to a solution of the original
                    480: *     system.
                    481: *
                    482:       IF( NOTRAN ) THEN
                    483:          IF( COLEQU ) THEN
                    484:             DO 110 J = 1, NRHS
                    485:                DO 100 I = 1, N
                    486:                   X( I, J ) = C( I )*X( I, J )
                    487:   100          CONTINUE
                    488:   110       CONTINUE
                    489:             DO 120 J = 1, NRHS
                    490:                FERR( J ) = FERR( J ) / COLCND
                    491:   120       CONTINUE
                    492:          END IF
                    493:       ELSE IF( ROWEQU ) THEN
                    494:          DO 140 J = 1, NRHS
                    495:             DO 130 I = 1, N
                    496:                X( I, J ) = R( I )*X( I, J )
                    497:   130       CONTINUE
                    498:   140    CONTINUE
                    499:          DO 150 J = 1, NRHS
                    500:             FERR( J ) = FERR( J ) / ROWCND
                    501:   150    CONTINUE
                    502:       END IF
                    503: *
                    504: *     Set INFO = N+1 if the matrix is singular to working precision.
                    505: *
                    506:       IF( RCOND.LT.DLAMCH( 'Epsilon' ) )
                    507:      $   INFO = N + 1
                    508: *
                    509:       WORK( 1 ) = RPVGRW
                    510:       RETURN
                    511: *
                    512: *     End of DGBSVX
                    513: *
                    514:       END

CVSweb interface <joel.bertrand@systella.fr>