NAG CL Interface
f08jsc (zsteqr)

Settings help

CL Name Style:


1 Purpose

f08jsc computes all the eigenvalues and, optionally, all the eigenvectors of a complex Hermitian matrix which has been reduced to tridiagonal form.

2 Specification

#include <nag.h>
void  f08jsc (Nag_OrderType order, Nag_ComputeZType compz, Integer n, double d[], double e[], Complex z[], Integer pdz, NagError *fail)
The function may be called by the names: f08jsc, nag_lapackeig_zsteqr or nag_zsteqr.

3 Description

f08jsc computes all the eigenvalues and, optionally, all the eigenvectors of a real symmetric tridiagonal matrix T. In other words, it can compute the spectral factorization of T as
T=ZΛZT,  
where Λ is a diagonal matrix whose diagonal elements are the eigenvalues λi, and Z is the orthogonal matrix whose columns are the eigenvectors zi. Thus
Tzi=λizi,  i=1,2,,n.  
The function stores the real orthogonal matrix Z in a complex array, so that it may also be used to compute all the eigenvalues and eigenvectors of a complex Hermitian matrix A which has been reduced to tridiagonal form T:
A =QTQH, where ​Q​ is unitary =(QZ)Λ(QZ)H.  
In this case, the matrix Q must be formed explicitly and passed to f08jsc, which must be called with compz=Nag_UpdateZ. The functions which must be called to perform the reduction to tridiagonal form and form Q are:
full matrix f08fsc and f08ftc
full matrix, packed storage f08gsc and f08gtc
band matrix f08hsc with vect=Nag_FormQ.
f08jsc uses the implicitly shifted QR algorithm, switching between the QR and QL variants in order to handle graded matrices effectively (see Greenbaum and Dongarra (1980)). The eigenvectors are normalized so that zi2=1, but are determined only to within a complex factor of absolute value 1.
If only the eigenvalues of T are required, it is more efficient to call f08jfc instead. If T is positive definite, small eigenvalues can be computed more accurately by f08juc.

4 References

Golub G H and Van Loan C F (1996) Matrix Computations (3rd Edition) Johns Hopkins University Press, Baltimore
Greenbaum A and Dongarra J J (1980) Experiments with QR/QL methods for the symmetric triangular eigenproblem LAPACK Working Note No. 17 (Technical Report CS-89-92) University of Tennessee, Knoxville https://www.netlib.org/lapack/lawnspdf/lawn17.pdf
Parlett B N (1998) The Symmetric Eigenvalue Problem SIAM, Philadelphia

5 Arguments

1: order Nag_OrderType Input
On entry: the order argument specifies the two-dimensional storage scheme being used, i.e., row-major ordering or column-major ordering. C language defined storage is specified by order=Nag_RowMajor. See Section 3.1.3 in the Introduction to the NAG Library CL Interface for a more detailed explanation of the use of this argument.
Constraint: order=Nag_RowMajor or Nag_ColMajor.
2: compz Nag_ComputeZType Input
On entry: indicates whether the eigenvectors are to be computed.
compz=Nag_NotZ
Only the eigenvalues are computed (and the array z is not referenced).
compz=Nag_UpdateZ
The eigenvalues and eigenvectors of A are computed (and the array z must contain the matrix Q on entry).
compz=Nag_InitZ
The eigenvalues and eigenvectors of T are computed (and the array z is initialized by the function).
Constraint: compz=Nag_NotZ, Nag_UpdateZ or Nag_InitZ.
3: n Integer Input
On entry: n, the order of the matrix T.
Constraint: n0.
4: d[dim] double Input/Output
Note: the dimension, dim, of the array d must be at least max(1,n).
On entry: the diagonal elements of the tridiagonal matrix T.
On exit: the n eigenvalues in ascending order, unless fail.code= NE_CONVERGENCE (in which case see Section 6).
5: e[dim] double Input/Output
Note: the dimension, dim, of the array e must be at least max(1,n-1).
On entry: the off-diagonal elements of the tridiagonal matrix T.
On exit: e is overwritten.
6: z[dim] Complex Input/Output
Note: the dimension, dim, of the array z must be at least
  • max(1,pdz×n) when compz=Nag_UpdateZ or Nag_InitZ;
  • 1 when compz=Nag_NotZ.
The (i,j)th element of the matrix Z is stored in
  • z[(j-1)×pdz+i-1] when order=Nag_ColMajor;
  • z[(i-1)×pdz+j-1] when order=Nag_RowMajor.
On entry: if compz=Nag_UpdateZ, z must contain the unitary matrix Q from the reduction to tridiagonal form.
If compz=Nag_InitZ, z need not be set.
On exit: if compz=Nag_UpdateZ or Nag_InitZ, the n required orthonormal eigenvectors stored as columns of Z; the ith column corresponds to the ith eigenvalue, where i=1,2,,n, unless fail.code= NE_CONVERGENCE.
If compz=Nag_NotZ, z is not referenced.
7: pdz Integer Input
On entry: the stride separating row or column elements (depending on the value of order) in the array z.
Constraints:
  • if compz=Nag_UpdateZ or Nag_InitZ, pdz max(1,n) ;
  • if compz=Nag_NotZ, pdz1.
8: fail NagError * Input/Output
The NAG error argument (see Section 7 in the Introduction to the NAG Library CL Interface).

6 Error Indicators and Warnings

NE_ALLOC_FAIL
Dynamic memory allocation failed.
See Section 3.1.2 in the Introduction to the NAG Library CL Interface for further information.
NE_BAD_PARAM
On entry, argument value had an illegal value.
NE_CONVERGENCE
The algorithm has failed to find all the eigenvalues after a total of 30×n iterations. In this case, d and e contain on exit the diagonal and off-diagonal elements, respectively, of a tridiagonal matrix unitarily similar to T. value off-diagonal elements have not converged to zero.
NE_ENUM_INT_2
On entry, compz=value, pdz=value and n=value.
Constraint: if compz=Nag_UpdateZ or Nag_InitZ, pdz max(1,n) ;
if compz=Nag_NotZ, pdz1.
NE_INT
On entry, n=value.
Constraint: n0.
On entry, pdz=value.
Constraint: pdz>0.
NE_INTERNAL_ERROR
An internal error has occurred in this function. Check the function call and any array sizes. If the call is correct then please contact NAG for assistance.
See Section 7.5 in the Introduction to the NAG Library CL Interface for further information.
NE_NO_LICENCE
Your licence key may have expired or may not have been installed correctly.
See Section 8 in the Introduction to the NAG Library CL Interface for further information.

7 Accuracy

The computed eigenvalues and eigenvectors are exact for a nearby matrix (T+E), where
E2 = O(ε) T2 ,  
and ε is the machine precision.
If λi is an exact eigenvalue and λ~i is the corresponding computed value, then
|λ~i-λi| c (n) ε T2 ,  
where c(n) is a modestly increasing function of n.
If zi is the corresponding exact eigenvector, and z~i is the corresponding computed eigenvector, then the angle θ(z~i,zi) between them is bounded as follows:
θ (z~i,zi) c(n)εT2 minij|λi-λj| .  
Thus the accuracy of a computed eigenvector depends on the gap between its eigenvalue and all the other eigenvalues.

8 Parallelism and Performance

Background information to multithreading can be found in the Multithreading documentation.
f08jsc is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f08jsc 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 function. Please also consult the Users' Note for your implementation for any additional implementation-specific information.

9 Further Comments

The total number of real floating-point operations is typically about 24n2 if compz=Nag_NotZ and about 14n3 if compz=Nag_UpdateZ or Nag_InitZ, but depends on how rapidly the algorithm converges. When compz=Nag_NotZ, the operations are all performed in scalar mode; the additional operations to compute the eigenvectors when compz=Nag_UpdateZ or Nag_InitZ can be vectorized and on some machines may be performed much faster.
The real analogue of this function is f08jec.

10 Example

See Section 10 in f08ftc, f08gtc or f08hsc, which illustrate the use of this function to compute the eigenvalues and eigenvectors of a full or band Hermitian matrix.