Developer Reference for Intel® oneAPI Math Kernel Library for Fortran

ID 766686
Date 6/24/2024
Public

A newer version of this document is available. Customers should click here to go to the newest version.

Document Table of Contents

?gghrd

Reduces a pair of matrices to generalized upper Hessenberg form using orthogonal/unitary transformations.

Syntax

call sgghrd(compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, info)

call dgghrd(compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, info)

call cgghrd(compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, info)

call zgghrd(compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, info)

call gghrd(a, b [,ilo] [,ihi] [,q] [,z] [,compq] [,compz] [,info])

Include Files

  • mkl.fi, lapack.f90

Description

The routine reduces a pair of real/complex matrices (A,B) to generalized upper Hessenberg form using orthogonal/unitary transformations, where A is a general matrix and B is upper triangular. The form of the generalized eigenvalue problem is A*x = λ*B*x, and B is typically made upper triangular by computing its QR factorization and moving the orthogonal matrix Q to the left side of the equation.

This routine simultaneously reduces A to a Hessenberg matrix H:

QH*A*Z = H

and transforms B to another upper triangular matrix T:

QH*B*Z = T

in order to reduce the problem to its standard form H*y = λ*T*y, where y = ZH*x.

The orthogonal/unitary matrices Q and Z are 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

Q1*A*Z1H = (Q1*Q)*H*(Z1*Z)H

Q1*B*Z1H = (Q1*Q)*T*(Z1*Z)H

If Q1 is the orthogonal/unitary matrix from the QR factorization of B in the original equation A*x = λ*B*x, then the routine ?gghrd reduces the original problem to generalized Hessenberg form.

Input Parameters

compq

CHARACTER*1. Must be 'N', 'I', or 'V'.

If compq = 'N', matrix Q is not computed.

If compq = 'I', Q is initialized to the unit matrix, and the orthogonal/unitary matrix Q is returned;

If compq = 'V', Q must contain an orthogonal/unitary matrix Q1 on entry, and the product Q1*Q is returned.

compz

CHARACTER*1. Must be 'N', 'I', or 'V'.

If compz = 'N', matrix Z is not computed.

If compz = 'I', Z is initialized to the unit matrix, and the orthogonal/unitary matrix Z is returned;

If compz = 'V', Z must contain an orthogonal/unitary matrix Z1 on entry, and the product Z1*Z is returned.

n

INTEGER. The order of the matrices A and B (n 0).

ilo, ihi

INTEGER. ilo and ihi mark the rows and columns of A which are to be reduced. It is assumed that A is already upper triangular in rows and columns 1:ilo-1 and ihi+1:n. Values of ilo and ihi are normally set by a previous call to ggbal; otherwise they should be set to 1 and n respectively.

Constraint:

If n > 0, then 1 iloihin;

if n = 0, then ilo = 1 and ihi = 0.

a, b, q, z

REAL for sgghrd

DOUBLE PRECISION for dgghrd

COMPLEX for cgghrd

DOUBLE COMPLEX for zgghrd.

Arrays:

a(lda,*) contains the n-by-n general matrix A.

The second dimension of a must be at least max(1, n).

b(ldb,*) contains the n-by-n upper triangular matrix B.

The second dimension of b must be at least max(1, n).

q(ldq,*)

If compq = 'N', then q is not referenced.

If compq = 'V', then q must contain the orthogonal/unitary matrix Q1, typically from the QR factorization of B.

The second dimension of q must be at least max(1, n).

z(ldz,*)

If compz = 'N', then z is not referenced.

If compz = 'V', then z must contain the orthogonal/unitary matrix Z1.

The second dimension of z must be at least max(1, n).

lda

INTEGER. The leading dimension of a; at least max(1, n).

ldb

INTEGER. The leading dimension of b; at least max(1, n).

ldq

INTEGER. The leading dimension of q;

If compq = 'N', then ldq 1.

If compq = 'I'or 'V', then ldq max(1, n).

ldz

INTEGER. The leading dimension of z;

If compz = 'N', then ldz 1.

If compz = 'I'or 'V', then ldz max(1, n).

Output Parameters

a

On exit, the upper triangle and the first subdiagonal of A are overwritten with the upper Hessenberg matrix H, and the rest is set to zero.

b

On exit, overwritten by the upper triangular matrix T = QH*B*Z. The elements below the diagonal are set to zero.

q

If compq = 'I', then q contains the orthogonal/unitary matrix Q, ;

If compq = 'V', then q is overwritten by the product Q1*Q.

z

If compz = 'I', then z contains the orthogonal/unitary matrix Z;

If compz = 'V', then z is overwritten by the product Z1*Z.

info

INTEGER.

If info = 0, the execution is successful.

If info = -i, the i-th parameter had an illegal value.

LAPACK 95 Interface Notes

Routines in Fortran 95 interface have fewer arguments in the calling sequence than their FORTRAN 77 counterparts. For general conventions applied to skip redundant or restorable arguments, see LAPACK 95 Interface Conventions.

Specific details for the routine gghrd interface are the following:

a

Holds the matrix A of size (n,n).

b

Holds the matrix B of size (n,n).

q

Holds the matrix Q of size (n,n).

z

Holds the matrix Z of size (n,n).

ilo

Default value for this argument is ilo = 1.

ihi

Default value for this argument is ihi = n.

compq

If omitted, this argument is restored based on the presence of argument q as follows: compq = 'I', if q is present, compq = 'N', if q is omitted.

If present, compq must be equal to 'I' or 'V' and the argument q must also be present. Note that there will be an error condition if compq is present and q omitted.

compz

If omitted, this argument is restored based on the presence of argument z as follows: compz = 'I', if z is present, compz = 'N', if z is omitted.

If present, compz must be equal to 'I' or 'V' and the argument z must also be present. Note that there will be an error condition if compz is present and z omitted.