NAG FL Interface
f08tbf (dspgvx)

Settings help

FL Name Style:


FL Specification Language:


1 Purpose

f08tbf computes selected eigenvalues and, optionally, eigenvectors of a real generalized symmetric-definite eigenproblem, of the form
Az=λBz ,   ABz=λz   or   BAz=λz ,  
where A and B are symmetric, stored in packed storage, and B is also positive definite. Eigenvalues and eigenvectors can be selected by specifying either a range of values or a range of indices for the desired eigenvalues.

2 Specification

Fortran Interface
Subroutine f08tbf ( itype, jobz, range, uplo, n, ap, bp, vl, vu, il, iu, abstol, m, w, z, ldz, work, iwork, jfail, info)
Integer, Intent (In) :: itype, n, 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 (Inout) :: ap(*), bp(*), z(ldz,*)
Real (Kind=nag_wp), Intent (Out) :: w(n), work(8*n)
Character (1), Intent (In) :: jobz, range, uplo
C Header Interface
#include <nag.h>
void  f08tbf_ (const Integer *itype, const char *jobz, const char *range, const char *uplo, const Integer *n, double ap[], double bp[], const double *vl, const double *vu, const Integer *il, const Integer *iu, const double *abstol, Integer *m, double w[], double z[], const Integer *ldz, double work[], 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 f08tbf, nagf_lapackeig_dspgvx or its LAPACK name dspgvx.

3 Description

f08tbf first performs a Cholesky factorization of the matrix B as B=UTU , when uplo='U' or B=LLT , when uplo='L'. The generalized problem is then reduced to a standard symmetric eigenvalue problem
Cx=λx ,  
which is solved for the desired eigenvalues and eigenvectors; the eigenvectors are then backtransformed to give the eigenvectors of the original problem.
For the problem Az=λBz , the eigenvectors are normalized so that the matrix of eigenvectors, Z, satisfies
ZT A Z = Λ   and   ZT B Z = I ,  
where Λ is the diagonal matrix whose diagonal elements are the eigenvalues. For the problem A B z = λ z we correspondingly have
Z-1 A Z-T = Λ   and   ZT B Z = I ,  
and for B A z = λ z we have
ZT A Z = Λ   and   ZT B-1 Z = I .  

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
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

5 Arguments

1: itype Integer Input
On entry: specifies the problem type to be solved.
itype=1
Az=λBz.
itype=2
ABz=λz.
itype=3
BAz=λz.
Constraint: itype=1, 2 or 3.
2: 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'.
3: 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'.
4: 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'.
5: n Integer Input
On entry: n, the order of the matrices A and B.
Constraint: n0.
6: ap(*) Real (Kind=nag_wp) array Input/Output
Note: the dimension of the array ap must be at least max(1,n×(n+1)/2).
On entry: the upper or lower triangle of the n×n symmetric matrix A, packed by columns.
More precisely,
  • if uplo='U', the upper triangle of A must be stored with element Aij in ap(i+j(j-1)/2) for ij;
  • if uplo='L', the lower triangle of A must be stored with element Aij in ap(i+(2n-j)(j-1)/2) for ij.
On exit: the contents of ap are destroyed.
7: bp(*) Real (Kind=nag_wp) array Input/Output
Note: the dimension of the array bp must be at least max(1,n×(n+1)/2).
On entry: the upper or lower triangle of the n×n symmetric matrix B, packed by columns.
More precisely,
  • if uplo='U', the upper triangle of B must be stored with element Bij in bp(i+j(j-1)/2) for ij;
  • if uplo='L', the lower triangle of B must be stored with element Bij in bp(i+(2n-j)(j-1)/2) for ij.
On exit: the triangular factor U or L from the Cholesky factorization B=UTU or B=LLT, in the same storage format as B.
8: vl Real (Kind=nag_wp) Input
9: 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.
10: il Integer Input
11: 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 .
12: 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).
13: m Integer Output
On exit: the total number of eigenvalues found. 0mn.
If range='A', m=n.
If range='I', m=iu-il+1.
14: w(n) Real (Kind=nag_wp) array Output
On exit: the first m elements contain the selected eigenvalues in ascending order.
15: z(ldz,*) Real (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 orthonormal eigenvectors of the matrix A corresponding to the selected eigenvalues, with the ith column of Z holding the eigenvector associated with w(i). The eigenvectors are normalized as follows:
    • if itype=1 or 2, ZTBZ=I;
    • if itype=3, ZTB-1Z=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.
16: ldz Integer Input
On entry: the first dimension of the array z as declared in the (sub)program from which f08tbf is called.
Constraints:
  • if jobz='V', ldz max(1,n) ;
  • otherwise ldz1.
17: work(8×n) Real (Kind=nag_wp) array Workspace
18: iwork(5×n) Integer array Workspace
19: 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.
20: 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, then the leading minor of order value of 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

Background information to multithreading can be found in the Multithreading documentation.
f08tbf is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f08tbf 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 .
The complex analogue of this routine is f08tpf.

10 Example

This example finds the eigenvalues in the half-open interval (-1.0,1.0] , and corresponding eigenvectors, of the generalized symmetric eigenproblem Az = λ Bz , where
A = ( 0.24 0.39 0.42 -0.16 0.39 -0.11 0.79 0.63 0.42 0.79 -0.25 0.48 -0.16 0.63 0.48 -0.03 )   and   B = ( 4.16 -3.12 0.56 -0.10 -3.12 5.03 -0.83 1.09 0.56 -0.83 0.76 0.34 -0.10 1.09 0.34 1.18 ) .  
The example program for f08tcf illustrates solving a generalized symmetric eigenproblem of the form ABz=λz .

10.1 Program Text

Program Text (f08tbfe.f90)

10.2 Program Data

Program Data (f08tbfe.d)

10.3 Program Results

Program Results (f08tbfe.r)