NAG Library Routine Document

f08frf  (zheevr)

 Contents

    1  Purpose
    7  Accuracy

1
Purpose

f08frf (zheevr) 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

Fortran Interface
Subroutine f08frf ( jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol, m, w, z, ldz, isuppz, work, lwork, rwork, lrwork, iwork, liwork, info)
Integer, Intent (In):: n, lda, il, iu, ldz, lwork, lrwork, liwork
Integer, Intent (Inout):: isuppz(*)
Integer, Intent (Out):: m, iwork(max(1,liwork)), info
Real (Kind=nag_wp), Intent (In):: vl, vu, abstol
Real (Kind=nag_wp), Intent (Inout):: w(*)
Real (Kind=nag_wp), Intent (Out):: rwork(max(1,lrwork))
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 nagmk26.h
void  f08frf_ ( 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, Integer isuppz[], Complex work[], const Integer *lwork, double rwork[], const Integer *lrwork, Integer iwork[], const Integer *liwork, Integer *info, const Charlen length_jobz, const Charlen length_range, const Charlen length_uplo)
The routine may be called by its LAPACK name zheevr.

3
Description

The Hermitian matrix is first reduced to a real tridiagonal matrix T, using unitary similarity transformations. Then whenever possible, f08frf (zheevr) computes the eigenspectrum using Relatively Robust Representations. f08frf (zheevr) computes eigenvalues by the dqds algorithm, while orthogonal eigenvectors are computed from various ‘good’ LDLT representations (also known as Relatively Robust Representations). Gram–Schmidt orthogonalization is avoided as far as possible. More specifically, the various steps of the algorithm are as follows. For the ith unreduced block of T:
(a) compute T - σi I = Li Di LiT , such that Li Di LiT  is a relatively robust representation,
(b) compute the eigenvalues, λj, of Li Di LiT  to high relative accuracy by the dqds algorithm,
(c) if there is a cluster of close eigenvalues, ‘choose’ σi close to the cluster, and go to (a),
(d) given the approximate eigenvalue λj of Li Di LiT , compute the corresponding eigenvector by forming a rank-revealing twisted factorization.
The desired accuracy of the output can be specified by the argument abstol. For more details, see Dhillon (1997) and Parlett and Dhillon (2000).

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
Barlow J and Demmel J W (1990) Computing accurate eigensystems of scaled diagonally dominant matrices SIAM J. Numer. Anal. 27 762–791
Demmel J W and Kahan W (1990) Accurate singular values of bidiagonal matrices SIAM J. Sci. Statist. Comput. 11 873–912
Dhillon I (1997) A new On2 algorithm for the symmetric tridiagonal eigenvalue/eigenvector problem Computer Science Division Technical Report No. UCB//CSD-97-971 UC Berkeley
Golub G H and Van Loan C F (1996) Matrix Computations (3rd Edition) Johns Hopkins University Press, Baltimore
Parlett B N and Dhillon I S (2000) Relatively robust representations of symmetric tridiagonals Linear Algebra Appl. 309 121–151

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.
For range='V' or 'I' and iu-il<n-1, f08jjf (dstebz) and f08jxf (zstein) are called.
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:     alda* – 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 f08frf (zheevr) 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 real tridiagonal matrix obtained by reducing A to tridiagonal form. See Demmel and Kahan (1990).
If high relative accuracy is important, set abstol to x02amf  , although doing so does not currently guarantee that eigenvalues are computed to high relative accuracy. See Barlow and Demmel (1990) for a discussion of which matrices can define their eigenvalues to high relative accuracy.
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:   zldz* – 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', 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 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 f08frf (zheevr) is called.
Constraints:
  • if jobz='V', ldz max1,n ;
  • otherwise ldz1.
16:   isuppz* – Integer arrayOutput
Note: the dimension of the array isuppz must be at least max1,2×m.
On exit: the support of the eigenvectors in z, i.e., the indices indicating the nonzero elements in z. The ith eigenvector is nonzero only in elements isuppz2×i-1 through isuppz2×i. Implemented only for range='A' or 'I' and iu-il=n-1.
17:   workmax1,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.
18:   lwork – IntegerInput
On entry: the dimension of the array work as declared in the (sub)program from which f08frf (zheevr) is called.
If lwork=-1, a workspace query is assumed; the routine only calculates the optimal sizes of the work, rwork and iwork arrays, returns these values as the first entries of the work, rwork and iwork arrays, and no error message related to lwork, lrwork or liwork 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.
19:   rworkmax1,lrwork – Real (Kind=nag_wp) arrayWorkspace
On exit: if info=0, rwork1 returns the optimal (and minimal) lrwork.
20:   lrwork – IntegerInput
On entry: the dimension of the array rwork as declared in the (sub)program from which f08frf (zheevr) is called.
If lrwork=-1, a workspace query is assumed; the routine only calculates the optimal sizes of the work, rwork and iwork arrays, returns these values as the first entries of the work, rwork and iwork arrays, and no error message related to lwork, lrwork or liwork is issued.
Constraint: lrworkmax1,24×n.
21:   iworkmax1,liwork – Integer arrayWorkspace
On exit: if info=0, iwork1 returns the optimal (and minimal) liwork.
22:   liwork – IntegerInput
On entry: the dimension of the array iwork as declared in the (sub)program from which f08frf (zheevr) is called.
If liwork=-1, a workspace query is assumed; the routine only calculates the optimal sizes of the work, rwork and iwork arrays, returns these values as the first entries of the work, rwork and iwork arrays, and no error message related to lwork, lrwork or liwork is issued.
Constraint: liwork max1,10×n .
23:   info – IntegerOutput
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
f08frf (zheevr) failed to converge.

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

f08frf (zheevr) is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f08frf (zheevr) 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 f08fdf (dsyevr).

10
Example

This example finds the eigenvalues with indices in the range 2,3 , 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 .  
Information on required and provided workspace is also output.

10.1
Program Text

Program Text (f08frfe.f90)

10.2
Program Data

Program Data (f08frfe.d)

10.3
Program Results

Program Results (f08frfe.r)

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