hide long namesshow long names
hide short namesshow short names
Integer type:  int32  int64  nag_int  show int32  show int32  show int64  show int64  show nag_int  show nag_int

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

NAG Toolbox: nag_lapack_zpteqr (f08ju)

Purpose

nag_lapack_zpteqr (f08ju) computes all the eigenvalues and, optionally, all the eigenvectors of a complex Hermitian positive definite matrix which has been reduced to tridiagonal form.

Syntax

[d, e, z, info] = f08ju(compz, d, e, z, 'n', n)
[d, e, z, info] = nag_lapack_zpteqr(compz, d, e, z, 'n', n)

Description

nag_lapack_zpteqr (f08ju) computes all the eigenvalues and, optionally, all the eigenvectors of a real symmetric positive definite tridiagonal matrix TT. In other words, it can compute the spectral factorization of TT as
T = ZΛZT,
T=ZΛZT,
where ΛΛ is a diagonal matrix whose diagonal elements are the eigenvalues λiλi, and ZZ is the orthogonal matrix whose columns are the eigenvectors zizi. Thus
Tzi = λizi,  i = 1,2,,n.
Tzi=λizi,  i=1,2,,n.
The function stores the real orthogonal matrix ZZ in a complex array, so that it may be used to compute all the eigenvalues and eigenvectors of a complex Hermitian positive definite matrix AA which has been reduced to tridiagonal form TT:
A = QTQH, where ​Q​ is unitary
= (QZ)Λ(QZ)H.
A =QTQH, where ​Q​ is unitary =(QZ)Λ(QZ)H.
In this case, the matrix QQ must be formed explicitly and passed to nag_lapack_zpteqr (f08ju), which must be called with compz = 'V'compz='V'. The functions which must be called to perform the reduction to tridiagonal form and form QQ are:
full matrix nag_lapack_zhetrd (f08fs) and nag_lapack_zungtr (f08ft)
full matrix, packed storage nag_lapack_zhptrd (f08gs) and nag_lapack_zupgtr (f08gt)
band matrix nag_lapack_zhbtrd (f08hs) with vect = 'V'vect='V'.
nag_lapack_zpteqr (f08ju) first factorizes TT as LDLHLDLH where LL is unit lower bidiagonal and DD is diagonal. It forms the bidiagonal matrix B = LD(1/2)B=LD12, and then calls nag_lapack_zbdsqr (f08ms) to compute the singular values of BB which are the same as the eigenvalues of TT. The method used by the function allows high relative accuracy to be achieved in the small eigenvalues of TT. The eigenvectors are normalized so that zi2 = 1zi2=1, but are determined only to within a complex factor of absolute value 11.

References

Barlow J and Demmel J W (1990) Computing accurate eigensystems of scaled diagonally dominant matrices SIAM J. Numer. Anal. 27 762–791

Parameters

Compulsory Input Parameters

1:     compz – string (length ≥ 1)
Indicates whether the eigenvectors are to be computed.
compz = 'N'compz='N'
Only the eigenvalues are computed (and the array z is not referenced).
compz = 'I'compz='I'
The eigenvalues and eigenvectors of TT are computed (and the array z is initialized by the function).
compz = 'V'compz='V'
The eigenvalues and eigenvectors of AA are computed (and the array z must contain the matrix QQ on entry).
Constraint: compz = 'N'compz='N', 'V''V' or 'I''I'.
2:     d( : :) – double array
Note: the dimension of the array d must be at least max (1,n)max(1,n).
The diagonal elements of the tridiagonal matrix TT.
3:     e( : :) – double array
Note: the dimension of the array e must be at least max (1,n1)max(1,n-1).
The off-diagonal elements of the tridiagonal matrix TT.
4:     z(ldz, : :) – complex array
The first dimension, ldz, of the array z must satisfy
  • if compz = 'I'compz='I' or 'V''V', ldz max (1,n) ldz max(1,n) ;
  • if compz = 'N'compz='N', ldz1ldz1.
The second dimension of the array must be at least max (1,n)max(1,n) if compz = 'V'compz='V' or 'I''I' and at least 11 if compz = 'N'compz='N'
If compz = 'V'compz='V', z must contain the unitary matrix QQ from the reduction to tridiagonal form.
If compz = 'I'compz='I', z need not be set.

Optional Input Parameters

1:     n – int64int32nag_int scalar
Default: The first dimension of the array d and the second dimension of the array d. (An error is raised if these dimensions are not equal.)
nn, the order of the matrix TT.
Constraint: n0n0.

Input Parameters Omitted from the MATLAB Interface

ldz work

Output Parameters

1:     d( : :) – double array
Note: the dimension of the array d must be at least max (1,n)max(1,n).
The nn eigenvalues in descending order, unless INFO > 0INFO>0, in which case d is overwritten.
2:     e( : :) – double array
Note: the dimension of the array e must be at least max (1,n1)max(1,n-1).
3:     z(ldz, : :) – complex array
The first dimension, ldz, of the array z will be
  • if compz = 'I'compz='I' or 'V''V', ldz max (1,n) ldz max(1,n) ;
  • if compz = 'N'compz='N', ldz1ldz1.
The second dimension of the array will be max (1,n)max(1,n) if compz = 'V'compz='V' or 'I''I' and at least 11 if compz = 'N'compz='N'
If compz = 'I'compz='I' or 'V''V', the nn required orthonormal eigenvectors stored as columns of ZZ; the iith column corresponds to the iith eigenvalue, where i = 1,2,,ni=1,2,,n, unless INFO > 0INFO>0.
If compz = 'N'compz='N', z is not referenced.
4:     info – int64int32nag_int scalar
info = 0info=0 unless the function detects an error (see Section [Error Indicators and Warnings]).

Error Indicators and Warnings

  info = iinfo=-i
If info = iinfo=-i, parameter ii had an illegal value on entry. The parameters are numbered as follows:
1: compz, 2: n, 3: d, 4: e, 5: z, 6: ldz, 7: work, 8: info.
It is possible that info refers to a parameter that is omitted from the MATLAB interface. This usually indicates that an error in one of the other input parameters has caused an incorrect value to be inferred.
  INFO > 0INFO>0
If info = iinfo=i, the leading minor of order ii is not positive definite and the Cholesky factorization of TT could not be completed. Hence TT itself is not positive definite.
If info = n + iinfo=n+i, the algorithm to compute the singular values of the Cholesky factor BB failed to converge; ii off-diagonal elements did not converge to zero.

Accuracy

The eigenvalues and eigenvectors of TT 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 QRQR 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 HH be the tridiagonal matrix defined by H = DTDH=DTD, where DD is diagonal with dii = tii(1/2) dii = t ii -12 , and hii = 1 hii = 1  for all ii. If λiλi is an exact eigenvalue of TT and λ̃iλ~i is the corresponding computed value, then
|λ̃iλi| c (n) ε κ2 (H) λi
| λ~i - λi | c (n) ε κ2 (H) λi
where c(n)c(n) is a modestly increasing function of nn, εε is the machine precision, and κ2(H)κ2(H) is the condition number of HH with respect to inversion defined by: κ2(H) = H · H1κ2(H)=H·H-1.
If zizi is the corresponding exact eigenvector of TT, and iz~i is the corresponding computed eigenvector, then the angle θ(i,zi)θ(z~i,zi) between them is bounded as follows:
θ (i,zi) ( c (n) ε κ2 (H) )/(relgapi)
θ (z~i,zi) c (n) ε κ2 (H) relgapi
where relgapirelgapi is the relative gap between λiλi and the other eigenvalues, defined by
relgapi = min ( |λiλj| )/((λi + λj)).
ij
relgapi = min ij | λi - λj | ( λi + λj ) .

Further Comments

The total number of real floating point operations is typically about 30n230n2 if compz = 'N'compz='N' and about 12n312n3 if compz = 'V'compz='V' or 'I''I', but depends on how rapidly the algorithm converges. When compz = 'N'compz='N', the operations are all performed in scalar mode; the additional operations to compute the eigenvectors when compz = 'V'compz='V' or 'I''I' can be vectorized and on some machines may be performed much faster.
The real analogue of this function is nag_lapack_dpteqr (f08jg).

Example

function nag_lapack_zpteqr_example
compz = 'V';
d = [6.02;
     2.738844788384059;
     5.173556804164482;
     2.467598407451455];
e = [2.74238946905796;
     1.835961995070032;
     1.695211553772095];
z = [complex(1),  0 + 0i,  0 + 0i,  0 + 0i;
      0 + 0i,  -0.1640904784230299 - 0.09116137690168336i, ...
     0.04492226830902458 - 0.1991468061366732i,  -0.7606249187911637 - 0.5869720526411456i;
      0 + 0i,  -0.4740391598887533 - 0.6344831832357161i, ...
     -0.4067593168412005 + 0.4544041694574636i,  0.02193769252276673 + 0.01733238795915084i;
      0 + 0i,  0.5287359860297633 + 0.240666035020444i, ...
    -0.1787167294506699 + 0.7446116967739244i,  -0.2225496702687938 - 0.1631058324212738i];
[dOut, eOut, zOut, info] = nag_lapack_zpteqr(compz, d, e, z)
 

dOut =

    7.9995
    5.9976
    2.0003
    0.4026


eOut =

     0
     0
     0


zOut =

   0.7289 + 0.0000i  -0.5130 + 0.0000i   0.2606 + 0.0000i  -0.3709 + 0.0000i
  -0.1651 - 0.2067i  -0.2486 - 0.3726i  -0.5981 - 0.4200i  -0.4009 - 0.1860i
  -0.4170 - 0.1413i  -0.3086 + 0.3554i   0.2957 + 0.1501i  -0.1848 - 0.6637i
   0.1748 + 0.4175i  -0.2188 + 0.5166i  -0.3501 - 0.4068i   0.4001 - 0.1798i


info =

                    0


function f08ju_example
compz = 'V';
d = [6.02;
     2.738844788384059;
     5.173556804164482;
     2.467598407451455];
e = [2.74238946905796;
     1.835961995070032;
     1.695211553772095];
z = [complex(1),  0 + 0i,  0 + 0i,  0 + 0i;
      0 + 0i,  -0.1640904784230299 - 0.09116137690168336i, ...
     0.04492226830902458 - 0.1991468061366732i,  -0.7606249187911637 - 0.5869720526411456i;
      0 + 0i,  -0.4740391598887533 - 0.6344831832357161i, ...
     -0.4067593168412005 + 0.4544041694574636i,  0.02193769252276673 + 0.01733238795915084i;
      0 + 0i,  0.5287359860297633 + 0.240666035020444i, ...
    -0.1787167294506699 + 0.7446116967739244i,  -0.2225496702687938 - 0.1631058324212738i];
[dOut, eOut, zOut, info] = f08ju(compz, d, e, z)
 

dOut =

    7.9995
    5.9976
    2.0003
    0.4026


eOut =

     0
     0
     0


zOut =

   0.7289 + 0.0000i  -0.5130 + 0.0000i   0.2606 + 0.0000i  -0.3709 + 0.0000i
  -0.1651 - 0.2067i  -0.2486 - 0.3726i  -0.5981 - 0.4200i  -0.4009 - 0.1860i
  -0.4170 - 0.1413i  -0.3086 + 0.3554i   0.2957 + 0.1501i  -0.1848 - 0.6637i
   0.1748 + 0.4175i  -0.2188 + 0.5166i  -0.3501 - 0.4068i   0.4001 - 0.1798i


info =

                    0



PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2013