NAG FL Interface
f08fqf (zheevd)

Settings help

FL Name Style:


FL Specification Language:


1 Purpose

f08fqf computes all the eigenvalues and, optionally, all the eigenvectors of a complex Hermitian matrix. If the eigenvectors are requested, then it uses a divide-and-conquer algorithm to compute eigenvalues and eigenvectors. However, if only eigenvalues are required, then it uses the Pal–Walker–Kahan variant of the QL or QR algorithm.

2 Specification

Fortran Interface
Subroutine f08fqf ( job, uplo, n, a, lda, w, work, lwork, rwork, lrwork, iwork, liwork, info)
Integer, Intent (In) :: n, lda, lwork, lrwork, liwork
Integer, Intent (Out) :: iwork(max(1,liwork)), info
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,*)
Complex (Kind=nag_wp), Intent (Out) :: work(max(1,lwork))
Character (1), Intent (In) :: job, uplo
C Header Interface
#include <nag.h>
void  f08fqf_ (const char *job, const char *uplo, const Integer *n, Complex a[], const Integer *lda, double w[], Complex work[], const Integer *lwork, double rwork[], const Integer *lrwork, Integer iwork[], const Integer *liwork, Integer *info, const Charlen length_job, const Charlen length_uplo)
The routine may be called by the names f08fqf, nagf_lapackeig_zheevd or its LAPACK name zheevd.

3 Description

f08fqf computes all the eigenvalues and, optionally, all the eigenvectors of a complex Hermitian matrix A. In other words, it can compute the spectral factorization of A as
A=ZΛZH,  
where Λ is a real diagonal matrix whose diagonal elements are the eigenvalues λi, and Z is the (complex) unitary matrix whose columns are the eigenvectors zi. Thus
Azi=λizi,  i=1,2,,n.  

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
Golub G H and Van Loan C F (1996) Matrix Computations (3rd Edition) Johns Hopkins University Press, Baltimore

5 Arguments

1: job Character(1) Input
On entry: indicates whether eigenvectors are computed.
job='N'
Only eigenvalues are computed.
job='V'
Eigenvalues and eigenvectors are computed.
Constraint: job='N' or 'V'.
2: uplo Character(1) Input
On entry: indicates whether the upper or lower triangular part of A is stored.
uplo='U'
The upper triangular part of A is stored.
uplo='L'
The lower triangular part of A is stored.
Constraint: uplo='U' or 'L'.
3: n Integer Input
On entry: n, the order of the matrix A.
Constraint: n0.
4: 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: if job='V', a is overwritten by the unitary matrix Z which contains the eigenvectors of A.
5: lda Integer Input
On entry: the first dimension of the array a as declared in the (sub)program from which f08fqf is called.
Constraint: ldamax(1,n).
6: w(*) Real (Kind=nag_wp) array Output
Note: the dimension of the array w must be at least max(1,n).
On exit: the eigenvalues of the matrix A in ascending order.
7: work(max(1,lwork)) Complex (Kind=nag_wp) array Workspace
On exit: if info=0, the real part of work(1) contains the required minimal size of lwork.
8: lwork Integer Input
On entry: the dimension of the array work as declared in the (sub)program from which f08fqf is called.
If lwork=−1, a workspace query is assumed; the routine only calculates the minimum dimension of the work array, returns this value as the first entry of the work array, and no error message related to lwork is issued.
Constraints:
  • if n1, lwork1 or lwork=−1;
  • if job='N' and n>1, lworkn+1 or lwork=−1;
  • if job='V' and n>1, lworkn×(n+2) or lwork=−1.
9: rwork(max(1,lrwork)) Real (Kind=nag_wp) array Workspace
On exit: if info=0, rwork(1) contains the required minimal size of lrwork.
10: lrwork Integer Input
On entry: the dimension of the array rwork as declared in the (sub)program from which f08fqf is called.
If lrwork=−1, a workspace query is assumed; the routine only calculates the minimum dimension of the rwork array, returns this value as the first entry of the rwork array, and no error message related to lrwork is issued.
Constraints:
  • if n1, lrwork1 or lrwork=−1;
  • if job='N' and n>1, lrworkn or lrwork=−1;
  • if job='V' and n>1, lrwork2×n2+5×n+1 or lrwork=−1.
11: iwork(max(1,liwork)) Integer array Workspace
On exit: if info=0, iwork(1) contains the required minimal size of liwork.
12: liwork Integer Input
On entry: the dimension of the array iwork as declared in the (sub)program from which f08fqf is called.
If liwork=−1, a workspace query is assumed; the routine only calculates the minimum dimension of the iwork array, returns this value as the first entry of the iwork array, and no error message related to liwork is issued.
Constraints:
  • if n1, liwork1 or liwork=−1;
  • if job='N' and n>1, liwork1 or liwork=−1;
  • if job='V' and n>1, liwork5×n+3 or liwork=−1.
13: 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
If info=value and job='N', the algorithm failed to converge; value elements of an intermediate tridiagonal form did not converge to zero; if info=value and job='V', then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and column value/(n+1) through value mod (n+1).

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

Background information to multithreading can be found in the Multithreading documentation.
f08fqf is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f08fqf 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 real analogue of this routine is f08fcf.

10 Example

This example computes all the eigenvalues and eigenvectors of the Hermitian matrix A, where
A = ( 1.0+0.0i 2.0-1.0i 3.0-1.0i 4.0-1.0i 2.0+1.0i 2.0+0.0i 3.0-2.0i 4.0-2.0i 3.0+1.0i 3.0+2.0i 3.0+0.0i 4.0-3.0i 4.0+1.0i 4.0+2.0i 4.0+3.0i 4.0+0.0i ) .  
The example program for f08fqf illustrates the computation of error bounds for the eigenvalues and eigenvectors.

10.1 Program Text

Program Text (f08fqfe.f90)

10.2 Program Data

Program Data (f08fqfe.d)

10.3 Program Results

Program Results (f08fqfe.r)