LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine cpbsvx ( character  FACT,
character  UPLO,
integer  N,
integer  KD,
integer  NRHS,
complex, dimension( ldab, * )  AB,
integer  LDAB,
complex, dimension( ldafb, * )  AFB,
integer  LDAFB,
character  EQUED,
real, dimension( * )  S,
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 
)

CPBSVX computes the solution to system of linear equations A * X = B for OTHER matrices

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

Purpose:
 CPBSVX uses the Cholesky factorization A = U**H*U or A = L*L**H to
 compute the solution to a complex system of linear equations
    A * X = B,
 where A is an N-by-N Hermitian positive definite band 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:
       diag(S) * A * diag(S) * inv(diag(S)) * X = diag(S) * 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(S)*A*diag(S) and B by diag(S)*B.

 2. If FACT = 'N' or 'E', the Cholesky decomposition is used to
    factor the matrix A (after equilibration if FACT = 'E') as
       A = U**H * U,  if UPLO = 'U', or
       A = L * L**H,  if UPLO = 'L',
    where U is an upper triangular band matrix, and L is a lower
    triangular band matrix.

 3. If the leading i-by-i principal minor is not positive definite,
    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(S) 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 contains the factored form of A.
                  If EQUED = 'Y', the matrix A has been equilibrated
                  with scaling factors given by S.  AB and AFB will not
                  be 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]UPLO
          UPLO is CHARACTER*1
          = 'U':  Upper triangle of A is stored;
          = 'L':  Lower triangle of A is stored.
[in]N
          N is INTEGER
          The number of linear equations, i.e., the order of the
          matrix A.  N >= 0.
[in]KD
          KD is INTEGER
          The number of superdiagonals of the matrix A if UPLO = 'U',
          or the number of subdiagonals if UPLO = 'L'.  KD >= 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 upper or lower triangle of the Hermitian band
          matrix A, stored in the first KD+1 rows of the array, except
          if FACT = 'F' and EQUED = 'Y', then A must contain the
          equilibrated matrix diag(S)*A*diag(S).  The j-th column of A
          is stored in the j-th column of the array AB as follows:
          if UPLO = 'U', AB(KD+1+i-j,j) = A(i,j) for max(1,j-KD)<=i<=j;
          if UPLO = 'L', AB(1+i-j,j)    = A(i,j) for j<=i<=min(N,j+KD).
          See below for further details.

          On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by
          diag(S)*A*diag(S).
[in]LDAB
          LDAB is INTEGER
          The leading dimension of the array A.  LDAB >= KD+1.
[in,out]AFB
          AFB is COMPLEX array, dimension (LDAFB,N)
          If FACT = 'F', then AFB is an input argument and on entry
          contains the triangular factor U or L from the Cholesky
          factorization A = U**H*U or A = L*L**H of the band matrix
          A, in the same storage format as A (see AB).  If EQUED = 'Y',
          then AFB is the factored form of the equilibrated matrix A.

          If FACT = 'N', then AFB is an output argument and on exit
          returns the triangular factor U or L from the Cholesky
          factorization A = U**H*U or A = L*L**H.

          If FACT = 'E', then AFB is an output argument and on exit
          returns the triangular factor U or L from the Cholesky
          factorization A = U**H*U or A = L*L**H of the equilibrated
          matrix A (see the description of A for the form of the
          equilibrated matrix).
[in]LDAFB
          LDAFB is INTEGER
          The leading dimension of the array AFB.  LDAFB >= KD+1.
[in,out]EQUED
          EQUED is CHARACTER*1
          Specifies the form of equilibration that was done.
          = 'N':  No equilibration (always true if FACT = 'N').
          = 'Y':  Equilibration was done, i.e., A has been replaced by
                  diag(S) * A * diag(S).
          EQUED is an input argument if FACT = 'F'; otherwise, it is an
          output argument.
[in,out]S
          S is REAL array, dimension (N)
          The scale factors for A; not accessed if EQUED = 'N'.  S is
          an input argument if FACT = 'F'; otherwise, S is an output
          argument.  If FACT = 'F' and EQUED = 'Y', each element of S
          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 EQUED = 'Y',
          B is overwritten by diag(S) * 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 if EQUED = 'Y',
          A and B are modified on exit, and the solution to the
          equilibrated system is inv(diag(S))*X.
[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)
[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:  the leading minor of order i of A is
                       not positive definite, so the factorization
                       could not be completed, and the solution has not
                       been 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
Further Details:
  The band storage scheme is illustrated by the following example, when
  N = 6, KD = 2, and UPLO = 'U':

  Two-dimensional storage of the Hermitian matrix A:

     a11  a12  a13
          a22  a23  a24
               a33  a34  a35
                    a44  a45  a46
                         a55  a56
     (aij=conjg(aji))         a66

  Band storage of the upper triangle of A:

      *    *   a13  a24  a35  a46
      *   a12  a23  a34  a45  a56
     a11  a22  a33  a44  a55  a66

  Similarly, if UPLO = 'L' the format of A is as follows:

     a11  a22  a33  a44  a55  a66
     a21  a32  a43  a54  a65   *
     a31  a42  a53  a64   *    *

  Array elements marked * are not used by the routine.

Definition at line 344 of file cpbsvx.f.

344 *
345 * -- LAPACK driver routine (version 3.4.1) --
346 * -- LAPACK is a software package provided by Univ. of Tennessee, --
347 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
348 * April 2012
349 *
350 * .. Scalar Arguments ..
351  CHARACTER equed, fact, uplo
352  INTEGER info, kd, ldab, ldafb, ldb, ldx, n, nrhs
353  REAL rcond
354 * ..
355 * .. Array Arguments ..
356  REAL berr( * ), ferr( * ), rwork( * ), s( * )
357  COMPLEX ab( ldab, * ), afb( ldafb, * ), b( ldb, * ),
358  $ work( * ), x( ldx, * )
359 * ..
360 *
361 * =====================================================================
362 *
363 * .. Parameters ..
364  REAL zero, one
365  parameter ( zero = 0.0e+0, one = 1.0e+0 )
366 * ..
367 * .. Local Scalars ..
368  LOGICAL equil, nofact, rcequ, upper
369  INTEGER i, infequ, j, j1, j2
370  REAL amax, anorm, bignum, scond, smax, smin, smlnum
371 * ..
372 * .. External Functions ..
373  LOGICAL lsame
374  REAL clanhb, slamch
375  EXTERNAL lsame, clanhb, slamch
376 * ..
377 * .. External Subroutines ..
378  EXTERNAL ccopy, clacpy, claqhb, cpbcon, cpbequ, cpbrfs,
379  $ cpbtrf, cpbtrs, xerbla
380 * ..
381 * .. Intrinsic Functions ..
382  INTRINSIC max, min
383 * ..
384 * .. Executable Statements ..
385 *
386  info = 0
387  nofact = lsame( fact, 'N' )
388  equil = lsame( fact, 'E' )
389  upper = lsame( uplo, 'U' )
390  IF( nofact .OR. equil ) THEN
391  equed = 'N'
392  rcequ = .false.
393  ELSE
394  rcequ = lsame( equed, 'Y' )
395  smlnum = slamch( 'Safe minimum' )
396  bignum = one / smlnum
397  END IF
398 *
399 * Test the input parameters.
400 *
401  IF( .NOT.nofact .AND. .NOT.equil .AND. .NOT.lsame( fact, 'F' ) )
402  $ THEN
403  info = -1
404  ELSE IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
405  info = -2
406  ELSE IF( n.LT.0 ) THEN
407  info = -3
408  ELSE IF( kd.LT.0 ) THEN
409  info = -4
410  ELSE IF( nrhs.LT.0 ) THEN
411  info = -5
412  ELSE IF( ldab.LT.kd+1 ) THEN
413  info = -7
414  ELSE IF( ldafb.LT.kd+1 ) THEN
415  info = -9
416  ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
417  $ ( rcequ .OR. lsame( equed, 'N' ) ) ) THEN
418  info = -10
419  ELSE
420  IF( rcequ ) THEN
421  smin = bignum
422  smax = zero
423  DO 10 j = 1, n
424  smin = min( smin, s( j ) )
425  smax = max( smax, s( j ) )
426  10 CONTINUE
427  IF( smin.LE.zero ) THEN
428  info = -11
429  ELSE IF( n.GT.0 ) THEN
430  scond = max( smin, smlnum ) / min( smax, bignum )
431  ELSE
432  scond = one
433  END IF
434  END IF
435  IF( info.EQ.0 ) THEN
436  IF( ldb.LT.max( 1, n ) ) THEN
437  info = -13
438  ELSE IF( ldx.LT.max( 1, n ) ) THEN
439  info = -15
440  END IF
441  END IF
442  END IF
443 *
444  IF( info.NE.0 ) THEN
445  CALL xerbla( 'CPBSVX', -info )
446  RETURN
447  END IF
448 *
449  IF( equil ) THEN
450 *
451 * Compute row and column scalings to equilibrate the matrix A.
452 *
453  CALL cpbequ( uplo, n, kd, ab, ldab, s, scond, amax, infequ )
454  IF( infequ.EQ.0 ) THEN
455 *
456 * Equilibrate the matrix.
457 *
458  CALL claqhb( uplo, n, kd, ab, ldab, s, scond, amax, equed )
459  rcequ = lsame( equed, 'Y' )
460  END IF
461  END IF
462 *
463 * Scale the right-hand side.
464 *
465  IF( rcequ ) THEN
466  DO 30 j = 1, nrhs
467  DO 20 i = 1, n
468  b( i, j ) = s( i )*b( i, j )
469  20 CONTINUE
470  30 CONTINUE
471  END IF
472 *
473  IF( nofact .OR. equil ) THEN
474 *
475 * Compute the Cholesky factorization A = U**H *U or A = L*L**H.
476 *
477  IF( upper ) THEN
478  DO 40 j = 1, n
479  j1 = max( j-kd, 1 )
480  CALL ccopy( j-j1+1, ab( kd+1-j+j1, j ), 1,
481  $ afb( kd+1-j+j1, j ), 1 )
482  40 CONTINUE
483  ELSE
484  DO 50 j = 1, n
485  j2 = min( j+kd, n )
486  CALL ccopy( j2-j+1, ab( 1, j ), 1, afb( 1, j ), 1 )
487  50 CONTINUE
488  END IF
489 *
490  CALL cpbtrf( uplo, n, kd, afb, ldafb, info )
491 *
492 * Return if INFO is non-zero.
493 *
494  IF( info.GT.0 )THEN
495  rcond = zero
496  RETURN
497  END IF
498  END IF
499 *
500 * Compute the norm of the matrix A.
501 *
502  anorm = clanhb( '1', uplo, n, kd, ab, ldab, rwork )
503 *
504 * Compute the reciprocal of the condition number of A.
505 *
506  CALL cpbcon( uplo, n, kd, afb, ldafb, anorm, rcond, work, rwork,
507  $ info )
508 *
509 * Compute the solution matrix X.
510 *
511  CALL clacpy( 'Full', n, nrhs, b, ldb, x, ldx )
512  CALL cpbtrs( uplo, n, kd, nrhs, afb, ldafb, x, ldx, info )
513 *
514 * Use iterative refinement to improve the computed solution and
515 * compute error bounds and backward error estimates for it.
516 *
517  CALL cpbrfs( uplo, n, kd, nrhs, ab, ldab, afb, ldafb, b, ldb, x,
518  $ ldx, ferr, berr, work, rwork, info )
519 *
520 * Transform the solution matrix X to a solution of the original
521 * system.
522 *
523  IF( rcequ ) THEN
524  DO 70 j = 1, nrhs
525  DO 60 i = 1, n
526  x( i, j ) = s( i )*x( i, j )
527  60 CONTINUE
528  70 CONTINUE
529  DO 80 j = 1, nrhs
530  ferr( j ) = ferr( j ) / scond
531  80 CONTINUE
532  END IF
533 *
534 * Set INFO = N+1 if the matrix is singular to working precision.
535 *
536  IF( rcond.LT.slamch( 'Epsilon' ) )
537  $ info = n + 1
538 *
539  RETURN
540 *
541 * End of CPBSVX
542 *
subroutine cpbequ(UPLO, N, KD, AB, LDAB, S, SCOND, AMAX, INFO)
CPBEQU
Definition: cpbequ.f:132
subroutine cpbcon(UPLO, N, KD, AB, LDAB, ANORM, RCOND, WORK, RWORK, INFO)
CPBCON
Definition: cpbcon.f:135
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine cpbrfs(UPLO, N, KD, NRHS, AB, LDAB, AFB, LDAFB, B, LDB, X, LDX, FERR, BERR, WORK, RWORK, INFO)
CPBRFS
Definition: cpbrfs.f:191
subroutine claqhb(UPLO, N, KD, AB, LDAB, S, SCOND, AMAX, EQUED)
CLAQHB scales a Hermitian band matrix, using scaling factors computed by cpbequ.
Definition: claqhb.f:143
real function clanhb(NORM, UPLO, N, K, AB, LDAB, WORK)
CLANHB returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a Hermitian band matrix.
Definition: clanhb.f:134
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 cpbtrf(UPLO, N, KD, AB, LDAB, INFO)
CPBTRF
Definition: cpbtrf.f:144
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine cpbtrs(UPLO, N, KD, NRHS, AB, LDAB, B, LDB, INFO)
CPBTRS
Definition: cpbtrs.f:123

Here is the call graph for this function:

Here is the caller graph for this function: