NAG CL Interface
f08wfc (dgghd3)

Settings help

CL Name Style:


1 Purpose

f08wfc 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>
void  f08wfc (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)
The function may be called by the names: f08wfc, nag_lapackeig_dgghd3 or nag_dgghd3.

3 Description

f08wfc is the third step in the solution of the real generalized eigenvalue problem
Ax=λBx.  
The (optional) first step balances the two matrices using f08whc. In the second step, matrix B is reduced to upper triangular form using the QR factorization function f08aec and this orthogonal transformation Q is applied to matrix A by calling f08agc. The driver, f08wcc, solves the real generalized eigenvalue problem by combining all the required steps including those just listed.
f08wfc 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=(Q1Q)H(Z1Z)T, Q1BZ1T=(Q1Q)T(Z1Z)T.  

4 References

Golub G H and Van Loan C F (2012) Matrix Computations (4th 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_OrderType Input
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.1.3 in the Introduction to the NAG Library CL Interface for a more detailed explanation of the use of this argument.
Constraint: order=Nag_RowMajor or Nag_ColMajor.
2: compq Nag_ComputeQType Input
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_ComputeZType Input
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 Integer Input
On entry: n, the order of the matrices A and B.
Constraint: n0.
5: ilo Integer Input
6: ihi Integer Input
On entry: ilo and ihi as determined by a previous call to 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] double Input/Output
Note: the dimension, dim, of the array a must be at least max(1,pda×n).
The (i,j)th 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 f08agc.
On exit: a is overwritten by the upper Hessenberg matrix H.
8: pda Integer Input
On entry: the stride separating row or column elements (depending on the value of order) in the array a.
Constraint: pdamax(1,n).
9: b[dim] double Input/Output
Note: the dimension, dim, of the array b must be at least max(1,pdb×n).
The (i,j)th 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 f08aec.
On exit: b is overwritten by the upper triangular matrix T.
10: pdb Integer Input
On entry: the stride separating row or column elements (depending on the value of order) in the array b.
Constraint: pdbmax(1,n).
11: q[dim] double Input/Output
Note: the dimension, dim, of the array q must be at least
  • max(1,pdq×n) when compq=Nag_InitQ or Nag_UpdateSchur;
  • 1 when compq=Nag_NotQ.
The (i,j)th 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 Integer Input
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 max(1,n) ;
  • if compq=Nag_NotQ, pdq1.
13: z[dim] double Input/Output
Note: the dimension, dim, of the array z must be at least
  • max(1,pdz×n) when compz=Nag_UpdateZ or Nag_InitZ;
  • 1 when compz=Nag_NotZ.
The (i,j)th 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 Integer Input
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 max(1,n) ;
  • if compz=Nag_NotZ, pdz1.
15: fail NagError * Input/Output
The NAG error argument (see Section 7 in the Introduction to the NAG Library CL Interface).

6 Error Indicators and Warnings

NE_ALLOC_FAIL
Dynamic memory allocation failed.
See Section 3.1.2 in the Introduction to the NAG Library CL Interface 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 max(1,n) ;
if compq=Nag_NotQ, pdq1.
On entry, compz=value, pdz=value and n=value.
Constraint: if compz=Nag_UpdateZ or Nag_InitZ, pdz max(1,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: pdamax(1,n).
On entry, pdb=value and n=value.
Constraint: pdbmax(1,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.
See Section 7.5 in the Introduction to the NAG Library CL Interface for further information.
NE_NO_LICENCE
Your licence key may have expired or may not have been installed correctly.
See Section 8 in the Introduction to the NAG Library CL Interface 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

f08wfc 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 f08xec which implements the QZ algorithm for computing generalized eigenvalues of a reduced pair of matrices.
The complex analogue of this function is f08wtc.

10 Example

See Section 10 in f08xec and f08ykc.