PDF version (NAG web site
, 64-bit version, 64-bit version)
NAG Toolbox: nag_lapack_zhegvd (f08sq)
Purpose
nag_lapack_zhegvd (f08sq) computes all the eigenvalues and, optionally, the eigenvectors of a complex generalized Hermitian-definite eigenproblem, of the form
where
and
are Hermitian and
is also positive definite. If eigenvectors are desired, it uses a divide-and-conquer algorithm.
Syntax
Description
nag_lapack_zhegvd (f08sq) first performs a Cholesky factorization of the matrix
as
, when
or
, when
. The generalized problem is then reduced to a standard symmetric eigenvalue problem
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
, the eigenvectors are normalized so that the matrix of eigenvectors,
, satisfies
where
is the diagonal matrix whose diagonal elements are the eigenvalues. For the problem
we correspondingly have
and for
we have
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:
– int64int32nag_int scalar
-
Specifies the problem type to be solved.
- .
- .
- .
Constraint:
, or .
- 2:
– string (length ≥ 1)
-
Indicates whether eigenvectors are computed.
- Only eigenvalues are computed.
- Eigenvalues and eigenvectors are computed.
Constraint:
or .
- 3:
– string (length ≥ 1)
-
If
, the upper triangles of
and
are stored.
If , the lower triangles of and are stored.
Constraint:
or .
- 4:
– complex array
-
The first dimension of the array
a must be at least
.
The second dimension of the array
a must be at least
.
The
by
Hermitian matrix
.
- If , the upper triangular part of must be stored and the elements of the array below the diagonal are not referenced.
- If , the lower triangular part of must be stored and the elements of the array above the diagonal are not referenced.
- 5:
– complex array
-
The first dimension of the array
b must be at least
.
The second dimension of the array
b must be at least
.
The
by
Hermitian matrix
.
- If , the upper triangular part of must be stored and the elements of the array below the diagonal are not referenced.
- If , the lower triangular part of must be stored and the elements of the array above the diagonal are not referenced.
Optional Input Parameters
- 1:
– 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.)
, the order of the matrices and .
Constraint:
.
Output Parameters
- 1:
– complex array
-
The first dimension of the array
a will be
.
The second dimension of the array
a will be
.
If
,
a contains the matrix
of eigenvectors. The eigenvectors are normalized as follows:
- if or , ;
- if , .
If
, the upper triangle (if
) or the lower triangle (if
) of
a, including the diagonal, is overwritten.
- 2:
– complex array
-
The first dimension of the array
b will be
.
The second dimension of the array
b will be
.
The triangular factor or from the Cholesky factorization or .
- 3:
– double array
-
The eigenvalues in ascending order.
- 4:
– int64int32nag_int scalar
unless the function detects an error (see
Error Indicators and Warnings).
Error Indicators and Warnings
-
If , parameter 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.
-
-
If
,
nag_lapack_zheevd (f08fq) failed to converge;
off-diagonal elements of an intermediate tridiagonal form did not converge to zero.
-
-
nag_lapack_zpotrf (f07fr) returned an error code; i.e., if
, for
, then the leading minor of order
of
is not positive definite. The factorization of
could not be completed and no eigenvalues or eigenvectors were computed.
Accuracy
If
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
differ widely in magnitude the eigenvalues and eigenvectors may be less sensitive than the condition of
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 .
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
, where
and
together with an estimate of the condition number of
, 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
.
Open in the MATLAB editor:
f08sq_example
function f08sq_example
fprintf('f08sq example results\n\n');
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);
itype = int64(2);
jobz = 'Vectors';
[Z, L, w, info] = f08sq( ...
itype, jobz, uplo, a, b);
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);
[rcond, info] = f07tu( ...
'1', uplo, 'N', L);
rcondb = rcond^2;
fprintf('Estimate of reciprocal condition number for B\n%12.1e\n', rcondb);
[rcondz, info] = f08fl( ...
'Eigenvectors', n, n, w);
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)
© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015