NAG Library Routine Document

f08jcf (dstevd)

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.

1
Purpose

f08jcf (dstevd) computes all the eigenvalues and, optionally, all the eigenvectors of a real symmetric tridiagonal matrix. 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 f08jcf ( job, n, d, e, 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):: d(*), e(*), z(ldz,*)
Real (Kind=nag_wp), Intent (Out):: work(max(1,lwork))
Character (1), Intent (In):: job
C Header Interface
#include <nagmk26.h>
void  f08jcf_ (const char *job, const Integer *n, double d[], double e[], double z[], const Integer *ldz, double work[], const Integer *lwork, Integer iwork[], const Integer *liwork, Integer *info, const Charlen length_job)
The routine may be called by its LAPACK name dstevd.

3
Description

f08jcf (dstevd) computes all the eigenvalues and, optionally, all the eigenvectors of a real symmetric tridiagonal matrix T. In other words, it can compute the spectral factorization of T as
T=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
Tzi=λizi,  i=1,2,,n.  

4
References

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:     n – IntegerInput
On entry: n, the order of the matrix T.
Constraint: n0.
3:     d* – Real (Kind=nag_wp) arrayInput/Output
Note: the dimension of the array d must be at least max1,n.
On entry: the n diagonal elements of the tridiagonal matrix T.
On exit: the eigenvalues of the matrix T in ascending order.
4:     e* – Real (Kind=nag_wp) arrayInput/Output
Note: the dimension of the array e must be at least max1,n.
On entry: the n-1 off-diagonal elements of the tridiagonal matrix T. The nth element of this array is used as workspace.
On exit: e is overwritten with intermediate results.
5:     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 T.
If job='N', z is not referenced.
6:     ldz – IntegerInput
On entry: the first dimension of the array z as declared in the (sub)program from which f08jcf (dstevd) is called.
Constraints:
  • if job='V', ldz max1,n ;
  • if job='N', ldz1.
7:     workmax1,lwork – Real (Kind=nag_wp) arrayWorkspace
On exit: if info=0, work1 contains the required minimal size of lwork.
8:     lwork – IntegerInput
On entry: the dimension of the array work as declared in the (sub)program from which f08jcf (dstevd) 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 job='N' or n1, lwork1 or lwork=-1;
  • if job='V' and n>1, lworkn2+4×n+1 or lwork=-1.
9:     iworkmax1,liwork – Integer arrayWorkspace
On exit: if info=0, iwork1 contains the required minimal size of liwork.
10:   liwork – IntegerInput
On entry: the dimension of the array iwork as declared in the (sub)program from which f08jcf (dstevd) 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.
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>0
The algorithm failed to converge; value off-diagonal elements of e did not converge to zero.

7
Accuracy

The computed eigenvalues and eigenvectors are exact for a nearby matrix T+E, where
E2 = Oε T2 ,  
and ε is the machine precision.
If λi is an exact eigenvalue and λ~i is the corresponding computed value, then
λ~i - λi c n ε T2 ,  
where cn is a modestly increasing function of n.
If zi is the corresponding exact eigenvector, and z~i is the corresponding computed eigenvector, then the angle θz~i,zi between them is bounded as follows:
θ z~i,zi c n ε T2 min ij λi - λj .  
Thus the accuracy of a computed eigenvector depends on the gap between its eigenvalue and all the other eigenvalues.

8
Parallelism and Performance

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

There is no complex analogue of this routine.

10
Example

This example computes all the eigenvalues and eigenvectors of the symmetric tridiagonal matrix T, where
T = 1.0 1.0 0.0 0.0 1.0 4.0 2.0 0.0 0.0 2.0 9.0 3.0 0.0 0.0 3.0 16.0 .  

10.1
Program Text

Program Text (f08jcfe.f90)

10.2
Program Data

Program Data (f08jcfe.d)

10.3
Program Results

Program Results (f08jcfe.r)