NAG Library Routine Document
f08ucf (dsbgvd)
1
Purpose
f08ucf (dsbgvd) computes all the eigenvalues and, optionally, the eigenvectors of a real generalized symmetric-definite banded eigenproblem, of the form
where
and
are symmetric and banded, and
is also positive definite. If eigenvectors are desired, it uses a divide-and-conquer algorithm.
2
Specification
Fortran Interface
Subroutine f08ucf ( |
jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, w, z, ldz, work, lwork, iwork, liwork, info) |
Integer, Intent (In) | :: | n, ka, kb, ldab, ldbb, ldz, lwork, liwork | Integer, Intent (Out) | :: | iwork(max(1,liwork)), info | Real (Kind=nag_wp), Intent (Inout) | :: | ab(ldab,*), bb(ldbb,*), z(ldz,*) | Real (Kind=nag_wp), Intent (Out) | :: | w(n), work(max(1,lwork)) | Character (1), Intent (In) | :: | jobz, uplo |
|
C Header Interface
#include <nagmk26.h>
void |
f08ucf_ (const char *jobz, const char *uplo, const Integer *n, const Integer *ka, const Integer *kb, double ab[], const Integer *ldab, double bb[], const Integer *ldbb, double w[], double z[], const Integer *ldz, double work[], const Integer *lwork, Integer iwork[], const Integer *liwork, Integer *info, const Charlen length_jobz, const Charlen length_uplo) |
|
The routine may be called by its
LAPACK
name dsbgvd.
3
Description
The generalized symmetric-definite band problem
is first reduced to a standard band symmetric problem
where
is a symmetric band matrix, using Wilkinson's modification to Crawford's algorithm (see
Crawford (1973) and
Wilkinson (1977)).
The symmetric eigenvalue problem is then solved for the eigenvalues and the eigenvectors, if required, which are then backtransformed to the eigenvectors of the original problem.
The eigenvectors are normalized so that the matrix of eigenvectors,
, satisfies
where
is the diagonal matrix whose diagonal elements are the eigenvalues.
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
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
Arguments
- 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
, the upper triangles of
and
are stored.
If , the lower triangles of and are stored.
Constraint:
or .
- 3: – IntegerInput
-
On entry: , the order of the matrices and .
Constraint:
.
- 4: – IntegerInput
-
On entry: if
, the number of superdiagonals,
, of the matrix
.
If , the number of subdiagonals, , of the matrix .
Constraint:
.
- 5: – IntegerInput
-
On entry: if
, the number of superdiagonals,
, of the matrix
.
If , the number of subdiagonals, , of the matrix .
Constraint:
.
- 6: – Real (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
symmetric 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.
- 7: – IntegerInput
-
On entry: the first dimension of the array
ab as declared in the (sub)program from which
f08ucf (dsbgvd) is called.
Constraint:
.
- 8: – Real (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
symmetric 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
f08uff (dpbstf).
- 9: – IntegerInput
-
On entry: the first dimension of the array
bb as declared in the (sub)program from which
f08ucf (dsbgvd) is called.
Constraint:
.
- 10: – Real (Kind=nag_wp) arrayOutput
-
On exit: the eigenvalues in ascending order.
- 11: – Real (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.
- 12: – IntegerInput
-
On entry: the first dimension of the array
z as declared in the (sub)program from which
f08ucf (dsbgvd) is called.
Constraints:
- if , ;
- otherwise .
- 13: – Real (Kind=nag_wp) arrayWorkspace
-
On exit: if
,
contains the minimum value of
lwork required for optimal performance.
- 14: – IntegerInput
-
On entry: the dimension of the array
work as declared in the (sub)program from which
f08ucf (dsbgvd) is called.
If
, a workspace query is assumed; the routine only calculates the minimum sizes of the
work and
iwork arrays, returns these values as the first entries of the
work and
iwork arrays, and no error message related to
lwork or
liwork is issued.
Constraints:
- if , ;
- if and , ;
- if and , .
- 15: – Integer arrayWorkspace
-
On exit: if
,
returns the minimum
liwork.
- 16: – IntegerInput
-
On entry: the dimension of the array
iwork as declared in the (sub)program from which
f08ucf (dsbgvd) is called.
If
, a workspace query is assumed; the routine only calculates the minimum sizes of the
work and
iwork arrays, returns these values as the first entries of the
work and
iwork arrays, and no error message related to
lwork or
liwork is issued.
Constraints:
- if or , ;
- if and , .
- 17: – 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.
-
The algorithm failed to converge; off-diagonal elements of an intermediate tridiagonal form did not converge to zero.
-
If
, for
,
f08uff (dpbstf) returned
:
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
f08ucf (dsbgvd) is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f08ucf (dsbgvd) 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, assuming that , is approximately proportional to otherwise.
The complex analogue of this routine is
f08uqf (zhbgvd).
10
Example
This example finds all the eigenvalues of the generalized band symmetric eigenproblem
, where
10.1
Program Text
Program Text (f08ucfe.f90)
10.2
Program Data
Program Data (f08ucfe.d)
10.3
Program Results
Program Results (f08ucfe.r)