NAG Library Routine Document

f08usf (zhbgst)

1
Purpose

f08usf (zhbgst) reduces a complex Hermitian-definite generalized eigenproblem Az=λBz to the standard form Cy=λy, where A and B are band matrices, A is a complex Hermitian matrix, and B has been factorized by f08utf (zpbstf).

2
Specification

Fortran Interface
Subroutine f08usf ( vect, uplo, n, ka, kb, ab, ldab, bb, ldbb, x, ldx, work, rwork, info)
Integer, Intent (In):: n, ka, kb, ldab, ldbb, ldx
Integer, Intent (Out):: info
Real (Kind=nag_wp), Intent (Out):: rwork(n)
Complex (Kind=nag_wp), Intent (In):: bb(ldbb,*)
Complex (Kind=nag_wp), Intent (Inout):: ab(ldab,*), x(ldx,*)
Complex (Kind=nag_wp), Intent (Out):: work(n)
Character (1), Intent (In):: vect, uplo
C Header Interface
#include <nagmk26.h>
void  f08usf_ (const char *vect, const char *uplo, const Integer *n, const Integer *ka, const Integer *kb, Complex ab[], const Integer *ldab, const Complex bb[], const Integer *ldbb, Complex x[], const Integer *ldx, Complex work[], double rwork[], Integer *info, const Charlen length_vect, const Charlen length_uplo)
The routine may be called by its LAPACK name zhbgst.

3
Description

To reduce the complex Hermitian-definite generalized eigenproblem Az=λBz to the standard form Cy=λy, where A, B and C are banded, f08usf (zhbgst) must be preceded by a call to f08utf (zpbstf) which computes the split Cholesky factorization of the positive definite matrix B: B=SHS. The split Cholesky factorization, compared with the ordinary Cholesky factorization, allows the work to be approximately halved.
This routine overwrites A with C=XHAX, where X=S-1Q and Q is a unitary matrix chosen (implicitly) to preserve the bandwidth of A. The routine also has an option to allow the accumulation of X, and then, if z is an eigenvector of C, Xz is an eigenvector of the original system.

4
References

Crawford C R (1973) Reduction of a band-symmetric generalized eigenvalue problem Comm. ACM 16 41–44
Kaufman L (1984) Banded eigenvalue solvers on vector machines ACM Trans. Math. Software 10 73–86

5
Arguments

1:     vect – Character(1)Input
On entry: indicates whether X is to be returned.
vect='N'
X is not returned.
vect='V'
X is returned.
Constraint: vect='N' or 'V'.
2:     uplo – Character(1)Input
On entry: indicates whether the upper or lower triangular part of A is stored.
uplo='U'
The upper triangular part of A is stored.
uplo='L'
The lower triangular part of A is stored.
Constraint: uplo='U' or 'L'.
3:     n – IntegerInput
On entry: n, the order of the matrices A and B.
Constraint: n0.
4:     ka – IntegerInput
On entry: if uplo='U', the number of superdiagonals, ka, of the matrix A.
If uplo='L', the number of subdiagonals, ka, of the matrix A.
Constraint: ka0.
5:     kb – IntegerInput
On entry: if uplo='U', the number of superdiagonals, kb, of the matrix B.
If uplo='L', the number of subdiagonals, kb, of the matrix B.
Constraint: kakb0.
6:     abldab* – Complex (Kind=nag_wp) arrayInput/Output
Note: the second dimension of the array ab must be at least max1,n.
On entry: the upper or lower triangle of the n by n Hermitian band matrix A.
The matrix is stored in rows 1 to ka+1, more precisely,
  • if uplo='U', the elements of the upper triangle of A within the band must be stored with element Aij in abka+1+i-jj​ for ​max1,j-kaij;
  • if uplo='L', the elements of the lower triangle of A within the band must be stored with element Aij in ab1+i-jj​ for ​jiminn,j+ka.
On exit: the upper or lower triangle of ab is overwritten by the corresponding upper or lower triangle of C as specified by uplo.
7:     ldab – IntegerInput
On entry: the first dimension of the array ab as declared in the (sub)program from which f08usf (zhbgst) is called.
Constraint: ldabka+1.
8:     bbldbb* – Complex (Kind=nag_wp) arrayInput
Note: the second dimension of the array bb must be at least max1,n.
On entry: the banded split Cholesky factor of B as specified by uplo, n and kb and returned by f08utf (zpbstf).
9:     ldbb – IntegerInput
On entry: the first dimension of the array bb as declared in the (sub)program from which f08usf (zhbgst) is called.
Constraint: ldbbkb+1.
10:   xldx* – Complex (Kind=nag_wp) arrayOutput
Note: the second dimension of the array x must be at least max1,n if vect='V' and at least 1 if vect='N'.
On exit: the n by n matrix X=S-1Q, if vect='V'.
If vect='N', x is not referenced.
11:   ldx – IntegerInput
On entry: the first dimension of the array x as declared in the (sub)program from which f08usf (zhbgst) is called.
Constraints:
  • if vect='V', ldx max1,n ;
  • if vect='N', ldx1.
12:   workn – Complex (Kind=nag_wp) arrayWorkspace
13:   rworkn – Real (Kind=nag_wp) arrayWorkspace
14:   info – IntegerOutput
On exit: info=0 unless the routine detects an error (see Section 6).

6
Error Indicators and Warnings

info<0
If info=-i, argument i had an illegal value. An explanatory message is output, and execution of the program is terminated.

7
Accuracy

Forming the reduced matrix C is a stable procedure. However it involves implicit multiplication by B-1. When f08usf (zhbgst) is used as a step in the computation of eigenvalues and eigenvectors of the original problem, there may be a significant loss of accuracy if B is ill-conditioned with respect to inversion.

8
Parallelism and Performance

f08usf (zhbgst) 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.

9
Further Comments

The total number of real floating-point operations is approximately 20n2kB, when vect='N', assuming nkA,kB; there are an additional 5n3kB/kA operations when vect='V'.
The real analogue of this routine is f08uef (dsbgst).

10
Example

This example computes all the eigenvalues of Az=λBz, where
A = -1.13+0.00i 1.94-2.10i -1.40+0.25i 0.00+0.00i 1.94+2.10i -1.91+0.00i -0.82-0.89i -0.67+0.34i -1.40-0.25i -0.82+0.89i -1.87+0.00i -1.10-0.16i 0.00+0.00i -0.67-0.34i -1.10+0.16i 0.50+0.00i  
and
B = 9.89+0.00i 1.08-1.73i 0.00+0.00i 0.00+0.00i 1.08+1.73i 1.69+0.00i -0.04+0.29i 0.00+0.00i 0.00+0.00i -0.04-0.29i 2.65+0.00i -0.33+2.24i 0.00+0.00i 0.00+0.00i -0.33-2.24i 2.17+0.00i .  
Here A is Hermitian, B is Hermitian positive definite, and A and B are treated as band matrices. B must first be factorized by f08utf (zpbstf). The program calls f08usf (zhbgst) to reduce the problem to the standard form Cy=λy, then f08hsf (zhbtrd) to reduce C to tridiagonal form, and f08jff (dsterf) to compute the eigenvalues.

10.1
Program Text

Program Text (f08usfe.f90)

10.2
Program Data

Program Data (f08usfe.d)

10.3
Program Results

Program Results (f08usfe.r)