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

◆ dgesvx()

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 (MAX(1,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.

Definition at line 346 of file dgesvx.f.

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