LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine cgbsvx ( character  FACT,
character  TRANS,
integer  N,
integer  KL,
integer  KU,
integer  NRHS,
complex, dimension( ldab, * )  AB,
integer  LDAB,
complex, dimension( ldafb, * )  AFB,
integer  LDAFB,
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 
)

CGBSVX computes the solution to system of linear equations A * X = B for GB matrices

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

Purpose:
 CGBSVX uses the LU factorization to compute the solution to a complex
 system of linear equations A * X = B, A**T * X = B, or A**H * X = B,
 where A is a band matrix of order N with KL subdiagonals and KU
 superdiagonals, 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 by this subroutine:

 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 = L * U,
    where L is a product of permutation and unit lower triangular
    matrices with KL subdiagonals, and U is upper triangular with
    KL+KU superdiagonals.

 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, AFB 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.
                  AB, AFB, and IPIV are not modified.
          = 'N':  The matrix A will be copied to AFB and factored.
          = 'E':  The matrix A will be equilibrated if necessary, then
                  copied to AFB 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]KL
          KL is INTEGER
          The number of subdiagonals within the band of A.  KL >= 0.
[in]KU
          KU is INTEGER
          The number of superdiagonals within the band of A.  KU >= 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]AB
          AB is COMPLEX array, dimension (LDAB,N)
          On entry, the matrix A in band storage, in rows 1 to KL+KU+1.
          The j-th column of A is stored in the j-th column of the
          array AB as follows:
          AB(KU+1+i-j,j) = A(i,j) for max(1,j-KU)<=i<=min(N,j+kl)

          If FACT = 'F' and EQUED is not 'N', then A must have been
          equilibrated by the scaling factors in R and/or C.  AB 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]LDAB
          LDAB is INTEGER
          The leading dimension of the array AB.  LDAB >= KL+KU+1.
[in,out]AFB
          AFB is COMPLEX array, dimension (LDAFB,N)
          If FACT = 'F', then AFB is an input argument and on entry
          contains details of the LU factorization of the band matrix
          A, as computed by CGBTRF.  U is stored as an upper triangular
          band matrix with KL+KU superdiagonals in rows 1 to KL+KU+1,
          and the multipliers used during the factorization are stored
          in rows KL+KU+2 to 2*KL+KU+1.  If EQUED .ne. 'N', then AFB is
          the factored form of the equilibrated matrix A.

          If FACT = 'N', then AFB is an output argument and on exit
          returns details of the LU factorization of A.

          If FACT = 'E', then AFB is an output argument and on exit
          returns details of the LU factorization of the equilibrated
          matrix A (see the description of AB for the form of the
          equilibrated matrix).
[in]LDAFB
          LDAFB is INTEGER
          The leading dimension of the array AFB.  LDAFB >= 2*KL+KU+1.
[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 = L*U
          as computed by CGBTRF; 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 = 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 = 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 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 (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.
Date
April 2012

Definition at line 372 of file cgbsvx.f.

372 *
373 * -- LAPACK driver routine (version 3.4.1) --
374 * -- LAPACK is a software package provided by Univ. of Tennessee, --
375 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
376 * April 2012
377 *
378 * .. Scalar Arguments ..
379  CHARACTER equed, fact, trans
380  INTEGER info, kl, ku, ldab, ldafb, ldb, ldx, n, nrhs
381  REAL rcond
382 * ..
383 * .. Array Arguments ..
384  INTEGER ipiv( * )
385  REAL berr( * ), c( * ), ferr( * ), r( * ),
386  $ rwork( * )
387  COMPLEX ab( ldab, * ), afb( ldafb, * ), b( ldb, * ),
388  $ work( * ), x( ldx, * )
389 * ..
390 *
391 * =====================================================================
392 * Moved setting of INFO = N+1 so INFO does not subsequently get
393 * overwritten. Sven, 17 Mar 05.
394 * =====================================================================
395 *
396 * .. Parameters ..
397  REAL zero, one
398  parameter ( zero = 0.0e+0, one = 1.0e+0 )
399 * ..
400 * .. Local Scalars ..
401  LOGICAL colequ, equil, nofact, notran, rowequ
402  CHARACTER norm
403  INTEGER i, infequ, j, j1, j2
404  REAL amax, anorm, bignum, colcnd, rcmax, rcmin,
405  $ rowcnd, rpvgrw, smlnum
406 * ..
407 * .. External Functions ..
408  LOGICAL lsame
409  REAL clangb, clantb, slamch
410  EXTERNAL lsame, clangb, clantb, slamch
411 * ..
412 * .. External Subroutines ..
413  EXTERNAL ccopy, cgbcon, cgbequ, cgbrfs, cgbtrf, cgbtrs,
414  $ clacpy, claqgb, xerbla
415 * ..
416 * .. Intrinsic Functions ..
417  INTRINSIC abs, max, min
418 * ..
419 * .. Executable Statements ..
420 *
421  info = 0
422  nofact = lsame( fact, 'N' )
423  equil = lsame( fact, 'E' )
424  notran = lsame( trans, 'N' )
425  IF( nofact .OR. equil ) THEN
426  equed = 'N'
427  rowequ = .false.
428  colequ = .false.
429  ELSE
430  rowequ = lsame( equed, 'R' ) .OR. lsame( equed, 'B' )
431  colequ = lsame( equed, 'C' ) .OR. lsame( equed, 'B' )
432  smlnum = slamch( 'Safe minimum' )
433  bignum = one / smlnum
434  END IF
435 *
436 * Test the input parameters.
437 *
438  IF( .NOT.nofact .AND. .NOT.equil .AND. .NOT.lsame( fact, 'F' ) )
439  $ THEN
440  info = -1
441  ELSE IF( .NOT.notran .AND. .NOT.lsame( trans, 'T' ) .AND. .NOT.
442  $ lsame( trans, 'C' ) ) THEN
443  info = -2
444  ELSE IF( n.LT.0 ) THEN
445  info = -3
446  ELSE IF( kl.LT.0 ) THEN
447  info = -4
448  ELSE IF( ku.LT.0 ) THEN
449  info = -5
450  ELSE IF( nrhs.LT.0 ) THEN
451  info = -6
452  ELSE IF( ldab.LT.kl+ku+1 ) THEN
453  info = -8
454  ELSE IF( ldafb.LT.2*kl+ku+1 ) THEN
455  info = -10
456  ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
457  $ ( rowequ .OR. colequ .OR. lsame( equed, 'N' ) ) ) THEN
458  info = -12
459  ELSE
460  IF( rowequ ) THEN
461  rcmin = bignum
462  rcmax = zero
463  DO 10 j = 1, n
464  rcmin = min( rcmin, r( j ) )
465  rcmax = max( rcmax, r( j ) )
466  10 CONTINUE
467  IF( rcmin.LE.zero ) THEN
468  info = -13
469  ELSE IF( n.GT.0 ) THEN
470  rowcnd = max( rcmin, smlnum ) / min( rcmax, bignum )
471  ELSE
472  rowcnd = one
473  END IF
474  END IF
475  IF( colequ .AND. info.EQ.0 ) THEN
476  rcmin = bignum
477  rcmax = zero
478  DO 20 j = 1, n
479  rcmin = min( rcmin, c( j ) )
480  rcmax = max( rcmax, c( j ) )
481  20 CONTINUE
482  IF( rcmin.LE.zero ) THEN
483  info = -14
484  ELSE IF( n.GT.0 ) THEN
485  colcnd = max( rcmin, smlnum ) / min( rcmax, bignum )
486  ELSE
487  colcnd = one
488  END IF
489  END IF
490  IF( info.EQ.0 ) THEN
491  IF( ldb.LT.max( 1, n ) ) THEN
492  info = -16
493  ELSE IF( ldx.LT.max( 1, n ) ) THEN
494  info = -18
495  END IF
496  END IF
497  END IF
498 *
499  IF( info.NE.0 ) THEN
500  CALL xerbla( 'CGBSVX', -info )
501  RETURN
502  END IF
503 *
504  IF( equil ) THEN
505 *
506 * Compute row and column scalings to equilibrate the matrix A.
507 *
508  CALL cgbequ( n, n, kl, ku, ab, ldab, r, c, rowcnd, colcnd,
509  $ amax, infequ )
510  IF( infequ.EQ.0 ) THEN
511 *
512 * Equilibrate the matrix.
513 *
514  CALL claqgb( n, n, kl, ku, ab, ldab, r, c, rowcnd, colcnd,
515  $ amax, equed )
516  rowequ = lsame( equed, 'R' ) .OR. lsame( equed, 'B' )
517  colequ = lsame( equed, 'C' ) .OR. lsame( equed, 'B' )
518  END IF
519  END IF
520 *
521 * Scale the right hand side.
522 *
523  IF( notran ) THEN
524  IF( rowequ ) THEN
525  DO 40 j = 1, nrhs
526  DO 30 i = 1, n
527  b( i, j ) = r( i )*b( i, j )
528  30 CONTINUE
529  40 CONTINUE
530  END IF
531  ELSE IF( colequ ) THEN
532  DO 60 j = 1, nrhs
533  DO 50 i = 1, n
534  b( i, j ) = c( i )*b( i, j )
535  50 CONTINUE
536  60 CONTINUE
537  END IF
538 *
539  IF( nofact .OR. equil ) THEN
540 *
541 * Compute the LU factorization of the band matrix A.
542 *
543  DO 70 j = 1, n
544  j1 = max( j-ku, 1 )
545  j2 = min( j+kl, n )
546  CALL ccopy( j2-j1+1, ab( ku+1-j+j1, j ), 1,
547  $ afb( kl+ku+1-j+j1, j ), 1 )
548  70 CONTINUE
549 *
550  CALL cgbtrf( n, n, kl, ku, afb, ldafb, ipiv, info )
551 *
552 * Return if INFO is non-zero.
553 *
554  IF( info.GT.0 ) THEN
555 *
556 * Compute the reciprocal pivot growth factor of the
557 * leading rank-deficient INFO columns of A.
558 *
559  anorm = zero
560  DO 90 j = 1, info
561  DO 80 i = max( ku+2-j, 1 ), min( n+ku+1-j, kl+ku+1 )
562  anorm = max( anorm, abs( ab( i, j ) ) )
563  80 CONTINUE
564  90 CONTINUE
565  rpvgrw = clantb( 'M', 'U', 'N', info, min( info-1, kl+ku ),
566  $ afb( max( 1, kl+ku+2-info ), 1 ), ldafb,
567  $ rwork )
568  IF( rpvgrw.EQ.zero ) THEN
569  rpvgrw = one
570  ELSE
571  rpvgrw = anorm / rpvgrw
572  END IF
573  rwork( 1 ) = rpvgrw
574  rcond = zero
575  RETURN
576  END IF
577  END IF
578 *
579 * Compute the norm of the matrix A and the
580 * reciprocal pivot growth factor RPVGRW.
581 *
582  IF( notran ) THEN
583  norm = '1'
584  ELSE
585  norm = 'I'
586  END IF
587  anorm = clangb( norm, n, kl, ku, ab, ldab, rwork )
588  rpvgrw = clantb( 'M', 'U', 'N', n, kl+ku, afb, ldafb, rwork )
589  IF( rpvgrw.EQ.zero ) THEN
590  rpvgrw = one
591  ELSE
592  rpvgrw = clangb( 'M', n, kl, ku, ab, ldab, rwork ) / rpvgrw
593  END IF
594 *
595 * Compute the reciprocal of the condition number of A.
596 *
597  CALL cgbcon( norm, n, kl, ku, afb, ldafb, ipiv, anorm, rcond,
598  $ work, rwork, info )
599 *
600 * Compute the solution matrix X.
601 *
602  CALL clacpy( 'Full', n, nrhs, b, ldb, x, ldx )
603  CALL cgbtrs( trans, n, kl, ku, nrhs, afb, ldafb, ipiv, x, ldx,
604  $ info )
605 *
606 * Use iterative refinement to improve the computed solution and
607 * compute error bounds and backward error estimates for it.
608 *
609  CALL cgbrfs( trans, n, kl, ku, nrhs, ab, ldab, afb, ldafb, ipiv,
610  $ b, ldb, x, ldx, ferr, berr, work, rwork, info )
611 *
612 * Transform the solution matrix X to a solution of the original
613 * system.
614 *
615  IF( notran ) THEN
616  IF( colequ ) THEN
617  DO 110 j = 1, nrhs
618  DO 100 i = 1, n
619  x( i, j ) = c( i )*x( i, j )
620  100 CONTINUE
621  110 CONTINUE
622  DO 120 j = 1, nrhs
623  ferr( j ) = ferr( j ) / colcnd
624  120 CONTINUE
625  END IF
626  ELSE IF( rowequ ) THEN
627  DO 140 j = 1, nrhs
628  DO 130 i = 1, n
629  x( i, j ) = r( i )*x( i, j )
630  130 CONTINUE
631  140 CONTINUE
632  DO 150 j = 1, nrhs
633  ferr( j ) = ferr( j ) / rowcnd
634  150 CONTINUE
635  END IF
636 *
637 * Set INFO = N+1 if the matrix is singular to working precision.
638 *
639  IF( rcond.LT.slamch( 'Epsilon' ) )
640  $ info = n + 1
641 *
642  rwork( 1 ) = rpvgrw
643  RETURN
644 *
645 * End of CGBSVX
646 *
subroutine cgbtrf(M, N, KL, KU, AB, LDAB, IPIV, INFO)
CGBTRF
Definition: cgbtrf.f:146
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
real function clantb(NORM, UPLO, DIAG, N, K, AB, LDAB, WORK)
CLANTB returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a triangular band matrix.
Definition: clantb.f:143
subroutine claqgb(M, N, KL, KU, AB, LDAB, R, C, ROWCND, COLCND, AMAX, EQUED)
CLAQGB scales a general band matrix, using row and column scaling factors computed by sgbequ...
Definition: claqgb.f:162
subroutine cgbcon(NORM, N, KL, KU, AB, LDAB, IPIV, ANORM, RCOND, WORK, RWORK, INFO)
CGBCON
Definition: cgbcon.f:149
real function clangb(NORM, N, KL, KU, AB, LDAB, WORK)
CLANGB returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: clangb.f:127
subroutine cgbequ(M, N, KL, KU, AB, LDAB, R, C, ROWCND, COLCND, AMAX, INFO)
CGBEQU
Definition: cgbequ.f:156
subroutine clacpy(UPLO, M, N, A, LDA, B, LDB)
CLACPY copies all or part of one two-dimensional array to another.
Definition: clacpy.f:105
subroutine ccopy(N, CX, INCX, CY, INCY)
CCOPY
Definition: ccopy.f:52
subroutine cgbrfs(TRANS, N, KL, KU, NRHS, AB, LDAB, AFB, LDAFB, IPIV, B, LDB, X, LDX, FERR, BERR, WORK, RWORK, INFO)
CGBRFS
Definition: cgbrfs.f:208
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine cgbtrs(TRANS, N, KL, KU, NRHS, AB, LDAB, IPIV, B, LDB, INFO)
CGBTRS
Definition: cgbtrs.f:140

Here is the call graph for this function:

Here is the caller graph for this function: