Developer Reference for Intel® oneAPI Math Kernel Library for Fortran

ID 766686
Date 11/07/2023
Public

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

Document Table of Contents

?gghd3

Reduces a pair of matrices to generalized upper Hessenberg form.

Syntax

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

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

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

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

Include Files

  • mkl.fi

Description

?gghd3 reduces a pair of real or 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/unitary matrix Q to the left side of the equation.

This subroutine simultaneously reduces A to a Hessenberg matrix H:

QT*A*Z = H for real flavors

or

QT*A*Z = H for complex flavors

and transforms B to another upper triangular matrix T:

QT*B*Z = T for real flavors

or

QT*B*Z = T for complex flavors

in order to reduce the problem to its standard form

H*y = λ*T*y

where y = ZT*x for real flavors

or

y = ZT*x for complex flavors.

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

for real flavors:

Q1 * A * Z1T = (Q1*Q) * H * (Z1*Z)T

Q1 * B * Z1T = (Q1*Q) * T * (Z1*Z)T

for complex flavors:

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

Q1 * B * Z1T = (Q1*Q) * T * (Z1*Z)T

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

This is a blocked variant of ?gghrd, using matrix-matrix multiplications for parts of the computation to enhance performance.

Input Parameters

compq

CHARACTER*1. = 'N': do not compute q;

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

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

compz

CHARACTER*1. = 'N': do not compute z;

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

= '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. ilo and ihi are normally set by a previous call to ?ggbal; otherwise they should be set to 1 and n, respectively.

1 iloihin, if n > 0; ilo=1 and ihi=0, if n=0.

a

REAL for sgghd3

DOUBLE PRECISION for dgghd3

COMPLEX for cgghd3

DOUBLE COMPLEX for zgghd3

Array, size (lda, n).

On entry, the n-by-n general matrix to be reduced.

lda

INTEGER. The leading dimension of the array a.

lda max(1,n).

b

REAL for sgghd3

DOUBLE PRECISION for dgghd3

COMPLEX for cgghd3

DOUBLE COMPLEX for zgghd3

Array, (ldb, n).

On entry, then-by-n upper triangular matrix B.

ldb

INTEGER. The leading dimension of the array b.

ldb max(1,n).

q

REAL for sgghd3

DOUBLE PRECISION for dgghd3

COMPLEX for cgghd3

DOUBLE COMPLEX for zgghd3

Array, size (ldq, n).

On entry, if compq = 'V', the orthogonal/unitary matrix Q1, typically from the QR factorization of b.

ldq

INTEGER. The leading dimension of the array q.

ldqn if compq='V' or 'I'; ldq 1 otherwise.

z

REAL for sgghd3

DOUBLE PRECISION for dgghd3

COMPLEX for cgghd3

DOUBLE COMPLEX for zgghd3

Array, size (ldz, n).

On entry, if compz = 'V', the orthogonal/unitary matrix Z1.

Not referenced if compz='N'.

ldz

INTEGER. The leading dimension of the array z. ldzn if compz='V' or 'I'; ldz 1 otherwise.

work

REAL for sgghd3

DOUBLE PRECISION for dgghd3

COMPLEX for cgghd3

DOUBLE COMPLEX for zgghd3

Array, size (lwork)

lwork

INTEGER. The length of the array work.

lwork 1.

For optimum performance lwork 6*n*NB, where NB is the optimal blocksize.

If lwork = -1, then a workspace query is assumed; the routine only calculates the optimal size of the work array, returns this value as the first entry of the work array, and no error message related to lwork is issued by xerbla.

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, the upper triangular matrix T = QTBZ for real flavors or T = QHBZ for complex flavors. The elements below the diagonal are set to zero.

q

On exit, if compq='I', the orthogonal/unitary matrix Q, and if compq = 'V', the product Q1*Q.

Not referenced if compq='N'.

z

On exit, if compz='I', the orthogonal/unitary matrix Z, and if compz = 'V', the product Z1*Z.

Not referenced if compz='N'.

work

On exit, if info = 0, work(1) returns the optimal lwork.

info

INTEGER. = 0: successful exit.

< 0: if info = -i, the i-th argument had an illegal value.

Application Notes

This routine reduces A to Hessenberg form and maintains B in using a blocked variant of Moler and Stewart's original algorithm, as described by Kagstrom, Kressner, Quintana-Orti, and Quintana-Orti (BIT 2008).