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_zhegvd (f08sq)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_lapack_zhegvd (f08sq) computes all the eigenvalues and, optionally, the eigenvectors of a complex generalized Hermitian-definite eigenproblem, of the form
Az=λBz ,   ABz=λz   or   BAz=λz ,  
where A and B are Hermitian and B is also positive definite. If eigenvectors are desired, it uses a divide-and-conquer algorithm.

Syntax

[a, b, w, info] = f08sq(itype, jobz, uplo, a, b, 'n', n)
[a, b, w, info] = nag_lapack_zhegvd(itype, jobz, uplo, a, b, 'n', n)

Description

nag_lapack_zhegvd (f08sq) first performs a Cholesky factorization of the matrix B as B=UHU , when uplo='U' or B=LLH , when uplo='L'. The generalized problem is then reduced to a standard symmetric eigenvalue problem
Cx=λx ,  
which is solved for the eigenvalues and, optionally, the eigenvectors; the eigenvectors are then backtransformed to give the eigenvectors of the original problem.
For the problem Az=λBz , the eigenvectors are normalized so that the matrix of eigenvectors, z, satisfies
ZH A Z = Λ   and   ZH B Z = I ,  
where Λ  is the diagonal matrix whose diagonal elements are the eigenvalues. For the problem A B z = λ z  we correspondingly have
Z-1 A Z-H = Λ   and   ZH B Z = I ,  
and for B A z = λ z  we have
ZH A Z = Λ   and   ZH B-1 Z = I .  

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 http://www.netlib.org/lapack/lug
Golub G H and Van Loan C F (1996) Matrix Computations (3rd Edition) Johns Hopkins University Press, Baltimore

Parameters

Compulsory Input Parameters

1:     itype int64int32nag_int scalar
Specifies the problem type to be solved.
itype=1
Az=λBz.
itype=2
ABz=λz.
itype=3
BAz=λz.
Constraint: itype=1, 2 or 3.
2:     jobz – string (length ≥ 1)
Indicates whether eigenvectors are computed.
jobz='N'
Only eigenvalues are computed.
jobz='V'
Eigenvalues and eigenvectors are computed.
Constraint: jobz='N' or 'V'.
3:     uplo – string (length ≥ 1)
If uplo='U', the upper triangles of A and B are stored.
If uplo='L', the lower triangles of A and B are stored.
Constraint: uplo='U' or 'L'.
4:     alda: – complex array
The first dimension of the array a must be at least max1,n.
The second dimension of the array a must be at least max1,n.
The n by 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.
5:     bldb: – complex array
The first dimension of the array b must be at least max1,n.
The second dimension of the array b must be at least max1,n.
The n by n Hermitian matrix B.
  • If uplo='U', the upper triangular part of b must be stored and the elements of the array below the diagonal are not referenced.
  • If uplo='L', the lower triangular part of b must be stored and the elements of the array above the diagonal are not referenced.

Optional Input Parameters

1:     n int64int32nag_int scalar
Default: the first dimension of the arrays a, b and the second dimension of the arrays a, b. (An error is raised if these dimensions are not equal.)
n, the order of the matrices A and B.
Constraint: n0.

Output Parameters

1:     alda: – complex array
The first dimension of the array a will be max1,n.
The second dimension of the array a will be max1,n.
If jobz='V', a contains the matrix Z of eigenvectors. The eigenvectors are normalized as follows:
  • if itype=1 or 2, ZHBZ=I;
  • if itype=3, ZHB-1Z=I.
If jobz='N', the upper triangle (if uplo='U') or the lower triangle (if uplo='L') of a, including the diagonal, is overwritten.
2:     bldb: – complex array
The first dimension of the array b will be max1,n.
The second dimension of the array b will be max1,n.
The triangular factor U or L from the Cholesky factorization B=UHU or B=LLH.
3:     wn – double array
The eigenvalues in ascending order.
4:     info int64int32nag_int scalar
info=0 unless the function detects an error (see Error Indicators and Warnings).

Error Indicators and Warnings

   info=-i
If info=-i, parameter i had an illegal value on entry. The parameters are numbered as follows:
1: itype, 2: jobz, 3: uplo, 4: n, 5: a, 6: lda, 7: b, 8: ldb, 9: w, 10: work, 11: lwork, 12: rwork, 13: lrwork, 14: iwork, 15: liwork, 16: 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=1ton
If info=i, nag_lapack_zheevd (f08fq) failed to converge; i i off-diagonal elements of an intermediate tridiagonal form did not converge to zero.
   info>n
nag_lapack_zpotrf (f07fr) returned an error code; i.e., if info=n+i, for 1in, then the leading minor of order i of B is not positive definite. The factorization of B could not be completed and no eigenvalues or eigenvectors were computed.

Accuracy

If B is ill-conditioned with respect to inversion, then the error bounds for the computed eigenvalues and vectors may be large, although when the diagonal elements of B differ widely in magnitude the eigenvalues and eigenvectors may be less sensitive than the condition of B would suggest. See Section 4.10 of Anderson et al. (1999) for details of the error bounds.
The example program below illustrates the computation of approximate error bounds.

Further Comments

The total number of floating-point operations is proportional to n3 .
The real analogue of this function is nag_lapack_dsygvd (f08sc).

Example

This example finds all the eigenvalues and eigenvectors of the generalized Hermitian eigenproblem ABz=λz , where
A = -7.36i+0.00 0.77-0.43i -0.64-0.92i 3.01-6.97i 0.77+0.43i 3.49i+0.00 2.19+4.45i 1.90+3.73i -0.64+0.92i 2.19-4.45i 0.12i+0.00 2.88-3.17i 3.01+6.97i 1.90-3.73i 2.88+3.17i -2.54i+0.00  
and
B = 3.23i+0.00 1.51-1.92i 1.90+0.84i 0.42+2.50i 1.51+1.92i 3.58i+0.00 -0.23+1.11i -1.18+1.37i 1.90-0.84i -0.23-1.11i 4.09i+0.00 2.33-0.14i 0.42-2.50i -1.18-1.37i 2.33+0.14i 4.29i+0.00 ,  
together with an estimate of the condition number of B, and approximate error bounds for the computed eigenvalues and eigenvectors.
The example program for nag_lapack_zhegv (f08sn) illustrates solving a generalized Hermitian eigenproblem of the form Az=λBz .
function f08sq_example


fprintf('f08sq example results\n\n');

% Upper triangular parts of Hermitian matrix A and positive definite matrix B
uplo = 'Upper';
n = int64(4);
a = [-7.36 + 0i,  0.77 - 0.43i, -0.64 - 0.92i,  3.01 - 6.97i;
      0    + 0i,  3.49 + 0i,     2.19 + 4.45i,  1.90 + 3.73i;
      0    + 0i,  0    + 0i,     0.12 + 0i,     2.88 - 3.17i;
      0    + 0i,  0    + 0i,     0    + 0i,    -2.54 + 0i];
b = [ 3.23 + 0i,  1.51 - 1.92i,  1.90 + 0.84i,  0.42 + 2.50i;
      0    + 0i,  3.58 + 0i,    -0.23 + 1.11i, -1.18 + 1.37i;
      0    + 0i,  0    + 0i,     4.09 + 0i,     2.33 - 0.14i;
      0    + 0i,  0    + 0i,     0    + 0i,     4.29 + 0i];

anorm = norm(a+a'-diag(diag(a)),1);
bnorm = norm(b+b'-diag(diag(b)),1);

% Generalized eigenvalues and eigenvectors for problem ABz = lambda z.
itype = int64(2);
jobz = 'Vectors';
[Z, L, w, info] = f08sq( ...
                         itype, jobz, uplo, a, b);

% Normalize: largest elements are real (preserving Z^HBZ = I)
for i = 1:n
  [~,k] = max(abs(real(Z(:,i)))+abs(imag(Z(:,i))));
  Z(:,i) = Z(:,i)*conj(Z(k,i))/abs(Z(k,i));
end

disp('Eigenvalues');
disp(w');
disp('Eigenvectors');
disp(Z);

% Get reciprocal condition number
[rcond, info] = f07tu( ...
                       '1', uplo, 'N', L);

rcondb = rcond^2;
fprintf('Estimate of reciprocal condition number for B\n%12.1e\n', rcondb);

% Eigenvector condition numbers
[rcondz, info] = f08fl( ...
		        'Eigenvectors', n, n, w);

% Error bounds
t1 = 1/rcond;
t2 = x02aj*t1;
t3 = anorm*bnorm;
errbnd(1:n) = x02aj*t3 + (x02aj/rcondb).*abs(w);
zerrbd(1:n) = t1*t2 + (t2*t3)./rcondz;

fprintf('\n');
disp('Error estimate for the eigenvalues');
fprintf('%12.1e',errbnd);
fprintf('\n\n');
disp('Error estimates for the eigenvectors');
fprintf('%12.1e',zerrbd);
fprintf('\n');


f08sq example results

Eigenvalues
  -61.7321   -6.6195    0.0725   43.1883

Eigenvectors
   0.3903 + 0.0000i  -0.1560 - 0.0404i   2.2909 + 0.0000i  -0.1943 - 0.0690i
  -0.1814 + 0.0114i  -0.1552 - 0.3651i  -0.5042 - 0.7120i   0.3884 + 0.0000i
   0.0438 + 0.0338i   0.5364 + 0.0000i  -1.2701 - 0.4547i   0.0657 - 0.2095i
  -0.2221 - 0.2272i  -0.1298 - 0.1880i   0.5706 + 1.3132i   0.2924 - 0.0675i

Estimate of reciprocal condition number for B
     2.5e-03

Error estimate for the eigenvalues
     2.7e-12     3.1e-13     2.6e-14     1.9e-12

Error estimates for the eigenvectors
     5.2e-14     1.1e-13     1.1e-13     5.4e-14

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–2015