Visible to Intel only — GUID: GUID-D2395E03-C774-4630-BB51-5EDD4069A9B1
Visible to Intel only — GUID: GUID-D2395E03-C774-4630-BB51-5EDD4069A9B1
?hegvx
Computes selected eigenvalues and, optionally, eigenvectors of a complex generalized Hermitian positive-definite eigenproblem.
lapack_int LAPACKE_chegvx( int matrix_layout, lapack_int itype, char jobz, char range, char uplo, lapack_int n, lapack_complex_float* a, lapack_int lda, lapack_complex_float* b, lapack_int ldb, float vl, float vu, lapack_int il, lapack_int iu, float abstol, lapack_int* m, float* w, lapack_complex_float* z, lapack_int ldz, lapack_int* ifail );
lapack_int LAPACKE_zhegvx( int matrix_layout, lapack_int itype, char jobz, char range, char uplo, lapack_int n, lapack_complex_double* a, lapack_int lda, lapack_complex_double* b, lapack_int ldb, double vl, double vu, lapack_int il, lapack_int iu, double abstol, lapack_int* m, double* w, lapack_complex_double* z, lapack_int ldz, lapack_int* ifail );
- mkl.h
The routine computes selected eigenvalues, and optionally, the eigenvectors of a complex generalized Hermitian positive-definite eigenproblem, of the form
A*x = λ*B*x, A*B*x = λ*x, or B*A*x = λ*x.
Here A and B are assumed to be Hermitian 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.
- matrix_layout
-
Specifies whether matrix storage layout is row major (LAPACK_ROW_MAJOR) or column major (LAPACK_COL_MAJOR).
- itype
-
Must be 1 or 2 or 3. Specifies the problem type to be solved:
if itype = 1, the problem type is A*x = λ*B*x;
if itype = 2, the problem type is A*B*x = λ*x;
if itype = 3, the problem type is B*A*x = λ*x.
- jobz
-
Must be 'N' or 'V'.
If jobz = 'N', then compute eigenvalues only.
If jobz = 'V', then compute eigenvalues and eigenvectors.
- range
-
Must be 'A' or 'V' or 'I'.
If range = 'A', the routine computes all eigenvalues.
If range = 'V', the routine computes eigenvalues w[i] in the half-open interval:
vl<w[i]≤vu.
If range = 'I', the routine computes eigenvalues with indices il to iu.
- uplo
-
Must be 'U' or 'L'.
If uplo = 'U', arrays a and b store the upper triangles of A and B;
If uplo = 'L', arrays a and b store the lower triangles of A and B.
- n
-
The order of the matrices A and B (n≥ 0).
- a, b
-
Arrays:
a (size at least max(1, lda*n)) contains the upper or lower triangle of the Hermitian matrix A, as specified by uplo.
b (size at least max(1, ldb*n)) contains the upper or lower triangle of the Hermitian positive definite matrix B, as specified by uplo.
- lda
-
The leading dimension of a; at least max(1, n).
- ldb
-
The leading dimension of b; at least max(1, n).
- vl, vu
-
If range = 'V', the lower and upper bounds of the interval to be searched for eigenvalues.
Constraint: vl< vu.
If range = 'A' or 'I', vl and vu are not referenced.
- il, iu
-
If range = 'I', the indices in ascending order of the smallest and largest eigenvalues to be returned.
Constraint: 1 ≤il≤iu≤n, if n > 0; il=1 and iu=0
if n = 0.
If range = 'A' or 'V', il and iu are not referenced.
- abstol
-
The absolute error tolerance for the eigenvalues. See Application Notes for more information.
- ldz
-
The leading dimension of the output array z. Constraints:
ldz≥ 1; if jobz = 'V', ldz≥ max(1, n) for column major layout and ldz≥ max(1, m) for row major layout.
- a
-
On exit, the upper triangle (if uplo = 'U') or the lower triangle (if uplo = 'L') of A, including the diagonal, is overwritten.
- b
-
On exit, if info≤n, the part of b containing the matrix is overwritten by the triangular factor U or L from the Cholesky factorization B = UH*U or B = L*LH.
- m
-
The total number of eigenvalues found,
0 ≤m≤n. If range = 'A', m = n, and if range = 'I',
m = iu-il+1.
- w
-
Array, size at least max(1, n).
The first m elements of w contain the selected eigenvalues in ascending order.
- z
-
Array z(size at least max(1, ldz*m) for column major layout and max(1, ldz*n) for row major layout).
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 - 1]. The eigenvectors are normalized as follows:
if itype = 1 or 2, ZH*B*Z = I;
if itype = 3, ZH*inv(B)*Z = I;
If jobz = 'N', then z is not referenced.
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: you 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.
- ifail
-
Array, size at least max(1, n).
If jobz = 'V', then if info = 0, the first m elements of ifail are zero; if info > 0, the ifail contains the indices of the eigenvectors that failed to converge.
If jobz = 'N', then ifail is not referenced.
This function returns a value info.
If info=0, the execution is successful.
If info = -i, the i-th parameter had an illegal value.
If info > 0, cpotrf/zpotrf and cheevx/zheevx returned an error code:
If info = i≤n, cheevx/zheevx failed to converge, and i eigenvectors failed to converge. Their indices are stored in the array ifail;
If info = n + i, for 1 ≤i≤n, then the leading minor of order i of B is not positive-definite. The factorization of B could not be completed and no eigenvalues or eigenvectors were computed.
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+ε*max(|a|,|b|), where ε is the machine precision.
If abstol is less than or equal to zero, then ε*||T||1 will be used in its place, where T is the tridiagonal matrix obtained by reducing C to tridiagonal form, where C is the symmetric matrix of the standard symmetric problem to which the generalized problem is transformed. Eigenvalues will be computed most accurately when abstol is set to twice the underflow threshold 2*?lamch('S'), not zero.
If this routine returns with info > 0, indicating that some eigenvectors did not converge, try setting abstol to 2*?lamch('S').