LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ dgesvxx()

subroutine dgesvxx ( character  fact,
character  trans,
integer  n,
integer  nrhs,
double precision, dimension( lda, * )  a,
integer  lda,
double precision, dimension( ldaf, * )  af,
integer  ldaf,
integer, dimension( * )  ipiv,
character  equed,
double precision, dimension( * )  r,
double precision, dimension( * )  c,
double precision, dimension( ldb, * )  b,
integer  ldb,
double precision, dimension( ldx , * )  x,
integer  ldx,
double precision  rcond,
double precision  rpvgrw,
double precision, dimension( * )  berr,
integer  n_err_bnds,
double precision, dimension( nrhs, * )  err_bnds_norm,
double precision, dimension( nrhs, * )  err_bnds_comp,
integer  nparams,
double precision, dimension( * )  params,
double precision, dimension( * )  work,
integer, dimension( * )  iwork,
integer  info 
)

DGESVXX computes the solution to system of linear equations A * X = B for GE matrices

Download DGESVXX + dependencies [TGZ] [ZIP] [TXT]

Purpose:
    DGESVXX uses the LU factorization to compute the solution to a
    double precision system of linear equations  A * X = B,  where A is an
    N-by-N matrix and X and B are N-by-NRHS matrices.

    If requested, both normwise and maximum componentwise error bounds
    are returned. DGESVXX will return a solution with a tiny
    guaranteed error (O(eps) where eps is the working machine
    precision) unless the matrix is very ill-conditioned, in which
    case a warning is returned. Relevant condition numbers also are
    calculated and returned.

    DGESVXX accepts user-provided factorizations and equilibration
    factors; see the definitions of the FACT and EQUED options.
    Solving with refinement and using a factorization from a previous
    DGESVXX call will also produce a solution with either O(eps)
    errors or warnings, but we cannot make that claim for general
    user-provided factorizations and equilibration factors if they
    differ from what DGESVXX would itself produce.
Description:
    The following steps are performed:

    1. If FACT = 'E', double precision scaling factors are computed to equilibrate
    the system:

      TRANS = 'N':  diag(R)*A*diag(C)     *inv(diag(C))*X = diag(R)*B
      TRANS = 'T': (diag(R)*A*diag(C))**T *inv(diag(R))*X = diag(C)*B
      TRANS = 'C': (diag(R)*A*diag(C))**H *inv(diag(R))*X = diag(C)*B

    Whether or not the system will be equilibrated depends on the
    scaling of the matrix A, but if equilibration is used, A is
    overwritten by diag(R)*A*diag(C) and B by diag(R)*B (if TRANS='N')
    or diag(C)*B (if TRANS = 'T' or 'C').

    2. If FACT = 'N' or 'E', the LU decomposition is used to factor
    the matrix A (after equilibration if FACT = 'E') as

      A = P * L * U,

    where P is a permutation matrix, L is a unit lower triangular
    matrix, and U is upper triangular.

    3. If some U(i,i)=0, so that U is exactly singular, then the
    routine returns with INFO = i. Otherwise, the factored form of A
    is used to estimate the condition number of the matrix A (see
    argument RCOND). If the reciprocal of the condition number is less
    than machine precision, the routine still goes on to solve for X
    and compute error bounds as described below.

    4. The system of equations is solved for X using the factored form
    of A.

    5. By default (unless PARAMS(LA_LINRX_ITREF_I) is set to zero),
    the routine will use iterative refinement to try to get a small
    error and error bounds.  Refinement calculates the residual to at
    least twice the working precision.

    6. If equilibration was used, the matrix X is premultiplied by
    diag(C) (if TRANS = 'N') or diag(R) (if TRANS = 'T' or 'C') so
    that it solves the original system before equilibration.
     Some optional parameters are bundled in the PARAMS array.  These
     settings determine how refinement is performed, but often the
     defaults are acceptable.  If the defaults are acceptable, users
     can pass NPARAMS = 0 which prevents the source code from accessing
     the PARAMS argument.
Parameters
[in]FACT
          FACT is CHARACTER*1
     Specifies whether or not the factored form of the matrix A is
     supplied on entry, and if not, whether the matrix A should be
     equilibrated before it is factored.
       = 'F':  On entry, AF and IPIV contain the factored form of A.
               If EQUED is not 'N', the matrix A has been
               equilibrated with scaling factors given by R and C.
               A, AF, and IPIV are not modified.
       = 'N':  The matrix A will be copied to AF and factored.
       = 'E':  The matrix A will be equilibrated if necessary, then
               copied to AF and factored.
[in]TRANS
          TRANS is CHARACTER*1
     Specifies the form of the system of equations:
       = 'N':  A * X = B     (No transpose)
       = 'T':  A**T * X = B  (Transpose)
       = 'C':  A**H * X = B  (Conjugate Transpose = Transpose)
[in]N
          N is INTEGER
     The number of linear equations, i.e., the order of the
     matrix A.  N >= 0.
[in]NRHS
          NRHS is INTEGER
     The number of right hand sides, i.e., the number of columns
     of the matrices B and X.  NRHS >= 0.
[in,out]A
          A is DOUBLE PRECISION array, dimension (LDA,N)
     On entry, the N-by-N matrix A.  If FACT = 'F' and EQUED is
     not 'N', then A must have been equilibrated by the scaling
     factors in R and/or C.  A is not modified if FACT = 'F' or
     'N', or if FACT = 'E' and EQUED = 'N' on exit.

     On exit, if EQUED .ne. 'N', A is scaled as follows:
     EQUED = 'R':  A := diag(R) * A
     EQUED = 'C':  A := A * diag(C)
     EQUED = 'B':  A := diag(R) * A * diag(C).
[in]LDA
          LDA is INTEGER
     The leading dimension of the array A.  LDA >= max(1,N).
[in,out]AF
          AF is DOUBLE PRECISION array, dimension (LDAF,N)
     If FACT = 'F', then AF is an input argument and on entry
     contains the factors L and U from the factorization
     A = P*L*U as computed by DGETRF.  If EQUED .ne. 'N', then
     AF is the factored form of the equilibrated matrix A.

     If FACT = 'N', then AF is an output argument and on exit
     returns the factors L and U from the factorization A = P*L*U
     of the original matrix A.

     If FACT = 'E', then AF is an output argument and on exit
     returns the factors L and U from the factorization A = P*L*U
     of the equilibrated matrix A (see the description of A for
     the form of the equilibrated matrix).
[in]LDAF
          LDAF is INTEGER
     The leading dimension of the array AF.  LDAF >= max(1,N).
[in,out]IPIV
          IPIV is INTEGER array, dimension (N)
     If FACT = 'F', then IPIV is an input argument and on entry
     contains the pivot indices from the factorization A = P*L*U
     as computed by DGETRF; row i of the matrix was interchanged
     with row IPIV(i).

     If FACT = 'N', then IPIV is an output argument and on exit
     contains the pivot indices from the factorization A = P*L*U
     of the original matrix A.

     If FACT = 'E', then IPIV is an output argument and on exit
     contains the pivot indices from the factorization A = P*L*U
     of the equilibrated matrix A.
[in,out]EQUED
          EQUED is CHARACTER*1
     Specifies the form of equilibration that was done.
       = 'N':  No equilibration (always true if FACT = 'N').
       = 'R':  Row equilibration, i.e., A has been premultiplied by
               diag(R).
       = 'C':  Column equilibration, i.e., A has been postmultiplied
               by diag(C).
       = 'B':  Both row and column equilibration, i.e., A has been
               replaced by diag(R) * A * diag(C).
     EQUED is an input argument if FACT = 'F'; otherwise, it is an
     output argument.
[in,out]R
          R is DOUBLE PRECISION array, dimension (N)
     The row scale factors for A.  If EQUED = 'R' or 'B', A is
     multiplied on the left by diag(R); if EQUED = 'N' or 'C', R
     is not accessed.  R is an input argument if FACT = 'F';
     otherwise, R is an output argument.  If FACT = 'F' and
     EQUED = 'R' or 'B', each element of R must be positive.
     If R is output, each element of R is a power of the radix.
     If R is input, each element of R should be a power of the radix
     to ensure a reliable solution and error estimates. Scaling by
     powers of the radix does not cause rounding errors unless the
     result underflows or overflows. Rounding errors during scaling
     lead to refining with a matrix that is not equivalent to the
     input matrix, producing error estimates that may not be
     reliable.
[in,out]C
          C is DOUBLE PRECISION array, dimension (N)
     The column scale factors for A.  If EQUED = 'C' or 'B', A is
     multiplied on the right by diag(C); if EQUED = 'N' or 'R', C
     is not accessed.  C is an input argument if FACT = 'F';
     otherwise, C is an output argument.  If FACT = 'F' and
     EQUED = 'C' or 'B', each element of C must be positive.
     If C is output, each element of C is a power of the radix.
     If C is input, each element of C should be a power of the radix
     to ensure a reliable solution and error estimates. Scaling by
     powers of the radix does not cause rounding errors unless the
     result underflows or overflows. Rounding errors during scaling
     lead to refining with a matrix that is not equivalent to the
     input matrix, producing error estimates that may not be
     reliable.
[in,out]B
          B is DOUBLE PRECISION array, dimension (LDB,NRHS)
     On entry, the N-by-NRHS right hand side matrix B.
     On exit,
     if EQUED = 'N', B is not modified;
     if TRANS = 'N' and EQUED = 'R' or 'B', B is overwritten by
        diag(R)*B;
     if TRANS = 'T' or 'C' and EQUED = 'C' or 'B', B is
        overwritten by diag(C)*B.
[in]LDB
          LDB is INTEGER
     The leading dimension of the array B.  LDB >= max(1,N).
[out]X
          X is DOUBLE PRECISION array, dimension (LDX,NRHS)
     If INFO = 0, the N-by-NRHS solution matrix X to the original
     system of equations.  Note that A and B are modified on exit
     if EQUED .ne. 'N', and the solution to the equilibrated system is
     inv(diag(C))*X if TRANS = 'N' and EQUED = 'C' or 'B', or
     inv(diag(R))*X if TRANS = 'T' or 'C' and EQUED = 'R' or 'B'.
[in]LDX
          LDX is INTEGER
     The leading dimension of the array X.  LDX >= max(1,N).
[out]RCOND
          RCOND is DOUBLE PRECISION
     Reciprocal scaled condition number.  This is an estimate of the
     reciprocal Skeel condition number of the matrix A after
     equilibration (if done).  If this is less than the machine
     precision (in particular, if it is zero), the matrix is singular
     to working precision.  Note that the error may still be small even
     if this number is very small and the matrix appears ill-
     conditioned.
[out]RPVGRW
          RPVGRW is DOUBLE PRECISION
     Reciprocal pivot growth.  On exit, this contains the reciprocal
     pivot growth factor norm(A)/norm(U). The "max absolute element"
     norm is used.  If this is much less than 1, then the stability of
     the LU factorization of the (equilibrated) matrix A could be poor.
     This also means that the solution X, estimated condition numbers,
     and error bounds could be unreliable. If factorization fails with
     0<INFO<=N, then this contains the reciprocal pivot growth factor
     for the leading INFO columns of A.  In DGESVX, this quantity is
     returned in WORK(1).
[out]BERR
          BERR is DOUBLE PRECISION array, dimension (NRHS)
     Componentwise relative backward error.  This is the
     componentwise relative backward error of each solution vector X(j)
     (i.e., the smallest relative change in any element of A or B that
     makes X(j) an exact solution).
[in]N_ERR_BNDS
          N_ERR_BNDS is INTEGER
     Number of error bounds to return for each right hand side
     and each type (normwise or componentwise).  See ERR_BNDS_NORM and
     ERR_BNDS_COMP below.
[out]ERR_BNDS_NORM
          ERR_BNDS_NORM is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
     For each right-hand side, this array contains information about
     various error bounds and condition numbers corresponding to the
     normwise relative error, which is defined as follows:

     Normwise relative error in the ith solution vector:
             max_j (abs(XTRUE(j,i) - X(j,i)))
            ------------------------------
                  max_j abs(X(j,i))

     The array is indexed by the type of error information as described
     below. There currently are up to three pieces of information
     returned.

     The first index in ERR_BNDS_NORM(i,:) corresponds to the ith
     right-hand side.

     The second index in ERR_BNDS_NORM(:,err) contains the following
     three fields:
     err = 1 "Trust/don't trust" boolean. Trust the answer if the
              reciprocal condition number is less than the threshold
              sqrt(n) * dlamch('Epsilon').

     err = 2 "Guaranteed" error bound: The estimated forward error,
              almost certainly within a factor of 10 of the true error
              so long as the next entry is greater than the threshold
              sqrt(n) * dlamch('Epsilon'). This error bound should only
              be trusted if the previous boolean is true.

     err = 3  Reciprocal condition number: Estimated normwise
              reciprocal condition number.  Compared with the threshold
              sqrt(n) * dlamch('Epsilon') to determine if the error
              estimate is "guaranteed". These reciprocal condition
              numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
              appropriately scaled matrix Z.
              Let Z = S*A, where S scales each row by a power of the
              radix so all absolute row sums of Z are approximately 1.

     See Lapack Working Note 165 for further details and extra
     cautions.
[out]ERR_BNDS_COMP
          ERR_BNDS_COMP is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
     For each right-hand side, this array contains information about
     various error bounds and condition numbers corresponding to the
     componentwise relative error, which is defined as follows:

     Componentwise relative error in the ith solution vector:
                    abs(XTRUE(j,i) - X(j,i))
             max_j ----------------------
                         abs(X(j,i))

     The array is indexed by the right-hand side i (on which the
     componentwise relative error depends), and the type of error
     information as described below. There currently are up to three
     pieces of information returned for each right-hand side. If
     componentwise accuracy is not requested (PARAMS(3) = 0.0), then
     ERR_BNDS_COMP is not accessed.  If N_ERR_BNDS < 3, then at most
     the first (:,N_ERR_BNDS) entries are returned.

     The first index in ERR_BNDS_COMP(i,:) corresponds to the ith
     right-hand side.

     The second index in ERR_BNDS_COMP(:,err) contains the following
     three fields:
     err = 1 "Trust/don't trust" boolean. Trust the answer if the
              reciprocal condition number is less than the threshold
              sqrt(n) * dlamch('Epsilon').

     err = 2 "Guaranteed" error bound: The estimated forward error,
              almost certainly within a factor of 10 of the true error
              so long as the next entry is greater than the threshold
              sqrt(n) * dlamch('Epsilon'). This error bound should only
              be trusted if the previous boolean is true.

     err = 3  Reciprocal condition number: Estimated componentwise
              reciprocal condition number.  Compared with the threshold
              sqrt(n) * dlamch('Epsilon') to determine if the error
              estimate is "guaranteed". These reciprocal condition
              numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
              appropriately scaled matrix Z.
              Let Z = S*(A*diag(x)), where x is the solution for the
              current right-hand side and S scales each row of
              A*diag(x) by a power of the radix so all absolute row
              sums of Z are approximately 1.

     See Lapack Working Note 165 for further details and extra
     cautions.
[in]NPARAMS
          NPARAMS is INTEGER
     Specifies the number of parameters set in PARAMS.  If <= 0, the
     PARAMS array is never referenced and default values are used.
[in,out]PARAMS
          PARAMS is DOUBLE PRECISION array, dimension (NPARAMS)
     Specifies algorithm parameters.  If an entry is < 0.0, then
     that entry will be filled with default value used for that
     parameter.  Only positions up to NPARAMS are accessed; defaults
     are used for higher-numbered parameters.

       PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative
            refinement or not.
         Default: 1.0D+0
            = 0.0:  No refinement is performed, and no error bounds are
                    computed.
            = 1.0:  Use the extra-precise refinement algorithm.
              (other values are reserved for future use)

       PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual
            computations allowed for refinement.
         Default: 10
         Aggressive: Set to 100 to permit convergence using approximate
                     factorizations or factorizations other than LU. If
                     the factorization uses a technique other than
                     Gaussian elimination, the guarantees in
                     err_bnds_norm and err_bnds_comp may no longer be
                     trustworthy.

       PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code
            will attempt to find a solution with small componentwise
            relative error in the double-precision algorithm.  Positive
            is true, 0.0 is false.
         Default: 1.0 (attempt componentwise convergence)
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (4*N)
[out]IWORK
          IWORK is INTEGER array, dimension (N)
[out]INFO
          INFO is INTEGER
       = 0:  Successful exit. The solution to every right-hand side is
         guaranteed.
       < 0:  If INFO = -i, the i-th argument had an illegal value
       > 0 and <= N:  U(INFO,INFO) is exactly zero.  The factorization
         has been completed, but the factor U is exactly singular, so
         the solution and error bounds could not be computed. RCOND = 0
         is returned.
       = N+J: The solution corresponding to the Jth right-hand side is
         not guaranteed. The solutions corresponding to other right-
         hand sides K with K > J may not be guaranteed as well, but
         only the first such right-hand side is reported. If a small
         componentwise error is not requested (PARAMS(3) = 0.0) then
         the Jth right-hand side is the first with a normwise error
         bound that is not guaranteed (the smallest J such
         that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0)
         the Jth right-hand side is the first with either a normwise or
         componentwise error bound that is not guaranteed (the smallest
         J such that either ERR_BNDS_NORM(J,1) = 0.0 or
         ERR_BNDS_COMP(J,1) = 0.0). See the definition of
         ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information
         about all of the right-hand sides check ERR_BNDS_NORM or
         ERR_BNDS_COMP.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 535 of file dgesvxx.f.

540*
541* -- LAPACK driver routine --
542* -- LAPACK is a software package provided by Univ. of Tennessee, --
543* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
544*
545* .. Scalar Arguments ..
546 CHARACTER EQUED, FACT, TRANS
547 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
548 $ N_ERR_BNDS
549 DOUBLE PRECISION RCOND, RPVGRW
550* ..
551* .. Array Arguments ..
552 INTEGER IPIV( * ), IWORK( * )
553 DOUBLE PRECISION A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
554 $ X( LDX , * ),WORK( * )
555 DOUBLE PRECISION R( * ), C( * ), PARAMS( * ), BERR( * ),
556 $ ERR_BNDS_NORM( NRHS, * ),
557 $ ERR_BNDS_COMP( NRHS, * )
558* ..
559*
560* =====================================================================
561*
562* .. Parameters ..
563 DOUBLE PRECISION ZERO, ONE
564 parameter( zero = 0.0d+0, one = 1.0d+0 )
565 INTEGER FINAL_NRM_ERR_I, FINAL_CMP_ERR_I, BERR_I
566 INTEGER RCOND_I, NRM_RCOND_I, NRM_ERR_I, CMP_RCOND_I
567 INTEGER CMP_ERR_I, PIV_GROWTH_I
568 parameter( final_nrm_err_i = 1, final_cmp_err_i = 2,
569 $ berr_i = 3 )
570 parameter( rcond_i = 4, nrm_rcond_i = 5, nrm_err_i = 6 )
571 parameter( cmp_rcond_i = 7, cmp_err_i = 8,
572 $ piv_growth_i = 9 )
573* ..
574* .. Local Scalars ..
575 LOGICAL COLEQU, EQUIL, NOFACT, NOTRAN, ROWEQU
576 INTEGER INFEQU, J
577 DOUBLE PRECISION AMAX, BIGNUM, COLCND, RCMAX, RCMIN, ROWCND,
578 $ SMLNUM
579* ..
580* .. External Functions ..
581 EXTERNAL lsame, dlamch, dla_gerpvgrw
582 LOGICAL LSAME
583 DOUBLE PRECISION DLAMCH, DLA_GERPVGRW
584* ..
585* .. External Subroutines ..
586 EXTERNAL dgeequb, dgetrf, dgetrs, dlacpy, dlaqge,
588* ..
589* .. Intrinsic Functions ..
590 INTRINSIC max, min
591* ..
592* .. Executable Statements ..
593*
594 info = 0
595 nofact = lsame( fact, 'N' )
596 equil = lsame( fact, 'E' )
597 notran = lsame( trans, 'N' )
598 smlnum = dlamch( 'Safe minimum' )
599 bignum = one / smlnum
600 IF( nofact .OR. equil ) THEN
601 equed = 'N'
602 rowequ = .false.
603 colequ = .false.
604 ELSE
605 rowequ = lsame( equed, 'R' ) .OR. lsame( equed, 'B' )
606 colequ = lsame( equed, 'C' ) .OR. lsame( equed, 'B' )
607 END IF
608*
609* Default is failure. If an input parameter is wrong or
610* factorization fails, make everything look horrible. Only the
611* pivot growth is set here, the rest is initialized in DGERFSX.
612*
613 rpvgrw = zero
614*
615* Test the input parameters. PARAMS is not tested until DGERFSX.
616*
617 IF( .NOT.nofact .AND. .NOT.equil .AND. .NOT.
618 $ lsame( fact, 'F' ) ) THEN
619 info = -1
620 ELSE IF( .NOT.notran .AND. .NOT.lsame( trans, 'T' ) .AND. .NOT.
621 $ lsame( trans, 'C' ) ) THEN
622 info = -2
623 ELSE IF( n.LT.0 ) THEN
624 info = -3
625 ELSE IF( nrhs.LT.0 ) THEN
626 info = -4
627 ELSE IF( lda.LT.max( 1, n ) ) THEN
628 info = -6
629 ELSE IF( ldaf.LT.max( 1, n ) ) THEN
630 info = -8
631 ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
632 $ ( rowequ .OR. colequ .OR. lsame( equed, 'N' ) ) ) THEN
633 info = -10
634 ELSE
635 IF( rowequ ) THEN
636 rcmin = bignum
637 rcmax = zero
638 DO 10 j = 1, n
639 rcmin = min( rcmin, r( j ) )
640 rcmax = max( rcmax, r( j ) )
641 10 CONTINUE
642 IF( rcmin.LE.zero ) THEN
643 info = -11
644 ELSE IF( n.GT.0 ) THEN
645 rowcnd = max( rcmin, smlnum ) / min( rcmax, bignum )
646 ELSE
647 rowcnd = one
648 END IF
649 END IF
650 IF( colequ .AND. info.EQ.0 ) THEN
651 rcmin = bignum
652 rcmax = zero
653 DO 20 j = 1, n
654 rcmin = min( rcmin, c( j ) )
655 rcmax = max( rcmax, c( j ) )
656 20 CONTINUE
657 IF( rcmin.LE.zero ) THEN
658 info = -12
659 ELSE IF( n.GT.0 ) THEN
660 colcnd = max( rcmin, smlnum ) / min( rcmax, bignum )
661 ELSE
662 colcnd = one
663 END IF
664 END IF
665 IF( info.EQ.0 ) THEN
666 IF( ldb.LT.max( 1, n ) ) THEN
667 info = -14
668 ELSE IF( ldx.LT.max( 1, n ) ) THEN
669 info = -16
670 END IF
671 END IF
672 END IF
673*
674 IF( info.NE.0 ) THEN
675 CALL xerbla( 'DGESVXX', -info )
676 RETURN
677 END IF
678*
679 IF( equil ) THEN
680*
681* Compute row and column scalings to equilibrate the matrix A.
682*
683 CALL dgeequb( n, n, a, lda, r, c, rowcnd, colcnd, amax,
684 $ infequ )
685 IF( infequ.EQ.0 ) THEN
686*
687* Equilibrate the matrix.
688*
689 CALL dlaqge( n, n, a, lda, r, c, rowcnd, colcnd, amax,
690 $ equed )
691 rowequ = lsame( equed, 'R' ) .OR. lsame( equed, 'B' )
692 colequ = lsame( equed, 'C' ) .OR. lsame( equed, 'B' )
693 END IF
694*
695* If the scaling factors are not applied, set them to 1.0.
696*
697 IF ( .NOT.rowequ ) THEN
698 DO j = 1, n
699 r( j ) = 1.0d+0
700 END DO
701 END IF
702 IF ( .NOT.colequ ) THEN
703 DO j = 1, n
704 c( j ) = 1.0d+0
705 END DO
706 END IF
707 END IF
708*
709* Scale the right-hand side.
710*
711 IF( notran ) THEN
712 IF( rowequ ) CALL dlascl2( n, nrhs, r, b, ldb )
713 ELSE
714 IF( colequ ) CALL dlascl2( n, nrhs, c, b, ldb )
715 END IF
716*
717 IF( nofact .OR. equil ) THEN
718*
719* Compute the LU factorization of A.
720*
721 CALL dlacpy( 'Full', n, n, a, lda, af, ldaf )
722 CALL dgetrf( n, n, af, ldaf, ipiv, info )
723*
724* Return if INFO is non-zero.
725*
726 IF( info.GT.0 ) THEN
727*
728* Pivot in column INFO is exactly 0
729* Compute the reciprocal pivot growth factor of the
730* leading rank-deficient INFO columns of A.
731*
732 rpvgrw = dla_gerpvgrw( n, info, a, lda, af, ldaf )
733 RETURN
734 END IF
735 END IF
736*
737* Compute the reciprocal pivot growth factor RPVGRW.
738*
739 rpvgrw = dla_gerpvgrw( n, n, a, lda, af, ldaf )
740*
741* Compute the solution matrix X.
742*
743 CALL dlacpy( 'Full', n, nrhs, b, ldb, x, ldx )
744 CALL dgetrs( trans, n, nrhs, af, ldaf, ipiv, x, ldx, info )
745*
746* Use iterative refinement to improve the computed solution and
747* compute error bounds and backward error estimates for it.
748*
749 CALL dgerfsx( trans, equed, n, nrhs, a, lda, af, ldaf,
750 $ ipiv, r, c, b, ldb, x, ldx, rcond, berr,
751 $ n_err_bnds, err_bnds_norm, err_bnds_comp, nparams, params,
752 $ work, iwork, info )
753*
754* Scale solutions.
755*
756 IF ( colequ .AND. notran ) THEN
757 CALL dlascl2 ( n, nrhs, c, x, ldx )
758 ELSE IF ( rowequ .AND. .NOT.notran ) THEN
759 CALL dlascl2 ( n, nrhs, r, x, ldx )
760 END IF
761*
762 RETURN
763*
764* End of DGESVXX
765
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgeequb(m, n, a, lda, r, c, rowcnd, colcnd, amax, info)
DGEEQUB
Definition dgeequb.f:146
subroutine dgerfsx(trans, equed, n, nrhs, a, lda, af, ldaf, ipiv, r, c, b, ldb, x, ldx, rcond, berr, n_err_bnds, err_bnds_norm, err_bnds_comp, nparams, params, work, iwork, info)
DGERFSX
Definition dgerfsx.f:414
subroutine dgetrf(m, n, a, lda, ipiv, info)
DGETRF
Definition dgetrf.f:108
subroutine dgetrs(trans, n, nrhs, a, lda, ipiv, b, ldb, info)
DGETRS
Definition dgetrs.f:121
double precision function dla_gerpvgrw(n, ncols, a, lda, af, ldaf)
DLA_GERPVGRW
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
subroutine dlaqge(m, n, a, lda, r, c, rowcnd, colcnd, amax, equed)
DLAQGE scales a general rectangular matrix, using row and column scaling factors computed by sgeequ.
Definition dlaqge.f:142
subroutine dlascl2(m, n, d, x, ldx)
DLASCL2 performs diagonal scaling on a matrix.
Definition dlascl2.f:90
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
Here is the call graph for this function:
Here is the caller graph for this function: