Contents


NAME

     zhpgvx  -  compute  selected  eigenvalues  and,  optionally,
     eigenvectors  of  a  complex  generalized Hermitian-definite
     eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x,
     or B*A*x=(lambda)*x

SYNOPSIS

     SUBROUTINE ZHPGVX(ITYPE, JOBZ, RANGE, UPLO, N, AP, BP, VL, VU, IL,
           IU, ABSTOL, M, W, Z, LDZ, WORK, RWORK, IWORK, IFAIL, INFO)

     CHARACTER * 1 JOBZ, RANGE, UPLO
     DOUBLE COMPLEX AP(*), BP(*), Z(LDZ,*), WORK(*)
     INTEGER ITYPE, N, IL, IU, M, LDZ, INFO
     INTEGER IWORK(*), IFAIL(*)
     DOUBLE PRECISION VL, VU, ABSTOL
     DOUBLE PRECISION W(*), RWORK(*)

     SUBROUTINE ZHPGVX_64(ITYPE, JOBZ, RANGE, UPLO, N, AP, BP, VL, VU, IL,
           IU, ABSTOL, M, W, Z, LDZ, WORK, RWORK, IWORK, IFAIL, INFO)

     CHARACTER * 1 JOBZ, RANGE, UPLO
     DOUBLE COMPLEX AP(*), BP(*), Z(LDZ,*), WORK(*)
     INTEGER*8 ITYPE, N, IL, IU, M, LDZ, INFO
     INTEGER*8 IWORK(*), IFAIL(*)
     DOUBLE PRECISION VL, VU, ABSTOL
     DOUBLE PRECISION W(*), RWORK(*)

  F95 INTERFACE
     SUBROUTINE HPGVX(ITYPE, JOBZ, RANGE, UPLO, [N], AP, BP, VL, VU, IL,
            IU, ABSTOL, M, W, Z, [LDZ], [WORK], [RWORK], [IWORK], IFAIL,
            [INFO])

     CHARACTER(LEN=1) :: JOBZ, RANGE, UPLO
     COMPLEX(8), DIMENSION(:) :: AP, BP, WORK
     COMPLEX(8), DIMENSION(:,:) :: Z
     INTEGER :: ITYPE, N, IL, IU, M, LDZ, INFO
     INTEGER, DIMENSION(:) :: IWORK, IFAIL
     REAL(8) :: VL, VU, ABSTOL
     REAL(8), DIMENSION(:) :: W, RWORK

     SUBROUTINE HPGVX_64(ITYPE, JOBZ, RANGE, UPLO, [N], AP, BP, VL, VU,
            IL, IU, ABSTOL, M, W, Z, [LDZ], [WORK], [RWORK], [IWORK], IFAIL,
            [INFO])

     CHARACTER(LEN=1) :: JOBZ, RANGE, UPLO
     COMPLEX(8), DIMENSION(:) :: AP, BP, WORK
     COMPLEX(8), DIMENSION(:,:) :: Z
     INTEGER(8) :: ITYPE, N, 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 zhpgvx(int itype, char jobz, char range, char uplo, int
               n,  doublecomplex  *ap,  doublecomplex *bp, double
               vl, double vu, int il, int iu, double abstol,  int
               *m,  double  *w,  doublecomplex  *z,  int ldz, int
               *ifail, int *info);

     void zhpgvx_64(long itype, char jobz, char range, char uplo,
               long n, doublecomplex *ap, doublecomplex *bp, dou-
               ble vl,  double  vu,  long  il,  long  iu,  double
               abstol, long *m, double *w, doublecomplex *z, long
               ldz, long *ifail, long *info);

PURPOSE

     zhpgvx computes selected eigenvalues and, optionally, eigen-
     vectors  of  a complex generalized Hermitian-definite eigen-
     problem, of the form A*x=(lambda)*B*x,  A*Bx=(lambda)*x,  or
     B*A*x=(lambda)*x.  Here A and B are assumed to be Hermitian,
     stored in packed format, and B is  also  positive  definite.
     Eigenvalues  and  eigenvectors can be selected by specifying
     either a range of values or  a  range  of  indices  for  the
     desired eigenvalues.

ARGUMENTS

     ITYPE (input)
               Specifies the problem type to be solved:
               = 1:  A*x = (lambda)*B*x
               = 2:  A*B*x = (lambda)*x
               = 3:  B*A*x = (lambda)*x

     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.

     AP (input/output) COMPLEX*16 array, dimension (N*(N+1)/2)
               On entry, the upper or lower triangle of the  Her-
               mitian  matrix  A,  packed  columnwise in a linear
               array.  The j-th column of  A  is  stored  in  the
               array  AP  as  follows:  if UPLO = 'U', AP(i + (j-
               1)*j/2) = A(i,j) for 1<=i<=j; if UPLO = 'L',  AP(i
               + (j-1)*(2*n-j)/2) = A(i,j) for j<=i<=n.

               On exit, the contents of AP are destroyed.

     BP (input/output) COMPLEX*16 array, dimension (N*(N+1)/2)
               On entry, the upper or lower triangle of the  Her-
               mitian  matrix  B,  packed  columnwise in a linear
               array.  The j-th column of  B  is  stored  in  the
               array  BP  as  follows:  if UPLO = 'U', BP(i + (j-
               1)*j/2) = B(i,j) for 1<=i<=j; if UPLO = 'L',  BP(i
               + (j-1)*(2*n-j)/2) = B(i,j) for j<=i<=n.

               On exit, the triangular factor U  or  L  from  the
               Cholesky  factorization  B = U**H*U or B = L*L**H,
               in the same storage format as B.

     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) DOUBLE PRECISION array, dimension (N)
               On normal exit, the first M elements  contain  the
               selected eigenvalues in ascending order.

     Z (input) COMPLEX*16 array, dimension (LDZ, N)
               If JOBZ = 'N', then Z is not referenced.  If  JOBZ
               =  'V', then if INFO = 0, the first M 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 asso-
               ciated with W(i).  The eigenvectors are normalized
               as follows:  if ITYPE = 1 or 2, Z**H*B*Z =  I;  if
               ITYPE = 3, Z**H*inv(B)*Z = I.

               If an eigenvector fails  to  converge,  then  that
               column  of  Z contains the latest approximation to
               the eigenvector, and the index of the  eigenvector
               is  returned in IFAIL.  Note: the user must ensure
               that at least max(1,M) columns are supplied in the
               array  Z;  if RANGE = 'V', the exact value of M 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)
               COMPLEX*16 array, dimension(2*N)

     RWORK (workspace)
               DOUBLE PRECISION array, dimension(7*N)

     IWORK (workspace)
               INTEGER array, 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 ille-
               gal value
               > 0:  CPPTRF or ZHPEVX returned an error code:
               <= N:  if INFO = i, ZHPEVX failed to  converge;  i
               eigenvectors  failed  to  converge.  Their indices
               are stored in array IFAIL.  > N:   if INFO =  N  +
               i,  for  1  <=  i  <= n, then the leading minor of
               order i of B is not positive definite.   The  fac-
               torization  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