.TH "SRC/zhpgvx.f" 3 "Version 3.12.0" "LAPACK" \" -*- nroff -*- .ad l .nh .SH NAME SRC/zhpgvx.f .SH SYNOPSIS .br .PP .SS "Functions/Subroutines" .in +1c .ti -1c .RI "subroutine \fBzhpgvx\fP (itype, jobz, range, uplo, n, ap, bp, vl, vu, il, iu, abstol, m, w, z, ldz, work, rwork, iwork, ifail, info)" .br .RI "\fBZHPGVX\fP " .in -1c .SH "Function/Subroutine Documentation" .PP .SS "subroutine zhpgvx (integer itype, character jobz, character range, character uplo, integer n, complex*16, dimension( * ) ap, complex*16, dimension( * ) bp, double precision vl, double precision vu, integer il, integer iu, double precision abstol, integer m, double precision, dimension( * ) w, complex*16, dimension( ldz, * ) z, integer ldz, complex*16, dimension( * ) work, double precision, dimension( * ) rwork, integer, dimension( * ) iwork, integer, dimension( * ) ifail, integer info)" .PP \fBZHPGVX\fP .PP \fBPurpose:\fP .RS 4 .PP .nf !> !> ZHPGVX computes 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\&. 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\&. !> .fi .PP .RE .PP \fBParameters\fP .RS 4 \fIITYPE\fP .PP .nf !> ITYPE is INTEGER !> 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 !> .fi .PP .br \fIJOBZ\fP .PP .nf !> JOBZ is CHARACTER*1 !> = 'N': Compute eigenvalues only; !> = 'V': Compute eigenvalues and eigenvectors\&. !> .fi .PP .br \fIRANGE\fP .PP .nf !> RANGE is CHARACTER*1 !> = '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\&. !> .fi .PP .br \fIUPLO\fP .PP .nf !> UPLO is CHARACTER*1 !> = 'U': Upper triangles of A and B are stored; !> = 'L': Lower triangles of A and B are stored\&. !> .fi .PP .br \fIN\fP .PP .nf !> N is INTEGER !> The order of the matrices A and B\&. N >= 0\&. !> .fi .PP .br \fIAP\fP .PP .nf !> AP is COMPLEX*16 array, dimension (N*(N+1)/2) !> On entry, the upper or lower triangle of the Hermitian 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\&. !> .fi .PP .br \fIBP\fP .PP .nf !> BP is COMPLEX*16 array, dimension (N*(N+1)/2) !> On entry, the upper or lower triangle of the Hermitian 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\&. !> .fi .PP .br \fIVL\fP .PP .nf !> VL is DOUBLE PRECISION !> !> If RANGE='V', the lower bound of the interval to !> be searched for eigenvalues\&. VL < VU\&. !> Not referenced if RANGE = 'A' or 'I'\&. !> .fi .PP .br \fIVU\fP .PP .nf !> VU is DOUBLE PRECISION !> !> If RANGE='V', the upper bound of the interval to !> be searched for eigenvalues\&. VL < VU\&. !> Not referenced if RANGE = 'A' or 'I'\&. !> .fi .PP .br \fIIL\fP .PP .nf !> IL is INTEGER !> !> If RANGE='I', the index of the !> smallest eigenvalue 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'\&. !> .fi .PP .br \fIIU\fP .PP .nf !> IU is INTEGER !> !> If RANGE='I', the index of the !> largest eigenvalue 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'\&. !> .fi .PP .br \fIABSTOL\fP .PP .nf !> ABSTOL is DOUBLE PRECISION !> 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 tridiagonal form\&. !> !> Eigenvalues will be computed most accurately when ABSTOL is !> set to twice the underflow threshold 2*DLAMCH('S'), not zero\&. !> If this routine returns with INFO>0, indicating that some !> eigenvectors did not converge, try setting ABSTOL to !> 2*DLAMCH('S')\&. !> .fi .PP .br \fIM\fP .PP .nf !> M is INTEGER !> The total number of eigenvalues found\&. 0 <= M <= N\&. !> If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1\&. !> .fi .PP .br \fIW\fP .PP .nf !> W is DOUBLE PRECISION array, dimension (N) !> On normal exit, the first M elements contain the selected !> eigenvalues in ascending order\&. !> .fi .PP .br \fIZ\fP .PP .nf !> Z is 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 associated 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\&. !> .fi .PP .br \fILDZ\fP .PP .nf !> LDZ is INTEGER !> The leading dimension of the array Z\&. LDZ >= 1, and if !> JOBZ = 'V', LDZ >= max(1,N)\&. !> .fi .PP .br \fIWORK\fP .PP .nf !> WORK is COMPLEX*16 array, dimension (2*N) !> .fi .PP .br \fIRWORK\fP .PP .nf !> RWORK is DOUBLE PRECISION array, dimension (7*N) !> .fi .PP .br \fIIWORK\fP .PP .nf !> IWORK is INTEGER array, dimension (5*N) !> .fi .PP .br \fIIFAIL\fP .PP .nf !> IFAIL is INTEGER array, dimension (N) !> If JOBZ = 'V', then if INFO = 0, the first M 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\&. !> .fi .PP .br \fIINFO\fP .PP .nf !> INFO is INTEGER !> = 0: successful exit !> < 0: if INFO = -i, the i-th argument had an illegal value !> > 0: ZPPTRF 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 !> principal minor of order i of B is not positive\&. !> The factorization of B could not be completed and !> no eigenvalues or eigenvectors were computed\&. !> .fi .PP .RE .PP \fBAuthor\fP .RS 4 Univ\&. of Tennessee .PP Univ\&. of California Berkeley .PP Univ\&. of Colorado Denver .PP NAG Ltd\&. .RE .PP \fBContributors:\fP .RS 4 Mark Fahey, Department of Mathematics, Univ\&. of Kentucky, USA .RE .PP .PP Definition at line \fB274\fP of file \fBzhpgvx\&.f\fP\&. .SH "Author" .PP Generated automatically by Doxygen for LAPACK from the source code\&.