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

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

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

Purpose:
 DPBSVX uses the Cholesky factorization A = U**T*U or A = L*L**T to
 compute the solution to a real system of linear equations
    A * X = B,
 where A is an N-by-N symmetric 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**T * U,  if UPLO = 'U', or
       A = L * L**T,  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 DOUBLE PRECISION array, dimension (LDAB,N)
          On entry, the upper or lower triangle of the symmetric 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 DOUBLE PRECISION 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**T*U or A = L*L**T 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**T*U or A = L*L**T.

          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**T*U or A = L*L**T 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 DOUBLE PRECISION 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 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 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 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 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 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 (3*N)
[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:  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 symmetric 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 345 of file dpbsvx.f.

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