Contents


NAME

     zhbgvx - compute all the eigenvalues,  and  optionally,  the
     eigenvectors  of  a  complex  generalized Hermitian-definite
     banded eigenproblem, of the form A*x=(lambda)*B*x

SYNOPSIS

     SUBROUTINE ZHBGVX(JOBZ, RANGE, UPLO, N, KA, KB, AB, LDAB, BB, LDBB,
           Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, RWORK, IWORK,
           IFAIL, INFO)

     CHARACTER * 1 JOBZ, RANGE, UPLO
     DOUBLE COMPLEX AB(LDAB,*), BB(LDBB,*),  Q(LDQ,*),  Z(LDZ,*),
     WORK(*)
     INTEGER N, KA, KB, LDAB, LDBB, LDQ, IL, IU, M, LDZ, INFO
     INTEGER IWORK(*), IFAIL(*)
     DOUBLE PRECISION VL, VU, ABSTOL
     DOUBLE PRECISION W(*), RWORK(*)

     SUBROUTINE ZHBGVX_64(JOBZ, RANGE, UPLO, N, KA, KB, AB, LDAB, BB,
           LDBB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, RWORK,
           IWORK, IFAIL, INFO)

     CHARACTER * 1 JOBZ, RANGE, UPLO
     DOUBLE COMPLEX AB(LDAB,*), BB(LDBB,*),  Q(LDQ,*),  Z(LDZ,*),
     WORK(*)
     INTEGER*8 N, KA, KB, LDAB, LDBB, LDQ, IL, IU, M, LDZ, INFO
     INTEGER*8 IWORK(*), IFAIL(*)
     DOUBLE PRECISION VL, VU, ABSTOL
     DOUBLE PRECISION W(*), RWORK(*)

  F95 INTERFACE
     SUBROUTINE HBGVX(JOBZ, RANGE, UPLO, [N], KA, KB, AB, [LDAB], BB,
            [LDBB], Q, [LDQ], VL, VU, IL, IU, ABSTOL, M, W, Z, [LDZ], [WORK],
            [RWORK], [IWORK], IFAIL, [INFO])

     CHARACTER(LEN=1) :: JOBZ, RANGE, UPLO
     COMPLEX(8), DIMENSION(:) :: WORK
     COMPLEX(8), DIMENSION(:,:) :: AB, BB, Q, Z
     INTEGER :: N, KA, KB, LDAB, LDBB, LDQ, IL, IU, M, LDZ, INFO
     INTEGER, DIMENSION(:) :: IWORK, IFAIL
     REAL(8) :: VL, VU, ABSTOL
     REAL(8), DIMENSION(:) :: W, RWORK

     SUBROUTINE HBGVX_64(JOBZ, RANGE, UPLO, [N], KA, KB, AB, [LDAB], BB,
            [LDBB], Q, [LDQ], VL, VU, IL, IU, ABSTOL, M, W, Z, [LDZ], [WORK],
            [RWORK], [IWORK], IFAIL, [INFO])

     CHARACTER(LEN=1) :: JOBZ, RANGE, UPLO
     COMPLEX(8), DIMENSION(:) :: WORK
     COMPLEX(8), DIMENSION(:,:) :: AB, BB, Q, Z
     INTEGER(8) :: N, KA, KB, LDAB, LDBB, LDQ, IL,  IU,  M,  LDZ,
     INFO
     INTEGER(8), DIMENSION(:) :: IWORK, IFAIL
     REAL(8) :: VL, VU, ABSTOL
     REAL(8), DIMENSION(:) :: W, RWORK

  C INTERFACE
     #include <sunperf.h>

     void zhbgvx(char jobz, char range, char uplo, int n, int ka,
               int kb, doublecomplex *ab, int ldab, doublecomplex
               *bb, int ldbb, doublecomplex *q, int  ldq,  double
               vl,  double vu, int il, int iu, double abstol, int
               *m, double *w,  doublecomplex  *z,  int  ldz,  int
               *ifail, int *info);

     void zhbgvx_64(char jobz, char range,  char  uplo,  long  n,
               long  ka,  long  kb, doublecomplex *ab, long ldab,
               doublecomplex *bb, long  ldbb,  doublecomplex  *q,
               long  ldq, double vl, double vu, long il, long iu,
               double abstol, long *m, double  *w,  doublecomplex
               *z, long ldz, long *ifail, long *info);

PURPOSE

     zhbgvx computes all the  eigenvalues,  and  optionally,  the
     eigenvectors  of  a  complex  generalized Hermitian-definite
     banded eigenproblem, of the form  A*x=(lambda)*B*x.  Here  A
     and  B are assumed to be Hermitian and banded, and B is also
     positive definite.   Eigenvalues  and  eigenvectors  can  be
     selected  by  specifying  either all eigenvalues, a range of
     values or a range of indices for the desired eigenvalues.

ARGUMENTS

     JOBZ (input)
               = 'N':  Compute eigenvalues only;
               = 'V':  Compute eigenvalues and eigenvectors.

     RANGE (input)
               = 'A': all eigenvalues will be found;
               = 'V': all eigenvalues in the  half-open  interval
               (VL,VU]  will  be  found; = 'I': the IL-th through
               IU-th eigenvalues will be found.

     UPLO (input)
               = 'U':  Upper triangles of A and B are stored;
               = 'L':  Lower triangles of A and B are stored.

     N (input) The order of the matrices A and B.  N >= 0.

     KA (input)
               The number of superdiagonals of the  matrix  A  if
               UPLO  = 'U', or the number of subdiagonals if UPLO
               = 'L'. KA >= 0.

     KB (input)
               The number of superdiagonals of the  matrix  B  if
               UPLO  = 'U', or the number of subdiagonals if UPLO
               = 'L'. KB >= 0.

     AB (input/output)
               On entry, the upper or lower triangle of the  Her-
               mitian  band  matrix  A,  stored in the first ka+1
               rows of the array.  The j-th column of A is stored
               in the j-th column of the array AB as follows:  if
               UPLO = 'U', AB(ka+1+i-j,j) = A(i,j)  for  max(1,j-
               ka)<=i<=j;  if UPLO = 'L', AB(1+i-j,j)    = A(i,j)
               for j<=i<=min(n,j+ka).

               On exit, the contents of AB are destroyed.

     LDAB (input)
               The leading dimension of the array  AB.   LDAB  >=
               KA+1.

     BB (input/output)
               On entry, the upper or lower triangle of the  Her-
               mitian  band  matrix  B,  stored in the first kb+1
               rows of the array.  The j-th column of B is stored
               in the j-th column of the array BB as follows:  if
               UPLO = 'U', BB(kb+1+i-j,j) = B(i,j)  for  max(1,j-
               kb)<=i<=j;  if UPLO = 'L', BB(1+i-j,j)    = B(i,j)
               for j<=i<=min(n,j+kb).

               On exit, the factor S from the split Cholesky fac-
               torization B = S**H*S, as returned by CPBSTF.

     LDBB (input)
               The leading dimension of the array  BB.   LDBB  >=
               KB+1.
     Q (output)
               If JOBZ = 'V',  the  n-by-n  matrix  used  in  the
               reduction  of A*x = (lambda)*B*x to standard form,
               i.e. C*x = (lambda)*x, and consequently C to  tri-
               diagonal  form.  If JOBZ = 'N', the array Q is not
               referenced.

     LDQ (input)
               The leading dimension of the array Q.  If  JOBZ  =
               'N', LDQ >= 1. If JOBZ = 'V', LDQ >= max(1,N).

     VL (input)
               If RANGE='V', the lower and upper  bounds  of  the
               interval  to be searched for eigenvalues. VL < VU.
               Not referenced if RANGE = 'A' or 'I'.

     VU (input)
               If RANGE='V', the lower and upper  bounds  of  the
               interval  to be searched for eigenvalues. VL < VU.
               Not referenced if RANGE = 'A' or 'I'.

     IL (input)
               If RANGE='I', the indices (in ascending order)  of
               the   smallest   and  largest  eigenvalues  to  be
               returned.  1 <= IL <= IU <= N, if N > 0;  IL  =  1
               and  IU  =  0 if N = 0.  Not referenced if RANGE =
               'A' or 'V'.

     IU (input)
               If RANGE='I', the indices (in ascending order)  of
               the   smallest   and  largest  eigenvalues  to  be
               returned.  1 <= IL <= IU <= N, if N > 0;  IL  =  1
               and  IU  =  0 if N = 0.  Not referenced if RANGE =
               'A' or 'V'.

     ABSTOL (input)
               The absolute error tolerance for the  eigenvalues.
               An approximate eigenvalue is accepted as converged
               when it is determined to lie in an interval  [a,b]
               of width less than or equal to

               ABSTOL + EPS *   max( |a|,|b| ) ,

               where EPS is the machine precision.  If ABSTOL  is
               less than or equal to zero, then  EPS*|T|  will be
               used in its place, where |T| is the 1-norm of  the
               tridiagonal matrix obtained by reducing AP to tri-
               diagonal form.

               Eigenvalues will be computed most accurately  when
               ABSTOL  is  set  to  twice the underflow threshold
               2*SLAMCH('S'), not zero.  If this routine  returns
               with INFO>0, indicating that some eigenvectors did
               not converge, try setting ABSTOL to 2*SLAMCH('S').

     M (output)
               The total number of eigenvalues found.  0 <= M  <=
               N.  If RANGE = 'A', M = N, and if RANGE = 'I', M =
               IU-IL+1.

     W (output)
               If INFO = 0, the eigenvalues in ascending order.

     Z (input) If JOBZ = 'V', then if INFO = 0,  Z  contains  the
               matrix  Z of eigenvectors, with the i-th column of
               Z holding the eigenvector  associated  with  W(i).
               The eigenvectors are normalized so that Z**H*B*Z =
               I.  If JOBZ = 'N', then Z is not referenced.

     LDZ (input)
               The leading dimension of the array Z.  LDZ  >=  1,
               and if JOBZ = 'V', LDZ >= N.

     WORK (workspace)
               dimension(N)

     RWORK (workspace)
               dimension(7*N)

     IWORK (workspace)
               dimension(5*N)

     IFAIL (output)
               If JOBZ = 'V', then if INFO = 0, the first M  ele-
               ments  of IFAIL are zero.  If INFO > 0, then IFAIL
               contains the  indices  of  the  eigenvectors  that
               failed  to converge.  If JOBZ = 'N', then IFAIL is
               not referenced.

     INFO (output)
               = 0:  successful exit
               < 0:  if INFO =  -i,  the  i-th  argument  had  an
               illegal value
               > 0:  if INFO = i, and i is:
               <= N:  then i  eigenvectors  failed  to  converge.
               Their  indices  are  stored  in array IFAIL.  > N:
               if INFO = N + i, for 1 <= i <= N, then CPBSTF
               returned INFO = i: B  is  not  positive  definite.
               The  factorization of B could not be completed and
               no eigenvalues or eigenvectors were computed.

FURTHER DETAILS

     Based on contributions by
        Mark Fahey, Department of Mathematics, Univ. of Kentucky,
     USA