NAG Library Routine Document
F08UPF (ZHBGVX)
1 Purpose
F08UPF (ZHBGVX) computes selected the eigenvalues and, optionally, the eigenvectors of a complex generalized Hermitian-definite banded eigenproblem, of the form
where
and
are Hermitian and banded, and
is also positive definite. Eigenvalues and eigenvectors can be selected by specifying either all eigenvalues, a range of values or a range of indices for the desired eigenvalues.
2 Specification
SUBROUTINE F08UPF ( |
JOBZ, RANGE, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, RWORK, IWORK, JFAIL, INFO) |
INTEGER |
N, KA, KB, LDAB, LDBB, LDQ, IL, IU, M, LDZ, IWORK(5*N), JFAIL(*), INFO |
REAL (KIND=nag_wp) |
VL, VU, ABSTOL, W(N), RWORK(7*N) |
COMPLEX (KIND=nag_wp) |
AB(LDAB,*), BB(LDBB,*), Q(LDQ,*), Z(LDZ,*), WORK(N) |
CHARACTER(1) |
JOBZ, RANGE, UPLO |
|
The routine may be called by its
LAPACK
name zhbgvx.
3 Description
The generalized Hermitian-definite band problem
is first reduced to a standard band Hermitian problem
where
is a Hermitian band matrix, using Wilkinson's modification to Crawford's algorithm (see
Crawford (1973) and
Wilkinson (1977)). The Hermitian eigenvalue problem is then solved for the required eigenvalues and eigenvectors, and the eigenvectors are then backtransformed to the eigenvectors of the original problem.
The eigenvectors are normalized so that
4 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
Crawford C R (1973) Reduction of a band-symmetric generalized eigenvalue problem Comm. ACM 16 41–44
Demmel J W and Kahan W (1990) Accurate singular values of bidiagonal matrices SIAM J. Sci. Statist. Comput. 11 873–912
Golub G H and Van Loan C F (1996) Matrix Computations (3rd Edition) Johns Hopkins University Press, Baltimore
Wilkinson J H (1977) Some recent advances in numerical linear algebra The State of the Art in Numerical Analysis (ed D A H Jacobs) Academic Press
5 Parameters
- 1: – CHARACTER(1)Input
-
On entry: indicates whether eigenvectors are computed.
- Only eigenvalues are computed.
- Eigenvalues and eigenvectors are computed.
Constraint:
or .
- 2: – CHARACTER(1)Input
-
On entry: 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: – CHARACTER(1)Input
-
On entry: if
, the upper triangles of
and
are stored.
If , the lower triangles of and are stored.
Constraint:
or .
- 4: – INTEGERInput
-
On entry: , the order of the matrices and .
Constraint:
.
- 5: – INTEGERInput
-
On entry: if
, the number of superdiagonals,
, of the matrix
.
If , the number of subdiagonals, , of the matrix .
Constraint:
.
- 6: – INTEGERInput
-
On entry: if
, the number of superdiagonals,
, of the matrix
.
If , the number of subdiagonals, , of the matrix .
Constraint:
.
- 7: – COMPLEX (KIND=nag_wp) arrayInput/Output
-
Note: the second dimension of the array
AB
must be at least
.
On entry: the upper or lower triangle of the
by
Hermitian band matrix
.
The matrix is stored in rows
to
, more precisely,
- if , the elements of the upper triangle of within the band must be stored with element in ;
- if , the elements of the lower triangle of within the band must be stored with element in
On exit: the contents of
AB are overwritten.
- 8: – INTEGERInput
-
On entry: the first dimension of the array
AB as declared in the (sub)program from which F08UPF (ZHBGVX) is called.
Constraint:
.
- 9: – COMPLEX (KIND=nag_wp) arrayInput/Output
-
Note: the second dimension of the array
BB
must be at least
.
On entry: the upper or lower triangle of the
by
Hermitian positive definite band matrix
.
The matrix is stored in rows
to
, more precisely,
- if , the elements of the upper triangle of within the band must be stored with element in ;
- if , the elements of the lower triangle of within the band must be stored with element in
On exit: the factor
from the split Cholesky factorization
, as returned by
F08UTF (ZPBSTF).
- 10: – INTEGERInput
-
On entry: the first dimension of the array
BB as declared in the (sub)program from which F08UPF (ZHBGVX) is called.
Constraint:
.
- 11: – COMPLEX (KIND=nag_wp) arrayOutput
-
Note: the second dimension of the array
Q
must be at least
if
, and at least
otherwise.
On exit: if
, the
by
matrix,
used in the reduction of the standard form, i.e.,
, from symmetric banded to tridiagonal form.
If
,
Q is not referenced.
- 12: – INTEGERInput
-
On entry: the first dimension of the array
Q as declared in the (sub)program from which F08UPF (ZHBGVX) is called.
Constraints:
- if , ;
- otherwise .
- 13: – REAL (KIND=nag_wp)Input
- 14: – REAL (KIND=nag_wp)Input
-
On entry: if
, the lower and upper bounds of the interval to be searched for eigenvalues.
If
or
,
VL and
VU are not referenced.
Constraint:
if , .
- 15: – INTEGERInput
- 16: – INTEGERInput
-
On entry: 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 , .
- 17: – REAL (KIND=nag_wp)Input
-
On entry: 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, where
is the tridiagonal matrix obtained by reducing
to tridiagonal form. Eigenvalues will be computed most accurately when
ABSTOL is set to twice the underflow threshold
, not zero. If this routine returns with
, indicating that some eigenvectors did not converge, try setting
ABSTOL to
. See
Demmel and Kahan (1990).
- 18: – INTEGEROutput
-
On exit: the total number of eigenvalues found.
.
If , .
If , .
- 19: – REAL (KIND=nag_wp) arrayOutput
-
On exit: the eigenvalues in ascending order.
- 20: – COMPLEX (KIND=nag_wp) arrayOutput
-
Note: the second dimension of the array
Z
must be at least
if
, and at least
otherwise.
On exit: if
,
Z contains the matrix
of eigenvectors, with the
th column of
holding the eigenvector associated with
. The eigenvectors are normalized so that
.
If
,
Z is not referenced.
- 21: – INTEGERInput
-
On entry: the first dimension of the array
Z as declared in the (sub)program from which F08UPF (ZHBGVX) is called.
Constraints:
- if , ;
- otherwise .
- 22: – COMPLEX (KIND=nag_wp) arrayWorkspace
-
- 23: – REAL (KIND=nag_wp) arrayWorkspace
-
- 24: – INTEGER arrayWorkspace
-
- 25: – INTEGER arrayOutput
-
Note: the dimension of the array
JFAIL
must be at least
.
On exit: if
, then
- if , the first M elements of JFAIL are zero;
- if , JFAIL contains the indices of the eigenvectors that failed to converge.
If
,
JFAIL is not referenced.
- 26: – INTEGEROutput
On exit:
unless the routine detects an error (see
Section 6).
6 Error Indicators and Warnings
-
If , argument had an illegal value. An explanatory message is output, and execution of the program is terminated.
-
If
, then
eigenvectors failed to converge. Their indices are stored in array
JFAIL. Please see
ABSTOL.
-
F08UFF (DPBSTF) 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.
7 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.
8 Parallelism and Performance
F08UPF (ZHBGVX) is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
F08UPF (ZHBGVX) 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
X06 Chapter Introduction for information on how to control and interrogate the OpenMP environment used within this routine. Please also consult the
Users' Note for your implementation for any additional implementation-specific information.
The total number of floating-point operations is proportional to if and , and assuming that , is approximately proportional to if . Otherwise the number of floating-point operations depends upon the number of eigenvectors computed.
The real analogue of this routine is
F08UBF (DSBGVX).
10 Example
This example finds the eigenvalues in the half-open interval
, and corresponding eigenvectors, of the generalized band Hermitian eigenproblem
, where
and
10.1 Program Text
Program Text (f08upfe.f90)
10.2 Program Data
Program Data (f08upfe.d)
10.3 Program Results
Program Results (f08upfe.r)