Annotation of rpl/lapack/lapack/dsbgv.f, revision 1.4

1.1       bertrand    1:       SUBROUTINE DSBGV( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, W, Z,
                      2:      $                  LDZ, WORK, INFO )
                      3: *
                      4: *  -- LAPACK driver routine (version 3.2) --
                      5: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
                      6: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
                      7: *     November 2006
                      8: *
                      9: *     .. Scalar Arguments ..
                     10:       CHARACTER          JOBZ, UPLO
                     11:       INTEGER            INFO, KA, KB, LDAB, LDBB, LDZ, N
                     12: *     ..
                     13: *     .. Array Arguments ..
                     14:       DOUBLE PRECISION   AB( LDAB, * ), BB( LDBB, * ), W( * ),
                     15:      $                   WORK( * ), Z( LDZ, * )
                     16: *     ..
                     17: *
                     18: *  Purpose
                     19: *  =======
                     20: *
                     21: *  DSBGV computes all the eigenvalues, and optionally, the eigenvectors
                     22: *  of a real generalized symmetric-definite banded eigenproblem, of
                     23: *  the form A*x=(lambda)*B*x. Here A and B are assumed to be symmetric
                     24: *  and banded, and B is also positive definite.
                     25: *
                     26: *  Arguments
                     27: *  =========
                     28: *
                     29: *  JOBZ    (input) CHARACTER*1
                     30: *          = 'N':  Compute eigenvalues only;
                     31: *          = 'V':  Compute eigenvalues and eigenvectors.
                     32: *
                     33: *  UPLO    (input) CHARACTER*1
                     34: *          = 'U':  Upper triangles of A and B are stored;
                     35: *          = 'L':  Lower triangles of A and B are stored.
                     36: *
                     37: *  N       (input) INTEGER
                     38: *          The order of the matrices A and B.  N >= 0.
                     39: *
                     40: *  KA      (input) INTEGER
                     41: *          The number of superdiagonals of the matrix A if UPLO = 'U',
                     42: *          or the number of subdiagonals if UPLO = 'L'. KA >= 0.
                     43: *
                     44: *  KB      (input) INTEGER
                     45: *          The number of superdiagonals of the matrix B if UPLO = 'U',
                     46: *          or the number of subdiagonals if UPLO = 'L'. KB >= 0.
                     47: *
                     48: *  AB      (input/output) DOUBLE PRECISION array, dimension (LDAB, N)
                     49: *          On entry, the upper or lower triangle of the symmetric band
                     50: *          matrix A, stored in the first ka+1 rows of the array.  The
                     51: *          j-th column of A is stored in the j-th column of the array AB
                     52: *          as follows:
                     53: *          if UPLO = 'U', AB(ka+1+i-j,j) = A(i,j) for max(1,j-ka)<=i<=j;
                     54: *          if UPLO = 'L', AB(1+i-j,j)    = A(i,j) for j<=i<=min(n,j+ka).
                     55: *
                     56: *          On exit, the contents of AB are destroyed.
                     57: *
                     58: *  LDAB    (input) INTEGER
                     59: *          The leading dimension of the array AB.  LDAB >= KA+1.
                     60: *
                     61: *  BB      (input/output) DOUBLE PRECISION array, dimension (LDBB, N)
                     62: *          On entry, the upper or lower triangle of the symmetric band
                     63: *          matrix B, stored in the first kb+1 rows of the array.  The
                     64: *          j-th column of B is stored in the j-th column of the array BB
                     65: *          as follows:
                     66: *          if UPLO = 'U', BB(kb+1+i-j,j) = B(i,j) for max(1,j-kb)<=i<=j;
                     67: *          if UPLO = 'L', BB(1+i-j,j)    = B(i,j) for j<=i<=min(n,j+kb).
                     68: *
                     69: *          On exit, the factor S from the split Cholesky factorization
                     70: *          B = S**T*S, as returned by DPBSTF.
                     71: *
                     72: *  LDBB    (input) INTEGER
                     73: *          The leading dimension of the array BB.  LDBB >= KB+1.
                     74: *
                     75: *  W       (output) DOUBLE PRECISION array, dimension (N)
                     76: *          If INFO = 0, the eigenvalues in ascending order.
                     77: *
                     78: *  Z       (output) DOUBLE PRECISION array, dimension (LDZ, N)
                     79: *          If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
                     80: *          eigenvectors, with the i-th column of Z holding the
                     81: *          eigenvector associated with W(i). The eigenvectors are
                     82: *          normalized so that Z**T*B*Z = I.
                     83: *          If JOBZ = 'N', then Z is not referenced.
                     84: *
                     85: *  LDZ     (input) INTEGER
                     86: *          The leading dimension of the array Z.  LDZ >= 1, and if
                     87: *          JOBZ = 'V', LDZ >= N.
                     88: *
                     89: *  WORK    (workspace) DOUBLE PRECISION array, dimension (3*N)
                     90: *
                     91: *  INFO    (output) INTEGER
                     92: *          = 0:  successful exit
                     93: *          < 0:  if INFO = -i, the i-th argument had an illegal value
                     94: *          > 0:  if INFO = i, and i is:
                     95: *             <= N:  the algorithm failed to converge:
                     96: *                    i off-diagonal elements of an intermediate
                     97: *                    tridiagonal form did not converge to zero;
                     98: *             > N:   if INFO = N + i, for 1 <= i <= N, then DPBSTF
                     99: *                    returned INFO = i: B is not positive definite.
                    100: *                    The factorization of B could not be completed and
                    101: *                    no eigenvalues or eigenvectors were computed.
                    102: *
                    103: *  =====================================================================
                    104: *
                    105: *     .. Local Scalars ..
                    106:       LOGICAL            UPPER, WANTZ
                    107:       CHARACTER          VECT
                    108:       INTEGER            IINFO, INDE, INDWRK
                    109: *     ..
                    110: *     .. External Functions ..
                    111:       LOGICAL            LSAME
                    112:       EXTERNAL           LSAME
                    113: *     ..
                    114: *     .. External Subroutines ..
                    115:       EXTERNAL           DPBSTF, DSBGST, DSBTRD, DSTEQR, DSTERF, XERBLA
                    116: *     ..
                    117: *     .. Executable Statements ..
                    118: *
                    119: *     Test the input parameters.
                    120: *
                    121:       WANTZ = LSAME( JOBZ, 'V' )
                    122:       UPPER = LSAME( UPLO, 'U' )
                    123: *
                    124:       INFO = 0
                    125:       IF( .NOT.( WANTZ .OR. LSAME( JOBZ, 'N' ) ) ) THEN
                    126:          INFO = -1
                    127:       ELSE IF( .NOT.( UPPER .OR. LSAME( UPLO, 'L' ) ) ) THEN
                    128:          INFO = -2
                    129:       ELSE IF( N.LT.0 ) THEN
                    130:          INFO = -3
                    131:       ELSE IF( KA.LT.0 ) THEN
                    132:          INFO = -4
                    133:       ELSE IF( KB.LT.0 .OR. KB.GT.KA ) THEN
                    134:          INFO = -5
                    135:       ELSE IF( LDAB.LT.KA+1 ) THEN
                    136:          INFO = -7
                    137:       ELSE IF( LDBB.LT.KB+1 ) THEN
                    138:          INFO = -9
                    139:       ELSE IF( LDZ.LT.1 .OR. ( WANTZ .AND. LDZ.LT.N ) ) THEN
                    140:          INFO = -12
                    141:       END IF
                    142:       IF( INFO.NE.0 ) THEN
                    143:          CALL XERBLA( 'DSBGV ', -INFO )
                    144:          RETURN
                    145:       END IF
                    146: *
                    147: *     Quick return if possible
                    148: *
                    149:       IF( N.EQ.0 )
                    150:      $   RETURN
                    151: *
                    152: *     Form a split Cholesky factorization of B.
                    153: *
                    154:       CALL DPBSTF( UPLO, N, KB, BB, LDBB, INFO )
                    155:       IF( INFO.NE.0 ) THEN
                    156:          INFO = N + INFO
                    157:          RETURN
                    158:       END IF
                    159: *
                    160: *     Transform problem to standard eigenvalue problem.
                    161: *
                    162:       INDE = 1
                    163:       INDWRK = INDE + N
                    164:       CALL DSBGST( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, Z, LDZ,
                    165:      $             WORK( INDWRK ), IINFO )
                    166: *
                    167: *     Reduce to tridiagonal form.
                    168: *
                    169:       IF( WANTZ ) THEN
                    170:          VECT = 'U'
                    171:       ELSE
                    172:          VECT = 'N'
                    173:       END IF
                    174:       CALL DSBTRD( VECT, UPLO, N, KA, AB, LDAB, W, WORK( INDE ), Z, LDZ,
                    175:      $             WORK( INDWRK ), IINFO )
                    176: *
                    177: *     For eigenvalues only, call DSTERF.  For eigenvectors, call SSTEQR.
                    178: *
                    179:       IF( .NOT.WANTZ ) THEN
                    180:          CALL DSTERF( N, W, WORK( INDE ), INFO )
                    181:       ELSE
                    182:          CALL DSTEQR( JOBZ, N, W, WORK( INDE ), Z, LDZ, WORK( INDWRK ),
                    183:      $                INFO )
                    184:       END IF
                    185:       RETURN
                    186: *
                    187: *     End of DSBGV
                    188: *
                    189:       END

CVSweb interface <joel.bertrand@systella.fr>