NAG Library Routine Document

f08gcf  (dspevd)

Warning. The specification of the arguments lwork and liwork changed at Mark 20 in the case where job='V' and n>1: the minimum dimension of the array work has been reduced whereas the minimum dimension of the array iwork has been increased.

 Contents

    1  Purpose
    7  Accuracy

1
Purpose

f08gcf (dspevd) computes all the eigenvalues and, optionally, all the eigenvectors of a real symmetric matrix held in packed storage. If the eigenvectors are requested, then it uses a divide-and-conquer algorithm to compute eigenvalues and eigenvectors. However, if only eigenvalues are required, then it uses the Pal–Walker–Kahan variant of the QL or QR algorithm.

2
Specification

Fortran Interface
Subroutine f08gcf ( job, uplo, n, ap, w, z, ldz, work, lwork, iwork, liwork, info)
Integer, Intent (In):: n, ldz, lwork, liwork
Integer, Intent (Out):: iwork(max(1,liwork)), info
Real (Kind=nag_wp), Intent (Inout):: ap(*), w(*), z(ldz,*)
Real (Kind=nag_wp), Intent (Out):: work(max(1,lwork))
Character (1), Intent (In):: job, uplo
C Header Interface
#include nagmk26.h
void  f08gcf_ ( const char *job, const char *uplo, const Integer *n, double ap[], double w[], double z[], const Integer *ldz, double work[], const Integer *lwork, Integer iwork[], const Integer *liwork, Integer *info, const Charlen length_job, const Charlen length_uplo)
The routine may be called by its LAPACK name dspevd.

3
Description

f08gcf (dspevd) computes all the eigenvalues and, optionally, all the eigenvectors of a real symmetric matrix A (held in packed storage). In other words, it can compute the spectral factorization of A as
A=ZΛZT,  
where Λ is a diagonal matrix whose diagonal elements are the eigenvalues λi, and Z is the orthogonal matrix whose columns are the eigenvectors zi. Thus
Azi=λizi,  i=1,2,,n.  

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:     job – Character(1)Input
On entry: indicates whether eigenvectors are computed.
job='N'
Only eigenvalues are computed.
job='V'
Eigenvalues and eigenvectors are computed.
Constraint: job='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 matrix A.
Constraint: n0.
4:     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: ap is overwritten by the values generated during the reduction to tridiagonal form. The elements of the diagonal and the off-diagonal of the tridiagonal matrix overwrite the corresponding elements of A.
5:     w* – Real (Kind=nag_wp) arrayOutput
Note: the dimension of the array w must be at least max1,n.
On exit: the eigenvalues of the matrix A in ascending order.
6:     zldz* – Real (Kind=nag_wp) arrayOutput
Note: the second dimension of the array z must be at least max1,n if job='V' and at least 1 if job='N'.
On exit: if job='V', z is overwritten by the orthogonal matrix Z which contains the eigenvectors of A.
If job='N', z is not referenced.
7:     ldz – IntegerInput
On entry: the first dimension of the array z as declared in the (sub)program from which f08gcf (dspevd) is called.
Constraints:
  • if job='V', ldz max1,n ;
  • if job='N', ldz1.
8:     workmax1,lwork – Real (Kind=nag_wp) arrayWorkspace
On exit: if info=0, work1 contains the required minimal size of lwork.
9:     lwork – IntegerInput
On entry: the dimension of the array work as declared in the (sub)program from which f08gcf (dspevd) is called.
If lwork=-1, a workspace query is assumed; the routine only calculates the minimum dimension of the work array, returns this value as the first entry of the work array, and no error message related to lwork is issued.
Constraints:
  • if n1, lwork1 or lwork=-1;
  • if job='N' and n>1, lwork2×n or lwork=-1;
  • if job='V' and n>1, lworkn2+6×n+1 or lwork=-1.
10:   iworkmax1,liwork – Integer arrayWorkspace
On exit: if info=0, iwork1 contains the required minimal size of liwork.
11:   liwork – IntegerInput
On entry: the dimension of the array iwork as declared in the (sub)program from which f08gcf (dspevd) is called.
If liwork=-1, a workspace query is assumed; the routine only calculates the minimum dimension of the iwork array, returns this value as the first entry of the iwork array, and no error message related to liwork is issued.
Constraints:
  • if job='N' or n1, liwork1 or liwork=-1;
  • if job='V' and n>1, liwork5×n+3 or liwork=-1.
12:   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>0
if info=i and job='N', the algorithm failed to converge; i elements of an intermediate tridiagonal form did not converge to zero; if info=i and job='V', then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and column i/n+1 through i mod n+1.

7
Accuracy

The computed eigenvalues and eigenvectors are exact for a nearby matrix A+E, where
E2 = Oε A2 ,  
and ε is the machine precision. See Section 4.7 of Anderson et al. (1999) for further details.

8
Parallelism and Performance

f08gcf (dspevd) is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f08gcf (dspevd) 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 complex analogue of this routine is f08gqf (zhpevd).

10
Example

This example computes all the eigenvalues and eigenvectors of the symmetric matrix A, where
A = 1.0 2.0 3.0 4.0 2.0 2.0 3.0 4.0 3.0 3.0 3.0 4.0 4.0 4.0 4.0 4.0 .  

10.1
Program Text

Program Text (f08gcfe.f90)

10.2
Program Data

Program Data (f08gcfe.d)

10.3
Program Results

Program Results (f08gcfe.r)

© The Numerical Algorithms Group Ltd, Oxford, UK. 2017