nag_dgghrd (f08wec) (PDF version)
f08 Chapter Contents
f08 Chapter Introduction
NAG Library Manual

NAG Library Function Document

nag_dgghrd (f08wec)

 Contents

    1  Purpose
    7  Accuracy
    10  Example

1  Purpose

nag_dgghrd (f08wec) reduces a pair of real matrices A,B, where B is upper triangular, to the generalized upper Hessenberg form using orthogonal transformations.

2  Specification

#include <nag.h>
#include <nagf08.h>
void  nag_dgghrd (Nag_OrderType order, Nag_ComputeQType compq, Nag_ComputeZType compz, Integer n, Integer ilo, Integer ihi, double a[], Integer pda, double b[], Integer pdb, double q[], Integer pdq, double z[], Integer pdz, NagError *fail)

3  Description

nag_dgghrd (f08wec) is the third step in the solution of the real generalized eigenvalue problem
Ax=λBx.  
The (optional) first step balances the two matrices using nag_dggbal (f08whc). In the second step, matrix B is reduced to upper triangular form using the QR factorization function nag_dgeqrf (f08aec) and this orthogonal transformation Q is applied to matrix A by calling nag_dormqr (f08agc).
nag_dgghrd (f08wec) reduces a pair of real matrices A,B, where B is upper triangular, to the generalized upper Hessenberg form using orthogonal transformations. This two-sided transformation is of the form
QTAZ=H QTBZ=T  
where H is an upper Hessenberg matrix, T is an upper triangular matrix and Q and Z are orthogonal matrices determined as products of Givens rotations. They may either be formed explicitly, or they may be postmultiplied into input matrices Q1 and Z1, so that
Q1AZ1T=Q1QHZ1ZT, Q1BZ1T=Q1QTZ1ZT.  

4  References

Golub G H and Van Loan C F (1996) Matrix Computations (3rd Edition) Johns Hopkins University Press, Baltimore
Moler C B and Stewart G W (1973) An algorithm for generalized matrix eigenproblems SIAM J. Numer. Anal. 10 241–256

5  Arguments

1:     order Nag_OrderTypeInput
On entry: the order argument specifies the two-dimensional storage scheme being used, i.e., row-major ordering or column-major ordering. C language defined storage is specified by order=Nag_RowMajor. See Section 3.2.1.3 in the Essential Introduction for a more detailed explanation of the use of this argument.
Constraint: order=Nag_RowMajor or Nag_ColMajor.
2:     compq Nag_ComputeQTypeInput
On entry: specifies the form of the computed orthogonal matrix Q.
compq=Nag_NotQ
Do not compute Q.
compq=Nag_InitQ
The orthogonal matrix Q is returned.
compq=Nag_UpdateSchur
q must contain an orthogonal matrix Q1, and the product Q1Q is returned.
Constraint: compq=Nag_NotQ, Nag_InitQ or Nag_UpdateSchur.
3:     compz Nag_ComputeZTypeInput
On entry: specifies the form of the computed orthogonal matrix Z.
compz=Nag_NotZ
Do not compute Z.
compz=Nag_InitZ
The orthogonal matrix Z is returned.
compz=Nag_UpdateZ
z must contain an orthogonal matrix Z1, and the product Z1Z is returned.
Constraint: compz=Nag_NotZ, Nag_UpdateZ or Nag_InitZ.
4:     n IntegerInput
On entry: n, the order of the matrices A and B.
Constraint: n0.
5:     ilo IntegerInput
6:     ihi IntegerInput
On entry: ilo and ihi as determined by a previous call to nag_dggbal (f08whc). Otherwise, they should be set to 1 and n, respectively.
Constraints:
  • if n>0, 1 ilo ihi n ;
  • if n=0, ilo=1 and ihi=0.
7:     a[dim] doubleInput/Output
Note: the dimension, dim, of the array a must be at least max1,pda×n.
The i,jth element of the matrix A is stored in
  • a[j-1×pda+i-1] when order=Nag_ColMajor;
  • a[i-1×pda+j-1] when order=Nag_RowMajor.
On entry: the matrix A of the matrix pair A,B. Usually, this is the matrix A returned by nag_dormqr (f08agc).
On exit: a is overwritten by the upper Hessenberg matrix H.
8:     pda IntegerInput
On entry: the stride separating row or column elements (depending on the value of order) in the array a.
Constraint: pdamax1,n.
9:     b[dim] doubleInput/Output
Note: the dimension, dim, of the array b must be at least max1,pdb×n.
The i,jth element of the matrix B is stored in
  • b[j-1×pdb+i-1] when order=Nag_ColMajor;
  • b[i-1×pdb+j-1] when order=Nag_RowMajor.
On entry: the upper triangular matrix B of the matrix pair A,B. Usually, this is the matrix B returned by the QR factorization function nag_dgeqrf (f08aec).
On exit: b is overwritten by the upper triangular matrix T.
10:   pdb IntegerInput
On entry: the stride separating row or column elements (depending on the value of order) in the array b.
Constraint: pdbmax1,n.
11:   q[dim] doubleInput/Output
Note: the dimension, dim, of the array q must be at least
  • max1,pdq×n when compq=Nag_InitQ or Nag_UpdateSchur;
  • 1 when compq=Nag_NotQ.
The i,jth element of the matrix Q is stored in
  • q[j-1×pdq+i-1] when order=Nag_ColMajor;
  • q[i-1×pdq+j-1] when order=Nag_RowMajor.
On entry: if compq=Nag_UpdateSchur, q must contain an orthogonal matrix Q1.
If compq=Nag_NotQ, q is not referenced.
On exit: if compq=Nag_InitQ, q contains the orthogonal matrix Q.
If compq=Nag_UpdateSchur, q is overwritten by Q1Q.
12:   pdq IntegerInput
On entry: the stride separating row or column elements (depending on the value of order) in the array q.
Constraints:
  • if compq=Nag_InitQ or Nag_UpdateSchur, pdq max1,n ;
  • if compq=Nag_NotQ, pdq1.
13:   z[dim] doubleInput/Output
Note: the dimension, dim, of the array z must be at least max1,pdz×n when compz=Nag_UpdateZ or Nag_InitZ.
The i,jth element of the matrix Z is stored in
  • z[j-1×pdz+i-1] when order=Nag_ColMajor;
  • z[i-1×pdz+j-1] when order=Nag_RowMajor.
On entry: if compz=Nag_UpdateZ, z must contain an orthogonal matrix Z1.
If compz=Nag_NotZ, z is not referenced.
On exit: if compz=Nag_InitZ, z contains the orthogonal matrix Z.
If compz=Nag_UpdateZ, z is overwritten by Z1Z.
14:   pdz IntegerInput
On entry: the stride separating row or column elements (depending on the value of order) in the array z.
Constraints:
  • if compz=Nag_UpdateZ or Nag_InitZ, pdz max1,n ;
  • if compz=Nag_NotZ, pdz1.
15:   fail NagError *Input/Output
The NAG error argument (see Section 3.6 in the Essential Introduction).

6  Error Indicators and Warnings

NE_ALLOC_FAIL
Dynamic memory allocation failed.
See Section 3.2.1.2 in the Essential Introduction for further information.
NE_BAD_PARAM
On entry, argument value had an illegal value.
NE_ENUM_INT_2
On entry, compq=value, pdq=value and n=value.
Constraint: if compq=Nag_InitQ or Nag_UpdateSchur, pdq max1,n ;
if compq=Nag_NotQ, pdq1.
On entry, compz=value, pdz=value and n=value.
Constraint: if compz=Nag_UpdateZ or Nag_InitZ, pdz max1,n ;
if compz=Nag_NotZ, pdz1.
NE_INT
On entry, n=value.
Constraint: n0.
On entry, pda=value.
Constraint: pda>0.
On entry, pdb=value.
Constraint: pdb>0.
On entry, pdq=value.
Constraint: pdq>0.
On entry, pdz=value.
Constraint: pdz>0.
NE_INT_2
On entry, pda=value and n=value.
Constraint: pdamax1,n.
On entry, pdb=value and n=value.
Constraint: pdbmax1,n.
NE_INT_3
On entry, n=value, ilo=value and ihi=value.
Constraint: if n>0, 1 ilo ihi n ;
if n=0, ilo=1 and ihi=0.
NE_INTERNAL_ERROR
An internal error has occurred in this function. Check the function call and any array sizes. If the call is correct then please contact NAG for assistance.
An unexpected error has been triggered by this function. Please contact NAG.
See Section 3.6.6 in the Essential Introduction for further information.
NE_NO_LICENCE
Your licence key may have expired or may not have been installed correctly.
See Section 3.6.5 in the Essential Introduction for further information.

7  Accuracy

The reduction to the generalized Hessenberg form is implemented using orthogonal transformations which are backward stable.

8  Parallelism and Performance

nag_dgghrd (f08wec) is not threaded by NAG in any implementation.
nag_dgghrd (f08wec) 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 function. Please also consult the Users' Note for your implementation for any additional implementation-specific information.

9  Further Comments

This function is usually followed by nag_dhgeqz (f08xec) which implements the QZ algorithm for computing generalized eigenvalues of a reduced pair of matrices.
The complex analogue of this function is nag_zgghrd (f08wsc).

10  Example

See Section 10 in nag_dhgeqz (f08xec) and nag_dtgevc (f08ykc).

nag_dgghrd (f08wec) (PDF version)
f08 Chapter Contents
f08 Chapter Introduction
NAG Library Manual

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