NAG Library Function Document
nag_zpteqr (f08juc)
1 Purpose
nag_zpteqr (f08juc) computes all the eigenvalues and, optionally, all the eigenvectors of a complex Hermitian positive definite matrix which has been reduced to tridiagonal form.
2 Specification
#include <nag.h> |
#include <nagf08.h> |
void |
nag_zpteqr (Nag_OrderType order,
Nag_ComputeZType compz,
Integer n,
double d[],
double e[],
Complex z[],
Integer pdz,
NagError *fail) |
|
3 Description
nag_zpteqr (f08juc) computes all the eigenvalues and, optionally, all the eigenvectors of a real symmetric positive definite tridiagonal matrix
.
In other words, it can compute the spectral factorization of
as
where
is a diagonal matrix whose diagonal elements are the eigenvalues
, and
is the orthogonal matrix whose columns are the eigenvectors
. Thus
The function stores the real orthogonal matrix
in a complex array, so that it may be used to compute all the eigenvalues and eigenvectors of a complex Hermitian positive definite matrix
which has been reduced to tridiagonal form
:
In this case, the matrix
must be formed explicitly and passed to nag_zpteqr (f08juc), which must be called with
. The functions which must be called to perform the reduction to tridiagonal form and form
are:
nag_zpteqr (f08juc) first factorizes
as
where
is unit lower bidiagonal and
is diagonal. It forms the bidiagonal matrix
, and then calls
nag_zbdsqr (f08msc) to compute the singular values of
which are the same as the eigenvalues of
. The method used by the function allows high relative accuracy to be achieved in the small eigenvalues of
. The eigenvectors are normalized so that
, but are determined only to within a complex factor of absolute value
.
4 References
Barlow J and Demmel J W (1990) Computing accurate eigensystems of scaled diagonally dominant matrices SIAM J. Numer. Anal. 27 762–791
5 Arguments
- 1:
– Nag_OrderTypeInput
-
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
. See
Section 3.2.1.3 in the Essential Introduction for a more detailed explanation of the use of this argument.
Constraint:
or .
- 2:
– Nag_ComputeZTypeInput
-
On entry: indicates whether the eigenvectors are to be computed.
- Only the eigenvalues are computed (and the array z is not referenced).
- The eigenvalues and eigenvectors of are computed (and the array z must contain the matrix on entry).
- The eigenvalues and eigenvectors of are computed (and the array z is initialized by the function).
Constraint:
, or .
- 3:
– IntegerInput
-
On entry: , the order of the matrix .
Constraint:
.
- 4:
– doubleInput/Output
-
Note: the dimension,
dim, of the array
d
must be at least
.
On entry: the diagonal elements of the tridiagonal matrix .
On exit: the
eigenvalues in descending order, unless
NE_CONVERGENCE or
NE_POS_DEF, in which case
d is overwritten.
- 5:
– doubleInput/Output
-
Note: the dimension,
dim, of the array
e
must be at least
.
On entry: the off-diagonal elements of the tridiagonal matrix .
On exit:
e is overwritten.
- 6:
– ComplexInput/Output
-
Note: the dimension,
dim, of the array
z
must be at least
- when
or ;
- when
.
The
th element of the matrix
is stored in
- when ;
- when .
On entry: if
,
z must contain the unitary matrix
from the reduction to tridiagonal form.
If
,
z need not be set.
On exit: if
or
, the
required orthonormal eigenvectors stored as columns of
; the
th column corresponds to the
th eigenvalue, where
, unless
NE_CONVERGENCE or
NE_POS_DEF.
If
,
z is not referenced.
- 7:
– IntegerInput
-
On entry: the stride separating row or column elements (depending on the value of
order) in the array
z.
Constraints:
- if or , ;
- if , .
- 8:
– NagError *Input/Output
-
The NAG error argument (see
Section 3.6 in the Essential Introduction).
6 Error Indicators and Warnings
- NE_ALLOC_FAIL
-
Dynamic memory allocation failed.
See
Section 3.2.1.2 in the Essential Introduction for further information.
- NE_BAD_PARAM
-
On entry, argument had an illegal value.
- NE_CONVERGENCE
-
The algorithm to compute the singular values of the Cholesky factor failed to converge; off-diagonal elements did not converge to zero.
- NE_ENUM_INT_2
-
On entry, , and .
Constraint: if or , ;
if , .
- NE_INT
-
On entry, .
Constraint: .
On entry, .
Constraint: .
- 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.
An unexpected error has been triggered by this function. Please contact
NAG.
See
Section 3.6.6 in the Essential Introduction for further information.
- NE_NO_LICENCE
-
Your licence key may have expired or may not have been installed correctly.
See
Section 3.6.5 in the Essential Introduction for further information.
- NE_POS_DEF
-
The leading minor of order is not positive definite and the Cholesky factorization of could not be completed. Hence itself is not positive definite.
7 Accuracy
The eigenvalues and eigenvectors of are computed to high relative accuracy which means that if they vary widely in magnitude, then any small eigenvalues (and corresponding eigenvectors) will be computed more accurately than, for example, with the standard method. However, the reduction to tridiagonal form (prior to calling the function) may exclude the possibility of obtaining high relative accuracy in the small eigenvalues of the original matrix if its eigenvalues vary widely in magnitude.
To be more precise, let
be the tridiagonal matrix defined by
, where
is diagonal with
, and
for all
. If
is an exact eigenvalue of
and
is the corresponding computed value, then
where
is a modestly increasing function of
,
is the
machine precision, and
is the condition number of
with respect to inversion defined by:
.
If
is the corresponding exact eigenvector of
, and
is the corresponding computed eigenvector, then the angle
between them is bounded as follows:
where
is the relative gap between
and the other eigenvalues, defined by
8 Parallelism and Performance
nag_zpteqr (f08juc) is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
nag_zpteqr (f08juc) 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.
The total number of real floating-point operations is typically about if and about if or , but depends on how rapidly the algorithm converges. When , the operations are all performed in scalar mode; the additional operations to compute the eigenvectors when or can be vectorized and on some machines may be performed much faster.
The real analogue of this function is
nag_dpteqr (f08jgc).
10 Example
This example computes all the eigenvalues and eigenvectors of the complex Hermitian positive definite matrix
, where
10.1 Program Text
Program Text (f08juce.c)
10.2 Program Data
Program Data (f08juce.d)
10.3 Program Results
Program Results (f08juce.r)