PDF version (NAG web site
, 64-bit version, 64-bit version)
NAG Toolbox: nag_lapack_dstevr (f08jd)
Purpose
nag_lapack_dstevr (f08jd) computes selected eigenvalues and, optionally, eigenvectors of a real by symmetric tridiagonal matrix . Eigenvalues and eigenvectors can be selected by specifying either a range of values or a range of indices for the desired eigenvalues.
Syntax
[
d,
e,
m,
w,
z,
isuppz,
info] = f08jd(
jobz,
range,
d,
e,
vl,
vu,
il,
iu,
abstol, 'n',
n)
[
d,
e,
m,
w,
z,
isuppz,
info] = nag_lapack_dstevr(
jobz,
range,
d,
e,
vl,
vu,
il,
iu,
abstol, 'n',
n)
Description
Whenever possible
nag_lapack_dstevr (f08jd) computes the eigenspectrum using Relatively Robust Representations.
nag_lapack_dstevr (f08jd) computes eigenvalues by the dqds algorithm, while orthogonal eigenvectors are computed from various ‘good’
representations (also known as Relatively Robust Representations). Gram–Schmidt orthogonalization is avoided as far as possible. More specifically, the various steps of the algorithm are as follows. For the
th unreduced block of
:
(a) |
compute
, such that
is a relatively robust representation, |
(b) |
compute the eigenvalues, , of
to high relative accuracy by the dqds algorithm, |
(c) |
if there is a cluster of close eigenvalues, ‘choose’ close to the cluster, and go to (a), |
(d) |
given the approximate eigenvalue of
, compute the corresponding eigenvector by forming a rank-revealing twisted factorization. |
The desired accuracy of the output can be specified by the argument
abstol. For more details, see
Dhillon (1997) and
Parlett and Dhillon (2000).
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
Barlow J and Demmel J W (1990) Computing accurate eigensystems of scaled diagonally dominant matrices SIAM J. Numer. Anal. 27 762–791
Demmel J W and Kahan W (1990) Accurate singular values of bidiagonal matrices SIAM J. Sci. Statist. Comput. 11 873–912
Dhillon I (1997) A new algorithm for the symmetric tridiagonal eigenvalue/eigenvector problem Computer Science Division Technical Report No. UCB//CSD-97-971 UC Berkeley
Golub G H and Van Loan C F (1996) Matrix Computations (3rd Edition) Johns Hopkins University Press, Baltimore
Parlett B N and Dhillon I S (2000) Relatively robust representations of symmetric tridiagonals Linear Algebra Appl. 309 121–151
Parameters
Compulsory Input Parameters
- 1:
– string (length ≥ 1)
-
Indicates whether eigenvectors are computed.
- Only eigenvalues are computed.
- Eigenvalues and eigenvectors are computed.
Constraint:
or .
- 2:
– string (length ≥ 1)
-
If
, all eigenvalues will be found.
If , all eigenvalues in the half-open interval will be found.
If
, the
ilth to
iuth eigenvalues will be found.
Constraint:
, or .
- 3:
– double array
-
The dimension of the array
d
must be at least
The diagonal elements of the tridiagonal matrix .
- 4:
– double array
-
The dimension of the array
e
must be at least
The subdiagonal elements of the tridiagonal matrix .
- 5:
– double scalar
- 6:
– double scalar
-
If
, the lower and upper bounds of the interval to be searched for eigenvalues.
If
or
,
vl and
vu are not referenced.
Constraint:
if , .
- 7:
– int64int32nag_int scalar
- 8:
– int64int32nag_int scalar
-
If
, the indices (in ascending order) of the smallest and largest eigenvalues to be returned.
If
or
,
il and
iu are not referenced.
Constraints:
- if and , and ;
- if and , .
- 9:
– double scalar
-
The absolute error tolerance for the eigenvalues. An approximate eigenvalue is accepted as converged when it is determined to lie in an interval
of width less than or equal to
where
is the
machine precision. If
abstol is less than or equal to zero, then
will be used in its place. See
Demmel and Kahan (1990).
If high relative accuracy is important, set
abstol to
, although doing so does not currently guarantee that eigenvalues are computed to high relative accuracy. See
Barlow and Demmel (1990) for a discussion of which matrices can define their eigenvalues to high relative accuracy.
Optional Input Parameters
- 1:
– int64int32nag_int scalar
-
Default:
the dimension of the array
d.
, the order of the matrix.
Constraint:
.
Output Parameters
- 1:
– double array
-
The dimension of the array
d will be
May be multiplied by a constant factor chosen to avoid over/underflow in computing the eigenvalues.
- 2:
– double array
-
The dimension of the array
e will be
May be multiplied by a constant factor chosen to avoid over/underflow in computing the eigenvalues.
- 3:
– int64int32nag_int scalar
-
The total number of eigenvalues found.
.
If , .
If , .
- 4:
– double array
-
The dimension of the array
w will be
The first
m elements contain the selected eigenvalues in ascending order.
- 5:
– double array
-
The first dimension,
, of the array
z will be
- if , ;
- otherwise .
The second dimension of the array
z will be
if
and
otherwise.
If
, the first
m columns of
contain the orthonormal eigenvectors of the matrix
corresponding to the selected eigenvalues, with the
th column of
holding the eigenvector associated with
.
If
,
z is not referenced.
- 6:
– int64int32nag_int array
-
The dimension of the array
isuppz will be
The support of the eigenvectors in
z, i.e., the indices indicating the nonzero elements in
z. The
th eigenvector is nonzero only in elements
through
. Implemented only for
or
and
.
- 7:
– 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:
jobz, 2:
range, 3:
n, 4:
d, 5:
e, 6:
vl, 7:
vu, 8:
il, 9:
iu, 10:
abstol, 11:
m, 12:
w, 13:
z, 14:
ldz, 15:
isuppz, 16:
work, 17:
lwork, 18:
iwork, 19:
liwork, 20:
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.
-
-
An internal error has occurred in this function. Please refer to
info in
nag_lapack_dstebz (f08jj).
Accuracy
The computed eigenvalues and eigenvectors are exact for a nearby matrix
, where
and
is the
machine precision. See Section 4.7 of
Anderson et al. (1999) for further details.
Further Comments
The total number of floating-point operations is proportional to if and is proportional to if and , otherwise the number of floating-point operations will depend upon the number of computed eigenvectors.
Example
This example finds the eigenvalues with indices in the range
, and the corresponding eigenvectors, of the symmetric tridiagonal matrix
Open in the MATLAB editor:
f08jd_example
function f08jd_example
fprintf('f08jd example results\n\n');
d = [1; 4; 9; 16];
e = [1; 2; 3];
jobz = 'Vectors';
range = 'Indices';
vl = 0;
vu = 0;
il = int64(2);
iu = int64(3);
abstol = 0;
[~, ~, m, w, z, ~, info] = ...
f08jd(...
jobz, range, d, e, vl, vu, il, iu, abstol);
for j = 1:m
[~,k] = max(abs(z(:,j)));
if z(k,j) < 0;
z(:,j) = -z(:,j);
end
end
disp(' Eigenvalues of A in range:');
disp(w(1:m)');
disp(' Corresponding eigenvectors:');
disp(z);
f08jd example results
Eigenvalues of A in range:
3.5470 8.6578
Corresponding eigenvectors:
0.3388 0.0494
0.8628 0.3781
-0.3648 0.8558
0.0879 -0.3497
PDF version (NAG web site
, 64-bit version, 64-bit version)
© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015