NAG FL Interface
f08fpf (zheevx)

Settings help

FL Name Style:


FL Specification Language:


1 Purpose

f08fpf computes selected eigenvalues and, optionally, eigenvectors of a complex n×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

Fortran Interface
Subroutine f08fpf ( jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol, m, w, z, ldz, work, lwork, rwork, iwork, jfail, info)
Integer, Intent (In) :: n, lda, il, iu, ldz, lwork
Integer, Intent (Inout) :: iwork(*), jfail(*)
Integer, Intent (Out) :: m, info
Real (Kind=nag_wp), Intent (In) :: vl, vu, abstol
Real (Kind=nag_wp), Intent (Inout) :: w(*), rwork(*)
Complex (Kind=nag_wp), Intent (Inout) :: a(lda,*), z(ldz,*)
Complex (Kind=nag_wp), Intent (Out) :: work(max(1,lwork))
Character (1), Intent (In) :: jobz, range, uplo
C Header Interface
#include <nag.h>
void  f08fpf_ (const char *jobz, const char *range, const char *uplo, const Integer *n, Complex a[], const Integer *lda, const double *vl, const double *vu, const Integer *il, const Integer *iu, const double *abstol, Integer *m, double w[], Complex z[], const Integer *ldz, Complex work[], const Integer *lwork, double rwork[], Integer iwork[], Integer jfail[], Integer *info, const Charlen length_jobz, const Charlen length_range, const Charlen length_uplo)
The routine may be called by the names f08fpf, nagf_lapackeig_zheevx or 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 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

5 Arguments

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 Integer Input
On entry: n, the order of the matrix A.
Constraint: n0.
5: a(lda,*) Complex (Kind=nag_wp) array Input/Output
Note: the second dimension of the array a must be at least max(1,n).
On entry: the n×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 Integer Input
On entry: the first dimension of the array a as declared in the (sub)program from which f08fpf is called.
Constraint: ldamax(1,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 Integer Input
10: iu Integer Input
On entry: if range='I', il and iu specify the indices (in ascending order) of the smallest and largest eigenvalues to be returned, respectively.
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+ε max(|a|,|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 Integer Output
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) array Output
Note: the dimension of the array w must be at least max(1,n).
On exit: the first m elements contain the selected eigenvalues in ascending order.
14: z(ldz,*) Complex (Kind=nag_wp) array Output
Note: the second dimension of the array z must be at least max(1,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 w(i);
  • 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 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 of at least n must be used.
15: ldz Integer Input
On entry: the first dimension of the array z as declared in the (sub)program from which f08fpf is called.
Constraints:
  • if jobz='V', ldz max(1,n) ;
  • otherwise ldz1.
16: work(max(1,lwork)) Complex (Kind=nag_wp) array Workspace
On exit: if info=0, the real part of work(1) contains the minimum value of lwork required for optimal performance.
17: lwork Integer Input
On entry: the dimension of the array work as declared in the (sub)program from which f08fpf 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, lwork(nb+1)×n, where nb is the largest optimal block size for f08fsf and for f08fuf.
Constraint: lworkmax(1,2×n) or lwork=−1.
18: rwork(*) Real (Kind=nag_wp) array Workspace
Note: the dimension of the array rwork must be at least max(1,7×n).
19: iwork(*) Integer array Workspace
Note: the dimension of the array iwork must be at least max(1,5×n).
20: jfail(*) Integer array Output
Note: the dimension of the array jfail must be at least max(1,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 Integer Output
On exit: info=0 unless the routine detects an error (see Section 6).

6 Error Indicators and Warnings

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
The algorithm failed to converge; value eigenvectors did not converge. Their indices are stored in array jfail.

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 Parallelism and Performance

f08fpf is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f08fpf makes calls to BLAS and/or LAPACK routines, which may be threaded within the vendor library used by this implementation. Consult the documentation for the vendor library for further information.
Please consult the X06 Chapter Introduction for information on how to control and interrogate the OpenMP environment used within this routine. Please also consult the Users' Note for your implementation for any additional implementation-specific information.

9 Further Comments

The total number of floating-point operations is proportional to n3.
The real analogue of this routine is f08fbf.

10 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 ) .  

10.1 Program Text

Program Text (f08fpfe.f90)

10.2 Program Data

Program Data (f08fpfe.d)

10.3 Program Results

Program Results (f08fpfe.r)