Contents


NAME

     zhbevx  -  compute  selected  eigenvalues  and,  optionally,
     eigenvectors of a complex Hermitian band matrix A

SYNOPSIS

     SUBROUTINE ZHBEVX(JOBZ, RANGE, UPLO, N, KD, A, LDA, Q, LDQ, VL,
           VU, IL, IU, ABTOL, NFOUND, W, Z, LDZ, WORK, WORK2, IWORK3, IFAIL,
           INFO)

     CHARACTER * 1 JOBZ, RANGE, UPLO
     DOUBLE COMPLEX A(LDA,*), Q(LDQ,*), Z(LDZ,*), WORK(*)
     INTEGER N, KD, LDA, LDQ, IL, IU, NFOUND, LDZ, INFO
     INTEGER IWORK3(*), IFAIL(*)
     DOUBLE PRECISION VL, VU, ABTOL
     DOUBLE PRECISION W(*), WORK2(*)

     SUBROUTINE ZHBEVX_64(JOBZ, RANGE, UPLO, N, KD, A, LDA, Q, LDQ, VL,
           VU, IL, IU, ABTOL, NFOUND, W, Z, LDZ, WORK, WORK2, IWORK3, IFAIL,
           INFO)

     CHARACTER * 1 JOBZ, RANGE, UPLO
     DOUBLE COMPLEX A(LDA,*), Q(LDQ,*), Z(LDZ,*), WORK(*)
     INTEGER*8 N, KD, LDA, LDQ, IL, IU, NFOUND, LDZ, INFO
     INTEGER*8 IWORK3(*), IFAIL(*)
     DOUBLE PRECISION VL, VU, ABTOL
     DOUBLE PRECISION W(*), WORK2(*)

  F95 INTERFACE
     SUBROUTINE HBEVX(JOBZ, RANGE, UPLO, [N], KD, A, [LDA], Q, [LDQ],
            VL, VU, IL, IU, ABTOL, [NFOUND], W, Z, [LDZ], [WORK], [WORK2],
            [IWORK3], IFAIL, [INFO])

     CHARACTER(LEN=1) :: JOBZ, RANGE, UPLO
     COMPLEX(8), DIMENSION(:) :: WORK
     COMPLEX(8), DIMENSION(:,:) :: A, Q, Z
     INTEGER :: N, KD, LDA, LDQ, IL, IU, NFOUND, LDZ, INFO
     INTEGER, DIMENSION(:) :: IWORK3, IFAIL
     REAL(8) :: VL, VU, ABTOL
     REAL(8), DIMENSION(:) :: W, WORK2

     SUBROUTINE HBEVX_64(JOBZ, RANGE, UPLO, [N], KD, A, [LDA], Q, [LDQ],
            VL, VU, IL, IU, ABTOL, [NFOUND], W, Z, [LDZ], [WORK], [WORK2],
            [IWORK3], IFAIL, [INFO])

     CHARACTER(LEN=1) :: JOBZ, RANGE, UPLO
     COMPLEX(8), DIMENSION(:) :: WORK
     COMPLEX(8), DIMENSION(:,:) :: A, Q, Z
     INTEGER(8) :: N, KD, LDA, LDQ, IL, IU, NFOUND, LDZ, INFO
     INTEGER(8), DIMENSION(:) :: IWORK3, IFAIL
     REAL(8) :: VL, VU, ABTOL
     REAL(8), DIMENSION(:) :: W, WORK2

  C INTERFACE
     #include <sunperf.h>

     void zhbevx(char jobz, char range, char uplo, int n, int kd,
               doublecomplex  *a,  int lda, doublecomplex *q, int
               ldq, double vl, double vu, int il, int iu,  double
               abtol,  int  *nfound, double *w, doublecomplex *z,
               int ldz, int *ifail, int *info);

     void zhbevx_64(char jobz, char range,  char  uplo,  long  n,
               long kd, doublecomplex *a, long lda, doublecomplex
               *q, long ldq, double vl, double vu, long il,  long
               iu,  double  abtol, long *nfound, double *w, doub-
               lecomplex *z, long ldz, long *ifail, long *info);

PURPOSE

     zhbevx computes selected eigenvalues and, optionally, eigen-
     vectors  of  a complex Hermitian band matrix A.  Eigenvalues
     and eigenvectors can be  selected  by  specifying  either  a
     range of values or a range of indices for the desired eigen-
     values.

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 triangle of A is stored;
               = 'L':  Lower triangle of A is stored.

     N (input) The order of the matrix A.  N >= 0.

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

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

               On exit, A is overwritten by values generated dur-
               ing the reduction to tridiagonal form.

     LDA (input)
               The leading dimension of the array A.  LDA >= KD +
               1.

     Q (output)
               If JOBZ = 'V', the N-by-N unitary matrix  used  in
               the reduction to tridiagonal form.  If JOBZ = 'N',
               the array Q is not referenced.

     LDQ (input)
               The leading dimension of the array Q.  If  JOBZ  =
               'V', then 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'.

     ABTOL (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

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

               where EPS is the machine precision.  If  ABTOL  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 A to tri-
               diagonal form.

               Eigenvalues will be computed most accurately  when
               ABTOL  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 ABTOL to 2*SLAMCH('S').

               See "Computing Small Singular Values of Bidiagonal
               Matrices  with Guaranteed High Relative Accuracy,"
               by Demmel and Kahan, LAPACK Working Note #3.

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

     W (output)
               The first NFOUND  elements  contain  the  selected
               eigenvalues in ascending order.

     Z (input) If JOBZ = 'V', then if INFO = 0, the first  NFOUND
               columns  of Z contain the orthonormal eigenvectors
               of the matrix  A  corresponding  to  the  selected
               eigenvalues, with the i-th column of Z holding the
               eigenvector associated with W(i).  If an eigenvec-
               tor  fails to converge, then that column of Z con-
               tains the latest approximation to the eigenvector,
               and  the  index  of the eigenvector is returned in
               IFAIL.  If JOBZ = 'N', then Z is  not  referenced.
               Note:   the   user   must  ensure  that  at  least
               max(1,NFOUND) columns are supplied in the array Z;
               if  RANGE  = 'V', the exact value of NFOUND is not
               known in advance and an upper bound must be used.

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

     WORK (workspace)
               dimension(N)

     WORK2 (workspace)
               dimension(7*N)

     IWORK3 (workspace)
               dimension(5*N)

     IFAIL (output)
               If JOBZ = 'V', then if INFO = 0, the first  NFOUND
               elements  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 ille-
               gal value
               > 0:  if INFO = i, then i eigenvectors  failed  to
               converge.   Their  indices  are  stored  in  array
               IFAIL.