NAG FL Interface
f08upf (zhbgvx)

Settings help

FL Name Style:


FL Specification Language:


1 Purpose

f08upf computes selected the eigenvalues and, optionally, the eigenvectors of a complex generalized Hermitian-definite banded eigenproblem, of the form
Az=λBz ,  
where A and B are Hermitian and banded, and B 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

Fortran Interface
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, Intent (In) :: n, ka, kb, ldab, ldbb, ldq, il, iu, ldz
Integer, Intent (Inout) :: jfail(*)
Integer, Intent (Out) :: m, iwork(5*n), info
Real (Kind=nag_wp), Intent (In) :: vl, vu, abstol
Real (Kind=nag_wp), Intent (Out) :: w(n), rwork(7*n)
Complex (Kind=nag_wp), Intent (Inout) :: ab(ldab,*), bb(ldbb,*), q(ldq,*), z(ldz,*)
Complex (Kind=nag_wp), Intent (Out) :: work(n)
Character (1), Intent (In) :: jobz, range, uplo
C Header Interface
#include <nag.h>
void  f08upf_ (const char *jobz, const char *range, const char *uplo, const Integer *n, const Integer *ka, const Integer *kb, Complex ab[], const Integer *ldab, Complex bb[], const Integer *ldbb, Complex q[], const Integer *ldq, const double *vl, const double *vu, const Integer *il, const Integer *iu, const double *abstol, Integer *m, double w[], Complex z[], const Integer *ldz, Complex work[], double rwork[], Integer iwork[], Integer jfail[], Integer *info, const Charlen length_jobz, const Charlen length_range, const Charlen length_uplo)
The routine may be called by the names f08upf, nagf_lapackeig_zhbgvx or its LAPACK name zhbgvx.

3 Description

The generalized Hermitian-definite band problem
Az = λ Bz  
is first reduced to a standard band Hermitian problem
Cx = λx ,  
where C 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
zH A z = λ   and   zH B z = 1 .  

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 https://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 Arguments

1: jobz Character(1) Input
On entry: indicates whether eigenvectors are computed.
jobz='N'
Only eigenvalues are computed.
jobz='V'
Eigenvalues and eigenvectors are computed.
Constraint: jobz='N' or 'V'.
2: range Character(1) Input
On entry: if range='A', all eigenvalues will be found.
If range='V', all eigenvalues in the half-open interval (vl,vu] will be found.
If range='I', the ilth to iuth eigenvalues will be found.
Constraint: range='A', 'V' or 'I'.
3: uplo Character(1) Input
On entry: if uplo='U', the upper triangles of A and B are stored.
If uplo='L', the lower triangles of A and B are stored.
Constraint: uplo='U' or 'L'.
4: n Integer Input
On entry: n, the order of the matrices A and B.
Constraint: n0.
5: ka Integer Input
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.
6: kb Integer Input
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.
7: ab(ldab,*) Complex (Kind=nag_wp) array Input/Output
Note: the second dimension of the array ab must be at least max(1,n).
On entry: the upper or lower triangle of the n×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 ab(ka+1+i-j,j)​ for ​max(1,j-ka)ij;
  • if uplo='L', the elements of the lower triangle of A within the band must be stored with element Aij in ab(1+i-j,j)​ for ​jimin(n,j+ka).
On exit: the contents of ab are overwritten.
8: ldab Integer Input
On entry: the first dimension of the array ab as declared in the (sub)program from which f08upf is called.
Constraint: ldabka+1.
9: bb(ldbb,*) Complex (Kind=nag_wp) array Input/Output
Note: the second dimension of the array bb must be at least max(1,n).
On entry: the upper or lower triangle of the n×n Hermitian positive definite band matrix B.
The matrix is stored in rows 1 to kb+1, more precisely,
  • if uplo='U', the elements of the upper triangle of B within the band must be stored with element Bij in bb(kb+1+i-j,j)​ for ​max(1,j-kb)ij;
  • if uplo='L', the elements of the lower triangle of B within the band must be stored with element Bij in bb(1+i-j,j)​ for ​jimin(n,j+kb).
On exit: the factor S from the split Cholesky factorization B=SHS, as returned by f08utf.
10: ldbb Integer Input
On entry: the first dimension of the array bb as declared in the (sub)program from which f08upf is called.
Constraint: ldbbkb+1.
11: q(ldq,*) Complex (Kind=nag_wp) array Output
Note: the second dimension of the array q must be at least max(1,n) if jobz='V', and at least 1 otherwise.
On exit: if jobz='V', the n×n matrix, Q used in the reduction of the standard form, i.e., Cx=λx, from symmetric banded to tridiagonal form.
If jobz='N', q is not referenced.
12: ldq Integer Input
On entry: the first dimension of the array q as declared in the (sub)program from which f08upf is called.
Constraints:
  • if jobz='V', ldq max(1,n) ;
  • otherwise ldq1.
13: vl Real (Kind=nag_wp) Input
14: vu Real (Kind=nag_wp) Input
On entry: if range='V', the lower and upper bounds of the interval to be searched for eigenvalues.
If range='A' or 'I', vl and vu are not referenced.
Constraint: if range='V', vl<vu.
15: il Integer Input
16: iu Integer Input
On entry: if range='I', il and iu specify the indices (in ascending order) of the smallest and largest eigenvalues to be returned, respectively.
If range='A' or 'V', il and iu are not referenced.
Constraints:
  • if range='I' and n=0, il=1 and iu=0;
  • if range='I' and n>0, 1 il iu n .
17: abstol 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 [a,b] of width less than or equal to
abstol+ε max(|a|,|b|) ,  
where ε is the machine precision. If abstol is less than or equal to zero, then ε T1 will be used in its place, where T is the tridiagonal matrix obtained by reducing C to tridiagonal form. Eigenvalues will be computed most accurately when abstol is set to twice the underflow threshold 2 × x02amf ( ) , not zero. If this routine returns with info=1,,n, indicating that some eigenvectors did not converge, try setting abstol to 2 × x02amf ( ) . See Demmel and Kahan (1990).
18: m Integer Output
On exit: the total number of eigenvalues found. 0mn.
If range='A', m=n.
If range='I', m=iu-il+1.
19: w(n) Real (Kind=nag_wp) array Output
On exit: the eigenvalues in ascending order.
20: z(ldz,*) Complex (Kind=nag_wp) array Output
Note: the second dimension of the array z must be at least max(1,m) if jobz='V', and at least 1 otherwise.
On exit: if jobz='V', then
  • if info=0, the first m columns of Z contain the eigenvectors corresponding to the selected eigenvalues, with the ith column of Z holding the eigenvector associated with w(i). The eigenvectors are normalized so that ZHBZ=I;
  • if an eigenvector fails to converge (info=1,,n), then that column of Z contains the latest approximation to the eigenvector, and the index of the eigenvector is returned in jfail.
If jobz='N', z is not referenced.
Note:  you must ensure that at least max(1,m) columns are supplied in the array z; if range='V', the exact value of m is not known in advance and an upper bound of at least n must be used.
21: ldz Integer Input
On entry: the first dimension of the array z as declared in the (sub)program from which f08upf is called.
Constraints:
  • if jobz='V', ldz max(1,n) ;
  • otherwise ldz1.
22: work(n) Complex (Kind=nag_wp) array Workspace
23: rwork(7×n) Real (Kind=nag_wp) array Workspace
24: iwork(5×n) Integer array Workspace
25: jfail(*) Integer array Output
Note: the dimension of the array jfail must be at least max(1,n).
On exit: if jobz='V', then
  • if info=0, the first m elements of jfail are zero;
  • if info=1,,n, the first info elements of jfail contains the indices of the eigenvectors that failed to converge.
If jobz='N', jfail is not referenced.
26: info Integer Output
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.
info=1,,n
The algorithm failed to converge; value eigenvectors did not converge. Their indices are stored in array jfail.
info>n
If info=n+value, for 1valuen, f08utf returned info=value: B is not positive definite. The factorization of B could not be completed and no eigenvalues or eigenvectors were computed.

7 Accuracy

If B 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 B differ widely in magnitude the eigenvalues and eigenvectors may be less sensitive than the condition of B would suggest. See Section 4.10 of Anderson et al. (1999) for details of the error bounds.

8 Parallelism and Performance

f08upf is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f08upf 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 floating-point operations is proportional to n3 if jobz='V' and range='A', and assuming that nka , is approximately proportional to n2 ka if jobz='N'. Otherwise the number of floating-point operations depends upon the number of eigenvectors computed.
The real analogue of this routine is f08ubf.

10 Example

This example finds the eigenvalues in the half-open interval (0.0,2.0] , and corresponding eigenvectors, of the generalized band Hermitian eigenproblem Az = λ Bz , where
A = ( -1.13i+0.00 1.94-2.10i -1.40+0.25i 0.00i+0.00 1.94+2.10i -1.91i+0.00 -0.82-0.89i -0.67+0.34i -1.40-0.25i -0.82+0.89i -1.87i+0.00 -1.10-0.16i 0.00i+0.00 -0.67-0.34i -1.10+0.16i 0.50i+0.00 )  
and
B = ( 9.89i+0.00 1.08-1.73i 0.00i+0.00 0.00i+0.00 1.08+1.73i 1.69i+0.00 -0.04+0.29i 0.00i+0.00 0.00i+0.00 -0.04-0.29i 2.65i+0.00 -0.33+2.24i 0.00i+0.00 0.00i+0.00 -0.33-2.24i 2.17i+0.00 ) .  

10.1 Program Text

Program Text (f08upfe.f90)

10.2 Program Data

Program Data (f08upfe.d)

10.3 Program Results

Program Results (f08upfe.r)