LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dgesvx ( 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, dimension( * )  FERR,
double precision, dimension( * )  BERR,
double precision, dimension( * )  WORK,
integer, dimension( * )  IWORK,
integer  INFO 
)

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

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

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

 Error bounds on the solution and a condition estimate are also
 provided.
Description:
 The following steps are performed:

 1. If FACT = 'E', real 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.  If the
    reciprocal of the condition number is less than machine precision,
    INFO = N+1 is returned as a warning, but 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. Iterative refinement is applied to improve the computed solution
    matrix and calculate error bounds and backward error estimates
    for it.

 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.
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  (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.
[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.
[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 or INFO = N+1, 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
          The estimate of the reciprocal condition number of the matrix
          A after equilibration (if done).  If RCOND is less than the
          machine precision (in particular, if RCOND = 0), the matrix
          is singular to working precision.  This condition is
          indicated by a return code of INFO > 0.
[out]FERR
          FERR is DOUBLE PRECISION array, dimension (NRHS)
          The estimated forward error bound for each solution vector
          X(j) (the j-th column of the solution matrix X).
          If XTRUE is the true solution corresponding to X(j), FERR(j)
          is an estimated upper bound for the magnitude of the largest
          element in (X(j) - XTRUE) divided by the magnitude of the
          largest element in X(j).  The estimate is as reliable as
          the estimate for RCOND, and is almost always a slight
          overestimate of the true error.
[out]BERR
          BERR is DOUBLE PRECISION array, dimension (NRHS)
          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).
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (4*N)
          On exit, WORK(1) contains the reciprocal pivot growth
          factor norm(A)/norm(U). The "max absolute element" norm is
          used. If WORK(1) 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, condition
          estimator RCOND, and forward error bound FERR could be
          unreliable. If factorization fails with 0<INFO<=N, then
          WORK(1) contains the reciprocal pivot growth factor for the
          leading INFO columns of A.
[out]IWORK
          IWORK is INTEGER array, dimension (N)
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  if INFO = i, and i is
                <= N:  U(i,i) 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+1: U is nonsingular, but RCOND is less than machine
                       precision, meaning that the matrix is singular
                       to working precision.  Nevertheless, the
                       solution and error bounds are computed because
                       there are a number of situations where the
                       computed solution can be more accurate than the
                       value of RCOND would suggest.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
April 2012

Definition at line 351 of file dgesvx.f.

351 *
352 * -- LAPACK driver routine (version 3.4.1) --
353 * -- LAPACK is a software package provided by Univ. of Tennessee, --
354 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
355 * April 2012
356 *
357 * .. Scalar Arguments ..
358  CHARACTER equed, fact, trans
359  INTEGER info, lda, ldaf, ldb, ldx, n, nrhs
360  DOUBLE PRECISION rcond
361 * ..
362 * .. Array Arguments ..
363  INTEGER ipiv( * ), iwork( * )
364  DOUBLE PRECISION a( lda, * ), af( ldaf, * ), b( ldb, * ),
365  $ berr( * ), c( * ), ferr( * ), r( * ),
366  $ work( * ), x( ldx, * )
367 * ..
368 *
369 * =====================================================================
370 *
371 * .. Parameters ..
372  DOUBLE PRECISION zero, one
373  parameter ( zero = 0.0d+0, one = 1.0d+0 )
374 * ..
375 * .. Local Scalars ..
376  LOGICAL colequ, equil, nofact, notran, rowequ
377  CHARACTER norm
378  INTEGER i, infequ, j
379  DOUBLE PRECISION amax, anorm, bignum, colcnd, rcmax, rcmin,
380  $ rowcnd, rpvgrw, smlnum
381 * ..
382 * .. External Functions ..
383  LOGICAL lsame
384  DOUBLE PRECISION dlamch, dlange, dlantr
385  EXTERNAL lsame, dlamch, dlange, dlantr
386 * ..
387 * .. External Subroutines ..
388  EXTERNAL dgecon, dgeequ, dgerfs, dgetrf, dgetrs, dlacpy,
389  $ dlaqge, xerbla
390 * ..
391 * .. Intrinsic Functions ..
392  INTRINSIC max, min
393 * ..
394 * .. Executable Statements ..
395 *
396  info = 0
397  nofact = lsame( fact, 'N' )
398  equil = lsame( fact, 'E' )
399  notran = lsame( trans, 'N' )
400  IF( nofact .OR. equil ) THEN
401  equed = 'N'
402  rowequ = .false.
403  colequ = .false.
404  ELSE
405  rowequ = lsame( equed, 'R' ) .OR. lsame( equed, 'B' )
406  colequ = lsame( equed, 'C' ) .OR. lsame( equed, 'B' )
407  smlnum = dlamch( 'Safe minimum' )
408  bignum = one / smlnum
409  END IF
410 *
411 * Test the input parameters.
412 *
413  IF( .NOT.nofact .AND. .NOT.equil .AND. .NOT.lsame( fact, 'F' ) )
414  $ THEN
415  info = -1
416  ELSE IF( .NOT.notran .AND. .NOT.lsame( trans, 'T' ) .AND. .NOT.
417  $ lsame( trans, 'C' ) ) THEN
418  info = -2
419  ELSE IF( n.LT.0 ) THEN
420  info = -3
421  ELSE IF( nrhs.LT.0 ) THEN
422  info = -4
423  ELSE IF( lda.LT.max( 1, n ) ) THEN
424  info = -6
425  ELSE IF( ldaf.LT.max( 1, n ) ) THEN
426  info = -8
427  ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
428  $ ( rowequ .OR. colequ .OR. lsame( equed, 'N' ) ) ) THEN
429  info = -10
430  ELSE
431  IF( rowequ ) THEN
432  rcmin = bignum
433  rcmax = zero
434  DO 10 j = 1, n
435  rcmin = min( rcmin, r( j ) )
436  rcmax = max( rcmax, r( j ) )
437  10 CONTINUE
438  IF( rcmin.LE.zero ) THEN
439  info = -11
440  ELSE IF( n.GT.0 ) THEN
441  rowcnd = max( rcmin, smlnum ) / min( rcmax, bignum )
442  ELSE
443  rowcnd = one
444  END IF
445  END IF
446  IF( colequ .AND. info.EQ.0 ) THEN
447  rcmin = bignum
448  rcmax = zero
449  DO 20 j = 1, n
450  rcmin = min( rcmin, c( j ) )
451  rcmax = max( rcmax, c( j ) )
452  20 CONTINUE
453  IF( rcmin.LE.zero ) THEN
454  info = -12
455  ELSE IF( n.GT.0 ) THEN
456  colcnd = max( rcmin, smlnum ) / min( rcmax, bignum )
457  ELSE
458  colcnd = one
459  END IF
460  END IF
461  IF( info.EQ.0 ) THEN
462  IF( ldb.LT.max( 1, n ) ) THEN
463  info = -14
464  ELSE IF( ldx.LT.max( 1, n ) ) THEN
465  info = -16
466  END IF
467  END IF
468  END IF
469 *
470  IF( info.NE.0 ) THEN
471  CALL xerbla( 'DGESVX', -info )
472  RETURN
473  END IF
474 *
475  IF( equil ) THEN
476 *
477 * Compute row and column scalings to equilibrate the matrix A.
478 *
479  CALL dgeequ( n, n, a, lda, r, c, rowcnd, colcnd, amax, infequ )
480  IF( infequ.EQ.0 ) THEN
481 *
482 * Equilibrate the matrix.
483 *
484  CALL dlaqge( n, n, a, lda, r, c, rowcnd, colcnd, amax,
485  $ equed )
486  rowequ = lsame( equed, 'R' ) .OR. lsame( equed, 'B' )
487  colequ = lsame( equed, 'C' ) .OR. lsame( equed, 'B' )
488  END IF
489  END IF
490 *
491 * Scale the right hand side.
492 *
493  IF( notran ) THEN
494  IF( rowequ ) THEN
495  DO 40 j = 1, nrhs
496  DO 30 i = 1, n
497  b( i, j ) = r( i )*b( i, j )
498  30 CONTINUE
499  40 CONTINUE
500  END IF
501  ELSE IF( colequ ) THEN
502  DO 60 j = 1, nrhs
503  DO 50 i = 1, n
504  b( i, j ) = c( i )*b( i, j )
505  50 CONTINUE
506  60 CONTINUE
507  END IF
508 *
509  IF( nofact .OR. equil ) THEN
510 *
511 * Compute the LU factorization of A.
512 *
513  CALL dlacpy( 'Full', n, n, a, lda, af, ldaf )
514  CALL dgetrf( n, n, af, ldaf, ipiv, info )
515 *
516 * Return if INFO is non-zero.
517 *
518  IF( info.GT.0 ) THEN
519 *
520 * Compute the reciprocal pivot growth factor of the
521 * leading rank-deficient INFO columns of A.
522 *
523  rpvgrw = dlantr( 'M', 'U', 'N', info, info, af, ldaf,
524  $ work )
525  IF( rpvgrw.EQ.zero ) THEN
526  rpvgrw = one
527  ELSE
528  rpvgrw = dlange( 'M', n, info, a, lda, work ) / rpvgrw
529  END IF
530  work( 1 ) = rpvgrw
531  rcond = zero
532  RETURN
533  END IF
534  END IF
535 *
536 * Compute the norm of the matrix A and the
537 * reciprocal pivot growth factor RPVGRW.
538 *
539  IF( notran ) THEN
540  norm = '1'
541  ELSE
542  norm = 'I'
543  END IF
544  anorm = dlange( norm, n, n, a, lda, work )
545  rpvgrw = dlantr( 'M', 'U', 'N', n, n, af, ldaf, work )
546  IF( rpvgrw.EQ.zero ) THEN
547  rpvgrw = one
548  ELSE
549  rpvgrw = dlange( 'M', n, n, a, lda, work ) / rpvgrw
550  END IF
551 *
552 * Compute the reciprocal of the condition number of A.
553 *
554  CALL dgecon( norm, n, af, ldaf, anorm, rcond, work, iwork, info )
555 *
556 * Compute the solution matrix X.
557 *
558  CALL dlacpy( 'Full', n, nrhs, b, ldb, x, ldx )
559  CALL dgetrs( trans, n, nrhs, af, ldaf, ipiv, x, ldx, info )
560 *
561 * Use iterative refinement to improve the computed solution and
562 * compute error bounds and backward error estimates for it.
563 *
564  CALL dgerfs( trans, n, nrhs, a, lda, af, ldaf, ipiv, b, ldb, x,
565  $ ldx, ferr, berr, work, iwork, info )
566 *
567 * Transform the solution matrix X to a solution of the original
568 * system.
569 *
570  IF( notran ) THEN
571  IF( colequ ) THEN
572  DO 80 j = 1, nrhs
573  DO 70 i = 1, n
574  x( i, j ) = c( i )*x( i, j )
575  70 CONTINUE
576  80 CONTINUE
577  DO 90 j = 1, nrhs
578  ferr( j ) = ferr( j ) / colcnd
579  90 CONTINUE
580  END IF
581  ELSE IF( rowequ ) THEN
582  DO 110 j = 1, nrhs
583  DO 100 i = 1, n
584  x( i, j ) = r( i )*x( i, j )
585  100 CONTINUE
586  110 CONTINUE
587  DO 120 j = 1, nrhs
588  ferr( j ) = ferr( j ) / rowcnd
589  120 CONTINUE
590  END IF
591 *
592  work( 1 ) = rpvgrw
593 *
594 * Set INFO = N+1 if the matrix is singular to working precision.
595 *
596  IF( rcond.LT.dlamch( 'Epsilon' ) )
597  $ info = n + 1
598  RETURN
599 *
600 * End of DGESVX
601 *
subroutine dgetrf(M, N, A, LDA, IPIV, INFO)
DGETRF
Definition: dgetrf.f:110
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:144
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine dgerfs(TRANS, N, NRHS, A, LDA, AF, LDAF, IPIV, B, LDB, X, LDX, FERR, BERR, WORK, IWORK, INFO)
DGERFS
Definition: dgerfs.f:187
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
subroutine dgeequ(M, N, A, LDA, R, C, ROWCND, COLCND, AMAX, INFO)
DGEEQU
Definition: dgeequ.f:141
subroutine dgetrs(TRANS, N, NRHS, A, LDA, IPIV, B, LDB, INFO)
DGETRS
Definition: dgetrs.f:123
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
double precision function dlange(NORM, M, N, A, LDA, WORK)
DLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: dlange.f:116
double precision function dlantr(NORM, UPLO, DIAG, M, N, A, LDA, WORK)
DLANTR returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a trapezoidal or triangular matrix.
Definition: dlantr.f:143
subroutine dgecon(NORM, N, A, LDA, ANORM, RCOND, WORK, IWORK, INFO)
DGECON
Definition: dgecon.f:126
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: