nag_dpteqr (f08jgc) (PDF version)
f08 Chapter Contents
f08 Chapter Introduction
NAG Library Manual

NAG Library Function Document

nag_dpteqr (f08jgc)

+ Contents

    1  Purpose
    7  Accuracy

1  Purpose

nag_dpteqr (f08jgc) computes all the eigenvalues and, optionally, all the eigenvectors of a real symmetric positive definite tridiagonal matrix, or of a real symmetric positive definite matrix which has been reduced to tridiagonal form.

2  Specification

#include <nag.h>
#include <nagf08.h>
void  nag_dpteqr (Nag_OrderType order, Nag_ComputeZType compz, Integer n, double d[], double e[], double z[], Integer pdz, NagError *fail)

3  Description

nag_dpteqr (f08jgc) computes all the eigenvalues and, optionally, all the eigenvectors of a real symmetric positive definite 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 may also be used to compute all the eigenvalues and eigenvectors of a real symmetric positive definite matrix A which has been reduced to tridiagonal form T:
A =QTQT, where ​Q​ is orthogonal =QZΛQZT.
In this case, the matrix Q must be formed explicitly and passed to nag_dpteqr (f08jgc), 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 nag_dsytrd (f08fec) and nag_dorgtr (f08ffc)
full matrix, packed storage nag_dsptrd (f08gec) and nag_dopgtr (f08gfc)
band matrix nag_dsbtrd (f08hec) with vect=Nag_FormQ.
nag_dpteqr (f08jgc) first factorizes T as LDLT where L is unit lower bidiagonal and D is diagonal. It forms the bidiagonal matrix B=LD12, and then calls nag_dbdsqr (f08mec) to compute the singular values of B which are the same as the eigenvalues of T. The method used by the function allows high relative accuracy to be achieved in the small eigenvalues of T. The eigenvectors are normalized so that zi2=1, but are determined only to within a factor ±1.

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:     orderNag_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 order=Nag_RowMajor. See Section 3.2.1.3 in the Essential Introduction for a more detailed explanation of the use of this argument.
Constraint: order=Nag_RowMajor or Nag_ColMajor.
2:     compzNag_ComputeZTypeInput
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_InitZ
The eigenvalues and eigenvectors of T are computed (and the array z is initialized by the function).
compz=Nag_UpdateZ
The eigenvalues and eigenvectors of A are computed (and the array z must contain the matrix Q on entry).
Constraint: compz=Nag_NotZ, Nag_UpdateZ or Nag_InitZ.
3:     nIntegerInput
On entry: n, the order of the matrix T.
Constraint: n0.
4:     d[dim]doubleInput/Output
Note: the dimension, dim, of the array d must be at least max1,n.
On entry: the diagonal elements of the tridiagonal matrix T.
On exit: the n eigenvalues in descending order, unless fail.code= NE_CONVERGENCE or NE_POS_DEF, in which case d is overwritten.
5:     e[dim]doubleInput/Output
Note: the dimension, dim, of the array e must be at least max1,n-1.
On entry: the off-diagonal elements of the tridiagonal matrix T.
On exit: e is overwritten.
6:     z[dim]doubleInput/Output
Note: the dimension, dim, of the array z must be at least
  • max1,pdz×n when compz=Nag_UpdateZ or Nag_InitZ and order=Nag_ColMajor;
  • max1,×pdz when compz=Nag_UpdateZ or Nag_InitZ and order=Nag_RowMajor;
  • 1 when compz=Nag_NotZ.
The i,jth 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 orthogonal matrix Q from the reduction to tridiagonal form.
If compz=Nag_InitZ, z need not be set.
On exit: if compz=Nag_InitZ or Nag_UpdateZ, 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 or NE_POS_DEF.
If compz=Nag_NotZ, z is not referenced.
7:     pdzIntegerInput
On entry: the stride separating row or column elements (depending on the value of order) in the array z.
Constraints:
  • if order=Nag_ColMajor,
    • if compz=Nag_InitZ or Nag_UpdateZ, pdz max1,n ;
    • if compz=Nag_NotZ, pdz1;
  • if order=Nag_RowMajor,
    • if compz=Nag_UpdateZ or Nag_InitZ, pdzmax1,n;
    • if compz=Nag_NotZ, pdz1.
8:     failNagError *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.
NE_BAD_PARAM
On entry, argument value had an illegal value.
NE_CONVERGENCE
The algorithm to compute the singular values of the Cholesky factor B failed to converge; value off-diagonal elements did not converge to zero.
NE_ENUM_INT_2
On entry, compz=value, pdz=value and n=value.
Constraint: if compz=Nag_InitZ or Nag_UpdateZ, pdz max1,n ;
if compz=Nag_NotZ, pdz1.
On entry, compz=value, pdz=value, n=value.
Constraint: if compz=Nag_UpdateZ or Nag_InitZ, pdzmax1,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.
NE_POS_DEF
The leading minor of order value is not positive definite and the Cholesky factorization of T could not be completed. Hence T itself is not positive definite.

7  Accuracy

The eigenvalues and eigenvectors of T 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 QR 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 H be the tridiagonal matrix defined by H=DTD, where D is diagonal with dii = t ii -12 , and hii = 1  for all i. If λi is an exact eigenvalue of T and λ~i is the corresponding computed value, then
λ~i - λi c n ε κ2 H λi
where cn is a modestly increasing function of n, ε is the machine precision, and κ2H is the condition number of H with respect to inversion defined by: κ2H=H·H-1.
If zi is the corresponding exact eigenvector of T, 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 ε κ2 H relgapi
where relgapi is the relative gap between λi and the other eigenvalues, defined by
relgapi = min ij λi - λj λi + λj .

8  Parallelism and Performance

nag_dpteqr (f08jgc) is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
nag_dpteqr (f08jgc) 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 Users' Note for your implementation for any additional implementation-specific information.

9  Further Comments

The total number of floating-point operations is typically about 30n2 if compz=Nag_NotZ and about 6n3 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 complex analogue of this function is nag_zpteqr (f08juc).

10  Example

This example computes all the eigenvalues and eigenvectors of the symmetric positive definite tridiagonal matrix T, where
T = 4.16 3.17 0.00 0.00 3.17 5.25 -0.97 0.00 0.00 -0.97 1.09 0.55 0.00 0.00 0.55 0.62 .

10.1  Program Text

Program Text (f08jgce.c)

10.2  Program Data

Program Data (f08jgce.d)

10.3  Program Results

Program Results (f08jgce.r)


nag_dpteqr (f08jgc) (PDF version)
f08 Chapter Contents
f08 Chapter Introduction
NAG Library Manual

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