NAG Library Routine Document

f08taf (dspgv)

1
Purpose

f08taf (dspgv) computes all the eigenvalues and, optionally, all the 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 format, and B is also positive definite.

2
Specification

Fortran Interface
Subroutine f08taf ( itype, jobz, uplo, n, ap, bp, w, z, ldz, work, info)
Integer, Intent (In):: itype, n, ldz
Integer, Intent (Out):: info
Real (Kind=nag_wp), Intent (Inout):: ap(*), bp(*), z(ldz,*)
Real (Kind=nag_wp), Intent (Out):: w(n), work(3*n)
Character (1), Intent (In):: jobz, uplo
C Header Interface
#include <nagmk26.h>
void  f08taf_ (const Integer *itype, const char *jobz, const char *uplo, const Integer *n, double ap[], double bp[], double w[], double z[], const Integer *ldz, double work[], Integer *info, const Charlen length_jobz, const Charlen length_uplo)
The routine may be called by its LAPACK name dspgv.

3
Description

f08taf (dspgv) 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 eigenvalues and, optionally, the 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 http://www.netlib.org/lapack/lug
Golub G H and Van Loan C F (1996) Matrix Computations (3rd Edition) Johns Hopkins University Press, Baltimore

5
Arguments

1:     itype – IntegerInput
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:     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 – IntegerInput
On entry: n, the order of the matrices A and B.
Constraint: n0.
5:     ap* – Real (Kind=nag_wp) arrayInput/Output
Note: the dimension of the array ap must be at least max1,n×n+1/2.
On entry: the upper or lower triangle of the n by n symmetric matrix A, packed by columns.
More precisely,
  • if uplo='U', the upper triangle of A must be stored with element Aij in api+jj-1/2 for ij;
  • if uplo='L', the lower triangle of A must be stored with element Aij in api+2n-jj-1/2 for ij.
On exit: the contents of ap are destroyed.
6:     bp* – Real (Kind=nag_wp) arrayInput/Output
Note: the dimension of the array bp must be at least max1,n×n+1/2.
On entry: the upper or lower triangle of the n by n symmetric matrix B, packed by columns.
More precisely,
  • if uplo='U', the upper triangle of B must be stored with element Bij in bpi+jj-1/2 for ij;
  • if uplo='L', the lower triangle of B must be stored with element Bij in bpi+2n-jj-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.
7:     wn – Real (Kind=nag_wp) arrayOutput
On exit: the eigenvalues in ascending order.
8:     zldz* – Real (Kind=nag_wp) arrayOutput
Note: the second dimension of the array z must be at least max1,n if jobz='V', and at least 1 otherwise.
On exit: if jobz='V', z contains the matrix Z of eigenvectors. The eigenvectors are normalized as follows:
  • if itype=1 or 2, ZTBZ=I;
  • if itype=3, ZTB-1Z=I.
If jobz='N', z is not referenced.
9:     ldz – IntegerInput
On entry: the first dimension of the array z as declared in the (sub)program from which f08taf (dspgv) is called.
Constraints:
  • if jobz='V', ldz max1,n ;
  • otherwise ldz1.
10:   work3×n – Real (Kind=nag_wp) arrayWorkspace
11:   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.
info=1ton
The algorithm failed to converge; value off-diagonal elements of an intermediate tridiagonal form did not converge to zero.
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.
The example program below illustrates the computation of approximate error bounds.

8
Parallelism and Performance

f08taf (dspgv) is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f08taf (dspgv) 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 f08tnf (zhpgv).

10
Example

This example finds all the eigenvalues and 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 ,  
together with an estimate of the condition number of B, and approximate error bounds for the computed eigenvalues and eigenvectors.
The example program for f08tcf (dspgvd) illustrates solving a generalized symmetric eigenproblem of the form ABz=λz .

10.1
Program Text

Program Text (f08tafe.f90)

10.2
Program Data

Program Data (f08tafe.d)

10.3
Program Results

Program Results (f08tafe.r)