F08FPF (ZHEEVX) (PDF version)
F08 Chapter Contents
F08 Chapter Introduction
NAG Library Manual

NAG Library Routine Document

F08FPF (ZHEEVX)

Note:  before using this routine, please read the Users' Note for your implementation to check the interpretation of bold italicised terms and other implementation-dependent details.

+ Contents

    1  Purpose
    7  Accuracy

1  Purpose

F08FPF (ZHEEVX) computes selected eigenvalues and, optionally, eigenvectors of a complex n by n Hermitian matrix A. Eigenvalues and eigenvectors can be selected by specifying either a range of values or a range of indices for the desired eigenvalues.

2  Specification

SUBROUTINE F08FPF ( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, LWORK, RWORK, IWORK, JFAIL, INFO)
INTEGER  N, LDA, IL, IU, M, LDZ, LWORK, IWORK(*), JFAIL(*), INFO
REAL (KIND=nag_wp)  VL, VU, ABSTOL, W(*), RWORK(*)
COMPLEX (KIND=nag_wp)  A(LDA,*), Z(LDZ,*), WORK(max(1,LWORK))
CHARACTER(1)  JOBZ, RANGE, UPLO
The routine may be called by its LAPACK name zheevx.

3  Description

The Hermitian matrix A is first reduced to real tridiagonal form, using unitary similarity transformations. The required eigenvalues and eigenvectors are then computed from the tridiagonal matrix; the method used depends upon whether all, or selected, eigenvalues and eigenvectors are required.

4  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 http://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

5  Parameters

1:     JOBZ – CHARACTER(1)Input
On entry: indicates whether eigenvectors are computed.
JOBZ='N'
Only eigenvalues are computed.
JOBZ='V'
Eigenvalues and eigenvectors are computed.
Constraint: JOBZ='N' or 'V'.
2:     RANGE – CHARACTER(1)Input
On entry: if RANGE='A', all eigenvalues will be found.
If RANGE='V', all eigenvalues in the half-open interval VL,VU will be found.
If RANGE='I', the ILth to IUth eigenvalues will be found.
Constraint: RANGE='A', 'V' or 'I'.
3:     UPLO – CHARACTER(1)Input
On entry: if UPLO='U', the upper triangular part of A is stored.
If UPLO='L', the lower triangular part of A is stored.
Constraint: UPLO='U' or 'L'.
4:     N – INTEGERInput
On entry: n, the order of the matrix A.
Constraint: N0.
5:     A(LDA,*) – COMPLEX (KIND=nag_wp) arrayInput/Output
Note: the second dimension of the array A must be at least max1,N.
On entry: the n by n Hermitian matrix A.
  • If UPLO='U', the upper triangular part of A must be stored and the elements of the array below the diagonal are not referenced.
  • If UPLO='L', the lower triangular part of A must be stored and the elements of the array above the diagonal are not referenced.
On exit: the lower triangle (if UPLO='L') or the upper triangle (if UPLO='U') of A, including the diagonal, is overwritten.
6:     LDA – INTEGERInput
On entry: the first dimension of the array A as declared in the (sub)program from which F08FPF (ZHEEVX) is called.
Constraint: LDAmax1,N.
7:     VL – REAL (KIND=nag_wp)Input
8:     VU – REAL (KIND=nag_wp)Input
On entry: if RANGE='V', the lower and upper bounds of the interval to be searched for eigenvalues.
If RANGE='A' or 'I', VL and VU are not referenced.
Constraint: if RANGE='V', VL<VU.
9:     IL – INTEGERInput
10:   IU – INTEGERInput
On entry: if RANGE='I', the indices (in ascending order) of the smallest and largest eigenvalues to be returned.
If RANGE='A' or 'V', IL and IU are not referenced.
Constraints:
  • if RANGE='I' and N=0, IL=1 and IU=0;
  • if RANGE='I' and N>0, 1 IL IU N .
11:   ABSTOL – REAL (KIND=nag_wp)Input
On entry: 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+ε maxa,b ,
where ε  is the machine precision. If ABSTOL is less than or equal to zero, then ε T1  will be used in its place, where T is the tridiagonal matrix obtained by reducing A to tridiagonal form. Eigenvalues will be computed most accurately when ABSTOL is set to twice the underflow threshold 2 × X02AMF   , not zero. If this routine returns with INFO>0, indicating that some eigenvectors did not converge, try setting ABSTOL to 2 × X02AMF   . See Demmel and Kahan (1990).
12:   M – INTEGEROutput
On exit: the total number of eigenvalues found. 0MN.
If RANGE='A', M=N.
If RANGE='I', M=IU-IL+1.
13:   W(*) – REAL (KIND=nag_wp) arrayOutput
Note: the dimension of the array W must be at least max1,N.
On exit: the first M elements contain the selected eigenvalues in ascending order.
14:   Z(LDZ,*) – COMPLEX (KIND=nag_wp) arrayOutput
Note: the second dimension of the array Z must be at least max1,M if JOBZ='V', and at least 1 otherwise.
On exit: 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 ith column of Z holding the eigenvector associated with Wi;
  • if an eigenvector fails to converge (INFO>0), then that column of Z contains the latest approximation to the eigenvector, and the index of the eigenvector is returned in JFAIL.
If JOBZ='N', Z is not referenced.
Note:  you must ensure that at least max1,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 of at least N must be used.
15:   LDZ – INTEGERInput
On entry: the first dimension of the array Z as declared in the (sub)program from which F08FPF (ZHEEVX) is called.
Constraints:
  • if JOBZ='V', LDZ max1,N ;
  • otherwise LDZ1.
16:   WORK(max1,LWORK) – COMPLEX (KIND=nag_wp) arrayWorkspace
On exit: if INFO=0, the real part of WORK1 contains the minimum value of LWORK required for optimal performance.
17:   LWORK – INTEGERInput
On entry: the dimension of the array WORK as declared in the (sub)program from which F08FPF (ZHEEVX) is called.
If LWORK=-1, a workspace query is assumed; the routine only calculates the optimal size of the WORK array, returns this value as the first entry of the WORK array, and no error message related to LWORK is issued.
Suggested value: for optimal performance, LWORKnb+1×N, where nb is the largest optimal block size for F08FSF (ZHETRD) and for F08FUF (ZUNMTR).
Constraint: LWORKmax1,2×N.
18:   RWORK(*) – REAL (KIND=nag_wp) arrayWorkspace
Note: the dimension of the array RWORK must be at least max1,7×N.
19:   IWORK(*) – INTEGER arrayWorkspace
Note: the dimension of the array IWORK must be at least max1,5×N.
20:   JFAIL(*) – INTEGER arrayOutput
Note: the dimension of the array JFAIL must be at least max1,N.
On exit: if JOBZ='V', then
  • if INFO=0, the first M elements of JFAIL are zero;
  • if INFO>0, JFAIL contains the indices of the eigenvectors that failed to converge.
If JOBZ='N', JFAIL is not referenced.
21:   INFO – INTEGEROutput
On exit: INFO=0 unless the routine detects an error (see Section 6).

6  Error Indicators and Warnings

Errors or warnings detected by the routine:
INFO<0
If INFO=-i, argument i had an illegal value. An explanatory message is output, and execution of the program is terminated.
INFO>0
If INFO=i, then i eigenvectors failed to converge. Their indices are stored in array JFAIL. Please see ABSTOL.

7  Accuracy

The computed eigenvalues and eigenvectors are exact for a nearby matrix A+E, where
E2 = Oε A2 ,
and ε is the machine precision. See Section 4.7 of Anderson et al. (1999) for further details.

8  Further Comments

The total number of floating point operations is proportional to n3.
The real analogue of this routine is F08FBF (DSYEVX).

9  Example

This example finds the eigenvalues in the half-open interval -2,2 , and the corresponding eigenvectors, of the Hermitian matrix
A = 1 2-i 3-i 4-i 2+i 2 3-2i 4-2i 3+i 3+2i 3 4-3i 4+i 4+2i 4+3i 4 .

9.1  Program Text

Program Text (f08fpfe.f90)

9.2  Program Data

Program Data (f08fpfe.d)

9.3  Program Results

Program Results (f08fpfe.r)


F08FPF (ZHEEVX) (PDF version)
F08 Chapter Contents
F08 Chapter Introduction
NAG Library Manual

© The Numerical Algorithms Group Ltd, Oxford, UK. 2012