naginterfaces.library.lapackeig.zhegvx

naginterfaces.library.lapackeig.zhegvx(itype, jobz, erange, uplo, a, b, vl, vu, il, iu, abstol)[source]

zhegvx computes selected eigenvalues and, optionally, eigenvectors of a complex generalized Hermitian-definite eigenproblem, of the form

where and are Hermitian and 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.

For full information please refer to the NAG Library document for f08sp

https://www.nag.com/numeric/nl/nagdoc_29.3/flhtml/f08/f08spf.html

Parameters
itypeint

Specifies the problem type to be solved.

.

.

.

jobzstr, length 1

Indicates whether eigenvectors are computed.

Only eigenvalues are computed.

Eigenvalues and eigenvectors are computed.

erangestr, length 1

If , all eigenvalues will be found.

If , all eigenvalues in the half-open interval will be found.

If , the th to th eigenvalues will be found.

uplostr, length 1

If , the upper triangles of and are stored.

If , the lower triangles of and are stored.

acomplex, array-like, shape

The Hermitian matrix .

bcomplex, array-like, shape

The Hermitian matrix .

vlfloat

If , the lower and upper bounds of the interval to be searched for eigenvalues.

If or , and are not referenced.

vufloat

If , the lower and upper bounds of the interval to be searched for eigenvalues.

If or , and are not referenced.

ilint

If , and specify the indices (in ascending order) of the smallest and largest eigenvalues to be returned, respectively.

If or , and are not referenced.

iuint

If , and specify the indices (in ascending order) of the smallest and largest eigenvalues to be returned, respectively.

If or , and are not referenced.

abstolfloat

The absolute error tolerance for the eigenvalues. An approximate eigenvalue is accepted as converged when it is determined to lie in an interval of width less than or equal to

where is the machine precision. If is less than or equal to zero, then will be used in its place, where is the tridiagonal matrix obtained by reducing to tridiagonal form. Eigenvalues will be computed most accurately when is set to twice the underflow threshold , not zero. If this function returns with = 1 … , indicating that some eigenvectors did not converge, try setting to . See Demmel and Kahan (1990).

Returns
acomplex, ndarray, shape

The lower triangle (if ) or the upper triangle (if ) of , including the diagonal, is overwritten.

bcomplex, ndarray, shape

The triangular factor or from the Cholesky factorization or .

mint

The total number of eigenvalues found. .

If , .

If , .

wfloat, ndarray, shape

The first elements contain the selected eigenvalues in ascending order.

zcomplex, ndarray, shape

If , then

if no exception or warning is raised, the first columns of contain the orthonormal eigenvectors of the matrix corresponding to the selected eigenvalues, with the th column of holding the eigenvector associated with . The eigenvectors are normalized as follows:

if or , ;

if , ;

if an eigenvector fails to converge ( = 1 … ), then that column of contains the latest approximation to the eigenvector, and the index of the eigenvector is returned in .

If , is not referenced.

Note: you must ensure that at least columns are supplied in the array ; if , the exact value of is not known in advance and an upper bound of at least must be used.

jfailint, ndarray, shape

If , then

if no exception or warning is raised, the first elements of are zero;

if = 1 … , the first elements of contains the indices of the eigenvectors that failed to converge.

If , is not referenced.

Raises
NagValueError
(errno )

On entry, error in parameter .

Constraint: , or .

(errno )

On entry, error in parameter .

Constraint: or .

(errno )

On entry, error in parameter .

Constraint: , or .

(errno )

On entry, error in parameter .

Constraint: or .

(errno )

On entry, error in parameter .

Constraint: .

(errno )

On entry, error in parameter .

Constraint: .

(errno )

On entry, error in parameter .

(errno )

On entry, error in parameter .

(errno )

If , for , then the leading minor of order of is not positive definite. The factorization of could not be completed and no eigenvalues or eigenvectors were computed.

Warns
NagAlgorithmicWarning
(errno )

The algorithm failed to converge; eigenvectors failed to converge.

Notes

zhegvx first performs a Cholesky factorization of the matrix as , when or , when . The generalized problem is then reduced to a standard symmetric eigenvalue problem

which is solved for the desired eigenvalues and eigenvectors; the eigenvectors are then backtransformed to give the eigenvectors of the original problem.

For the problem , the eigenvectors are normalized so that the matrix of eigenvectors, , satisfies

where is the diagonal matrix whose diagonal elements are the eigenvalues. For the problem we correspondingly have

and for we have

References

Anderson, E, Bai, Z, Bischof, C, Blackford, S, Demmel, J, Dongarra, J J, Du Croz, J J, Greenbaum, A, Hammarling, S, McKenney, A and Sorensen, D, 1999, LAPACK Users’ Guide, (3rd Edition), SIAM, Philadelphia, https://www.netlib.org/lapack/lug

Demmel, J W and Kahan, W, 1990, Accurate singular values of bidiagonal matrices, SIAM J. Sci. Statist. Comput. (11), 873–912

Golub, G H and Van Loan, C F, 1996, Matrix Computations, (3rd Edition), Johns Hopkins University Press, Baltimore