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

◆ cgesvx()

subroutine cgesvx ( character  fact,
character  trans,
integer  n,
integer  nrhs,
complex, dimension( lda, * )  a,
integer  lda,
complex, dimension( ldaf, * )  af,
integer  ldaf,
integer, dimension( * )  ipiv,
character  equed,
real, dimension( * )  r,
real, dimension( * )  c,
complex, dimension( ldb, * )  b,
integer  ldb,
complex, dimension( ldx, * )  x,
integer  ldx,
real  rcond,
real, dimension( * )  ferr,
real, dimension( * )  berr,
complex, dimension( * )  work,
real, dimension( * )  rwork,
integer  info 
)

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

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

Purpose:
 CGESVX uses the LU factorization to compute the solution to a complex
 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  (Conjugate 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 COMPLEX 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 COMPLEX 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 CGETRF.  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 CGETRF; 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 REAL 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 REAL 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 COMPLEX 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 COMPLEX 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 REAL
          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 REAL 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 REAL 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 COMPLEX array, dimension (2*N)
[out]RWORK
          RWORK is REAL array, dimension (MAX(1,2*N))
          On exit, RWORK(1) contains the reciprocal pivot growth
          factor norm(A)/norm(U). The "max absolute element" norm is
          used. If RWORK(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
          RWORK(1) contains the reciprocal pivot growth factor for the
          leading INFO columns of A.
[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 347 of file cgesvx.f.

350*
351* -- LAPACK driver routine --
352* -- LAPACK is a software package provided by Univ. of Tennessee, --
353* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
354*
355* .. Scalar Arguments ..
356 CHARACTER EQUED, FACT, TRANS
357 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS
358 REAL RCOND
359* ..
360* .. Array Arguments ..
361 INTEGER IPIV( * )
362 REAL BERR( * ), C( * ), FERR( * ), R( * ),
363 $ RWORK( * )
364 COMPLEX A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
365 $ WORK( * ), X( LDX, * )
366* ..
367*
368* =====================================================================
369*
370* .. Parameters ..
371 REAL ZERO, ONE
372 parameter( zero = 0.0e+0, one = 1.0e+0 )
373* ..
374* .. Local Scalars ..
375 LOGICAL COLEQU, EQUIL, NOFACT, NOTRAN, ROWEQU
376 CHARACTER NORM
377 INTEGER I, INFEQU, J
378 REAL AMAX, ANORM, BIGNUM, COLCND, RCMAX, RCMIN,
379 $ ROWCND, RPVGRW, SMLNUM
380* ..
381* .. External Functions ..
382 LOGICAL LSAME
383 REAL CLANGE, CLANTR, SLAMCH
384 EXTERNAL lsame, clange, clantr, slamch
385* ..
386* .. External Subroutines ..
387 EXTERNAL cgecon, cgeequ, cgerfs, cgetrf, cgetrs, clacpy,
388 $ claqge, xerbla
389* ..
390* .. Intrinsic Functions ..
391 INTRINSIC max, min
392* ..
393* .. Executable Statements ..
394*
395 info = 0
396 nofact = lsame( fact, 'N' )
397 equil = lsame( fact, 'E' )
398 notran = lsame( trans, 'N' )
399 IF( nofact .OR. equil ) THEN
400 equed = 'N'
401 rowequ = .false.
402 colequ = .false.
403 ELSE
404 rowequ = lsame( equed, 'R' ) .OR. lsame( equed, 'B' )
405 colequ = lsame( equed, 'C' ) .OR. lsame( equed, 'B' )
406 smlnum = slamch( 'Safe minimum' )
407 bignum = one / smlnum
408 END IF
409*
410* Test the input parameters.
411*
412 IF( .NOT.nofact .AND. .NOT.equil .AND. .NOT.lsame( fact, 'F' ) )
413 $ THEN
414 info = -1
415 ELSE IF( .NOT.notran .AND. .NOT.lsame( trans, 'T' ) .AND. .NOT.
416 $ lsame( trans, 'C' ) ) THEN
417 info = -2
418 ELSE IF( n.LT.0 ) THEN
419 info = -3
420 ELSE IF( nrhs.LT.0 ) THEN
421 info = -4
422 ELSE IF( lda.LT.max( 1, n ) ) THEN
423 info = -6
424 ELSE IF( ldaf.LT.max( 1, n ) ) THEN
425 info = -8
426 ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
427 $ ( rowequ .OR. colequ .OR. lsame( equed, 'N' ) ) ) THEN
428 info = -10
429 ELSE
430 IF( rowequ ) THEN
431 rcmin = bignum
432 rcmax = zero
433 DO 10 j = 1, n
434 rcmin = min( rcmin, r( j ) )
435 rcmax = max( rcmax, r( j ) )
436 10 CONTINUE
437 IF( rcmin.LE.zero ) THEN
438 info = -11
439 ELSE IF( n.GT.0 ) THEN
440 rowcnd = max( rcmin, smlnum ) / min( rcmax, bignum )
441 ELSE
442 rowcnd = one
443 END IF
444 END IF
445 IF( colequ .AND. info.EQ.0 ) THEN
446 rcmin = bignum
447 rcmax = zero
448 DO 20 j = 1, n
449 rcmin = min( rcmin, c( j ) )
450 rcmax = max( rcmax, c( j ) )
451 20 CONTINUE
452 IF( rcmin.LE.zero ) THEN
453 info = -12
454 ELSE IF( n.GT.0 ) THEN
455 colcnd = max( rcmin, smlnum ) / min( rcmax, bignum )
456 ELSE
457 colcnd = one
458 END IF
459 END IF
460 IF( info.EQ.0 ) THEN
461 IF( ldb.LT.max( 1, n ) ) THEN
462 info = -14
463 ELSE IF( ldx.LT.max( 1, n ) ) THEN
464 info = -16
465 END IF
466 END IF
467 END IF
468*
469 IF( info.NE.0 ) THEN
470 CALL xerbla( 'CGESVX', -info )
471 RETURN
472 END IF
473*
474 IF( equil ) THEN
475*
476* Compute row and column scalings to equilibrate the matrix A.
477*
478 CALL cgeequ( n, n, a, lda, r, c, rowcnd, colcnd, amax, infequ )
479 IF( infequ.EQ.0 ) THEN
480*
481* Equilibrate the matrix.
482*
483 CALL claqge( n, n, a, lda, r, c, rowcnd, colcnd, amax,
484 $ equed )
485 rowequ = lsame( equed, 'R' ) .OR. lsame( equed, 'B' )
486 colequ = lsame( equed, 'C' ) .OR. lsame( equed, 'B' )
487 END IF
488 END IF
489*
490* Scale the right hand side.
491*
492 IF( notran ) THEN
493 IF( rowequ ) THEN
494 DO 40 j = 1, nrhs
495 DO 30 i = 1, n
496 b( i, j ) = r( i )*b( i, j )
497 30 CONTINUE
498 40 CONTINUE
499 END IF
500 ELSE IF( colequ ) THEN
501 DO 60 j = 1, nrhs
502 DO 50 i = 1, n
503 b( i, j ) = c( i )*b( i, j )
504 50 CONTINUE
505 60 CONTINUE
506 END IF
507*
508 IF( nofact .OR. equil ) THEN
509*
510* Compute the LU factorization of A.
511*
512 CALL clacpy( 'Full', n, n, a, lda, af, ldaf )
513 CALL cgetrf( n, n, af, ldaf, ipiv, info )
514*
515* Return if INFO is non-zero.
516*
517 IF( info.GT.0 ) THEN
518*
519* Compute the reciprocal pivot growth factor of the
520* leading rank-deficient INFO columns of A.
521*
522 rpvgrw = clantr( 'M', 'U', 'N', info, info, af, ldaf,
523 $ rwork )
524 IF( rpvgrw.EQ.zero ) THEN
525 rpvgrw = one
526 ELSE
527 rpvgrw = clange( 'M', n, info, a, lda, rwork ) /
528 $ rpvgrw
529 END IF
530 rwork( 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 = clange( norm, n, n, a, lda, rwork )
545 rpvgrw = clantr( 'M', 'U', 'N', n, n, af, ldaf, rwork )
546 IF( rpvgrw.EQ.zero ) THEN
547 rpvgrw = one
548 ELSE
549 rpvgrw = clange( 'M', n, n, a, lda, rwork ) / rpvgrw
550 END IF
551*
552* Compute the reciprocal of the condition number of A.
553*
554 CALL cgecon( norm, n, af, ldaf, anorm, rcond, work, rwork, info )
555*
556* Compute the solution matrix X.
557*
558 CALL clacpy( 'Full', n, nrhs, b, ldb, x, ldx )
559 CALL cgetrs( 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 cgerfs( trans, n, nrhs, a, lda, af, ldaf, ipiv, b, ldb, x,
565 $ ldx, ferr, berr, work, rwork, 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* Set INFO = N+1 if the matrix is singular to working precision.
593*
594 IF( rcond.LT.slamch( 'Epsilon' ) )
595 $ info = n + 1
596*
597 rwork( 1 ) = rpvgrw
598 RETURN
599*
600* End of CGESVX
601*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine cgecon(norm, n, a, lda, anorm, rcond, work, rwork, info)
CGECON
Definition cgecon.f:132
subroutine cgeequ(m, n, a, lda, r, c, rowcnd, colcnd, amax, info)
CGEEQU
Definition cgeequ.f:140
subroutine cgerfs(trans, n, nrhs, a, lda, af, ldaf, ipiv, b, ldb, x, ldx, ferr, berr, work, rwork, info)
CGERFS
Definition cgerfs.f:186
subroutine cgetrf(m, n, a, lda, ipiv, info)
CGETRF
Definition cgetrf.f:108
subroutine cgetrs(trans, n, nrhs, a, lda, ipiv, b, ldb, info)
CGETRS
Definition cgetrs.f:121
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:103
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
real function clange(norm, m, n, a, lda, work)
CLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition clange.f:115
real function clantr(norm, uplo, diag, m, n, a, lda, work)
CLANTR returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition clantr.f:142
subroutine claqge(m, n, a, lda, r, c, rowcnd, colcnd, amax, equed)
CLAQGE scales a general rectangular matrix, using row and column scaling factors computed by sgeequ.
Definition claqge.f:143
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: