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

◆ cpbsvx()

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 principal minor of order i is not positive,
    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 principal minor of order i of A
                       is not positive, 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.
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 339 of file cpbsvx.f.

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