Scroll to navigation

gghd3(3) LAPACK gghd3(3)

NAME

gghd3 - gghd3: reduction to Hessenberg, level 3

SYNOPSIS

Functions


subroutine cgghd3 (compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, work, lwork, info)
CGGHD3 subroutine dgghd3 (compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, work, lwork, info)
DGGHD3 subroutine sgghd3 (compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, work, lwork, info)
SGGHD3 subroutine zgghd3 (compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, work, lwork, info)
ZGGHD3

Detailed Description

Function Documentation

subroutine cgghd3 (character compq, character compz, integer n, integer ilo, integer ihi, complex, dimension( lda, * ) a, integer lda, complex, dimension( ldb, * ) b, integer ldb, complex, dimension( ldq, * ) q, integer ldq, complex, dimension( ldz, * ) z, integer ldz, complex, dimension( * ) work, integer lwork, integer info)

CGGHD3

Purpose:


CGGHD3 reduces a pair of complex matrices (A,B) to generalized upper
Hessenberg form using unitary transformations, where A is a
general matrix and B is upper triangular. The form of the
generalized eigenvalue problem is
A*x = lambda*B*x,
and B is typically made upper triangular by computing its QR
factorization and moving the unitary matrix Q to the left side
of the equation.
This subroutine simultaneously reduces A to a Hessenberg matrix H:
Q**H*A*Z = H
and transforms B to another upper triangular matrix T:
Q**H*B*Z = T
in order to reduce the problem to its standard form
H*y = lambda*T*y
where y = Z**H*x.
The 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 * Z1**H = (Q1*Q) * H * (Z1*Z)**H
Q1 * B * Z1**H = (Q1*Q) * T * (Z1*Z)**H
If Q1 is the unitary matrix from the QR factorization of B in the
original equation A*x = lambda*B*x, then CGGHD3 reduces the original
problem to generalized Hessenberg form.
This is a blocked variant of CGGHRD, using matrix-matrix
multiplications for parts of the computation to enhance performance.

Parameters

COMPQ


COMPQ is CHARACTER*1
= 'N': do not compute Q;
= 'I': Q is initialized to the unit matrix, and the
unitary matrix Q is returned;
= 'V': Q must contain a unitary matrix Q1 on entry,
and the product Q1*Q is returned.

COMPZ


COMPZ is CHARACTER*1
= 'N': do not compute Z;
= 'I': Z is initialized to the unit matrix, and the
unitary matrix Z is returned;
= 'V': Z must contain a unitary matrix Z1 on entry,
and the product Z1*Z is returned.

N


N is INTEGER
The order of the matrices A and B. N >= 0.

ILO


ILO is INTEGER

IHI


IHI is 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 CGGBAL; otherwise they
should be set to 1 and N respectively.
1 <= ILO <= IHI <= N, if N > 0; ILO=1 and IHI=0, if N=0.

A


A is COMPLEX array, dimension (LDA, N)
On entry, the N-by-N general matrix to be reduced.
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.

LDA


LDA is INTEGER
The leading dimension of the array A. LDA >= max(1,N).

B


B is COMPLEX array, dimension (LDB, N)
On entry, the N-by-N upper triangular matrix B.
On exit, the upper triangular matrix T = Q**H B Z. The
elements below the diagonal are set to zero.

LDB


LDB is INTEGER
The leading dimension of the array B. LDB >= max(1,N).

Q


Q is COMPLEX array, dimension (LDQ, N)
On entry, if COMPQ = 'V', the unitary matrix Q1, typically
from the QR factorization of B.
On exit, if COMPQ='I', the unitary matrix Q, and if
COMPQ = 'V', the product Q1*Q.
Not referenced if COMPQ='N'.

LDQ


LDQ is INTEGER
The leading dimension of the array Q.
LDQ >= N if COMPQ='V' or 'I'; LDQ >= 1 otherwise.

Z


Z is COMPLEX array, dimension (LDZ, N)
On entry, if COMPZ = 'V', the unitary matrix Z1.
On exit, if COMPZ='I', the unitary matrix Z, and if
COMPZ = 'V', the product Z1*Z.
Not referenced if COMPZ='N'.

LDZ


LDZ is INTEGER
The leading dimension of the array Z.
LDZ >= N if COMPZ='V' or 'I'; LDZ >= 1 otherwise.

WORK


WORK is COMPLEX array, dimension (LWORK)
On exit, if INFO = 0, WORK(1) returns the optimal LWORK.

LWORK


LWORK is 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.

INFO


INFO is INTEGER
= 0: successful exit.
< 0: if INFO = -i, the i-th argument had an illegal value.

Author

Univ. of Tennessee

Univ. of California Berkeley

Univ. of Colorado Denver

NAG Ltd.

Further Details:


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

subroutine dgghd3 (character compq, character compz, integer n, integer ilo, integer ihi, double precision, dimension( lda, * ) a, integer lda, double precision, dimension( ldb, * ) b, integer ldb, double precision, dimension( ldq, * ) q, integer ldq, double precision, dimension( ldz, * ) z, integer ldz, double precision, dimension( * ) work, integer lwork, integer info)

DGGHD3

Purpose:


DGGHD3 reduces a pair of real matrices (A,B) to generalized upper
Hessenberg form using orthogonal transformations, where A is a
general matrix and B is upper triangular. The form of the
generalized eigenvalue problem is
A*x = lambda*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 subroutine simultaneously reduces A to a Hessenberg matrix H:
Q**T*A*Z = H
and transforms B to another upper triangular matrix T:
Q**T*B*Z = T
in order to reduce the problem to its standard form
H*y = lambda*T*y
where y = Z**T*x.
The orthogonal 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 * Z1**T = (Q1*Q) * H * (Z1*Z)**T
Q1 * B * Z1**T = (Q1*Q) * T * (Z1*Z)**T
If Q1 is the orthogonal matrix from the QR factorization of B in the
original equation A*x = lambda*B*x, then DGGHD3 reduces the original
problem to generalized Hessenberg form.
This is a blocked variant of DGGHRD, using matrix-matrix
multiplications for parts of the computation to enhance performance.

Parameters

COMPQ


COMPQ is CHARACTER*1
= 'N': do not compute Q;
= 'I': Q is initialized to the unit matrix, and the
orthogonal matrix Q is returned;
= 'V': Q must contain an orthogonal matrix Q1 on entry,
and the product Q1*Q is returned.

COMPZ


COMPZ is CHARACTER*1
= 'N': do not compute Z;
= 'I': Z is initialized to the unit matrix, and the
orthogonal matrix Z is returned;
= 'V': Z must contain an orthogonal matrix Z1 on entry,
and the product Z1*Z is returned.

N


N is INTEGER
The order of the matrices A and B. N >= 0.

ILO


ILO is INTEGER

IHI


IHI is 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 DGGBAL; otherwise they
should be set to 1 and N respectively.
1 <= ILO <= IHI <= N, if N > 0; ILO=1 and IHI=0, if N=0.

A


A is DOUBLE PRECISION array, dimension (LDA, N)
On entry, the N-by-N general matrix to be reduced.
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.

LDA


LDA is INTEGER
The leading dimension of the array A. LDA >= max(1,N).

B


B is DOUBLE PRECISION array, dimension (LDB, N)
On entry, the N-by-N upper triangular matrix B.
On exit, the upper triangular matrix T = Q**T B Z. The
elements below the diagonal are set to zero.

LDB


LDB is INTEGER
The leading dimension of the array B. LDB >= max(1,N).

Q


Q is DOUBLE PRECISION array, dimension (LDQ, N)
On entry, if COMPQ = 'V', the orthogonal matrix Q1,
typically from the QR factorization of B.
On exit, if COMPQ='I', the orthogonal matrix Q, and if
COMPQ = 'V', the product Q1*Q.
Not referenced if COMPQ='N'.

LDQ


LDQ is INTEGER
The leading dimension of the array Q.
LDQ >= N if COMPQ='V' or 'I'; LDQ >= 1 otherwise.

Z


Z is DOUBLE PRECISION array, dimension (LDZ, N)
On entry, if COMPZ = 'V', the orthogonal matrix Z1.
On exit, if COMPZ='I', the orthogonal matrix Z, and if
COMPZ = 'V', the product Z1*Z.
Not referenced if COMPZ='N'.

LDZ


LDZ is INTEGER
The leading dimension of the array Z.
LDZ >= N if COMPZ='V' or 'I'; LDZ >= 1 otherwise.

WORK


WORK is DOUBLE PRECISION array, dimension (LWORK)
On exit, if INFO = 0, WORK(1) returns the optimal LWORK.

LWORK


LWORK is 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.

INFO


INFO is INTEGER
= 0: successful exit.
< 0: if INFO = -i, the i-th argument had an illegal value.

Author

Univ. of Tennessee

Univ. of California Berkeley

Univ. of Colorado Denver

NAG Ltd.

Further Details:


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

subroutine sgghd3 (character compq, character compz, integer n, integer ilo, integer ihi, real, dimension( lda, * ) a, integer lda, real, dimension( ldb, * ) b, integer ldb, real, dimension( ldq, * ) q, integer ldq, real, dimension( ldz, * ) z, integer ldz, real, dimension( * ) work, integer lwork, integer info)

SGGHD3

Purpose:


SGGHD3 reduces a pair of real matrices (A,B) to generalized upper
Hessenberg form using orthogonal transformations, where A is a
general matrix and B is upper triangular. The form of the
generalized eigenvalue problem is
A*x = lambda*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 subroutine simultaneously reduces A to a Hessenberg matrix H:
Q**T*A*Z = H
and transforms B to another upper triangular matrix T:
Q**T*B*Z = T
in order to reduce the problem to its standard form
H*y = lambda*T*y
where y = Z**T*x.
The orthogonal 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 * Z1**T = (Q1*Q) * H * (Z1*Z)**T
Q1 * B * Z1**T = (Q1*Q) * T * (Z1*Z)**T
If Q1 is the orthogonal matrix from the QR factorization of B in the
original equation A*x = lambda*B*x, then SGGHD3 reduces the original
problem to generalized Hessenberg form.
This is a blocked variant of SGGHRD, using matrix-matrix
multiplications for parts of the computation to enhance performance.

Parameters

COMPQ


COMPQ is CHARACTER*1
= 'N': do not compute Q;
= 'I': Q is initialized to the unit matrix, and the
orthogonal matrix Q is returned;
= 'V': Q must contain an orthogonal matrix Q1 on entry,
and the product Q1*Q is returned.

COMPZ


COMPZ is CHARACTER*1
= 'N': do not compute Z;
= 'I': Z is initialized to the unit matrix, and the
orthogonal matrix Z is returned;
= 'V': Z must contain an orthogonal matrix Z1 on entry,
and the product Z1*Z is returned.

N


N is INTEGER
The order of the matrices A and B. N >= 0.

ILO


ILO is INTEGER

IHI


IHI is 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 SGGBAL; otherwise they
should be set to 1 and N respectively.
1 <= ILO <= IHI <= N, if N > 0; ILO=1 and IHI=0, if N=0.

A


A is REAL array, dimension (LDA, N)
On entry, the N-by-N general matrix to be reduced.
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.

LDA


LDA is INTEGER
The leading dimension of the array A. LDA >= max(1,N).

B


B is REAL array, dimension (LDB, N)
On entry, the N-by-N upper triangular matrix B.
On exit, the upper triangular matrix T = Q**T B Z. The
elements below the diagonal are set to zero.

LDB


LDB is INTEGER
The leading dimension of the array B. LDB >= max(1,N).

Q


Q is REAL array, dimension (LDQ, N)
On entry, if COMPQ = 'V', the orthogonal matrix Q1,
typically from the QR factorization of B.
On exit, if COMPQ='I', the orthogonal matrix Q, and if
COMPQ = 'V', the product Q1*Q.
Not referenced if COMPQ='N'.

LDQ


LDQ is INTEGER
The leading dimension of the array Q.
LDQ >= N if COMPQ='V' or 'I'; LDQ >= 1 otherwise.

Z


Z is REAL array, dimension (LDZ, N)
On entry, if COMPZ = 'V', the orthogonal matrix Z1.
On exit, if COMPZ='I', the orthogonal matrix Z, and if
COMPZ = 'V', the product Z1*Z.
Not referenced if COMPZ='N'.

LDZ


LDZ is INTEGER
The leading dimension of the array Z.
LDZ >= N if COMPZ='V' or 'I'; LDZ >= 1 otherwise.

WORK


WORK is REAL array, dimension (LWORK)
On exit, if INFO = 0, WORK(1) returns the optimal LWORK.

LWORK


LWORK is 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.

INFO


INFO is INTEGER
= 0: successful exit.
< 0: if INFO = -i, the i-th argument had an illegal value.

Author

Univ. of Tennessee

Univ. of California Berkeley

Univ. of Colorado Denver

NAG Ltd.

Further Details:


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

subroutine zgghd3 (character compq, character compz, integer n, integer ilo, integer ihi, complex*16, dimension( lda, * ) a, integer lda, complex*16, dimension( ldb, * ) b, integer ldb, complex*16, dimension( ldq, * ) q, integer ldq, complex*16, dimension( ldz, * ) z, integer ldz, complex*16, dimension( * ) work, integer lwork, integer info)

ZGGHD3

Purpose:


ZGGHD3 reduces a pair of complex matrices (A,B) to generalized upper
Hessenberg form using unitary transformations, where A is a
general matrix and B is upper triangular. The form of the
generalized eigenvalue problem is
A*x = lambda*B*x,
and B is typically made upper triangular by computing its QR
factorization and moving the unitary matrix Q to the left side
of the equation.
This subroutine simultaneously reduces A to a Hessenberg matrix H:
Q**H*A*Z = H
and transforms B to another upper triangular matrix T:
Q**H*B*Z = T
in order to reduce the problem to its standard form
H*y = lambda*T*y
where y = Z**H*x.
The 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 * Z1**H = (Q1*Q) * H * (Z1*Z)**H
Q1 * B * Z1**H = (Q1*Q) * T * (Z1*Z)**H
If Q1 is the unitary matrix from the QR factorization of B in the
original equation A*x = lambda*B*x, then ZGGHD3 reduces the original
problem to generalized Hessenberg form.
This is a blocked variant of CGGHRD, using matrix-matrix
multiplications for parts of the computation to enhance performance.

Parameters

COMPQ


COMPQ is CHARACTER*1
= 'N': do not compute Q;
= 'I': Q is initialized to the unit matrix, and the
unitary matrix Q is returned;
= 'V': Q must contain a unitary matrix Q1 on entry,
and the product Q1*Q is returned.

COMPZ


COMPZ is CHARACTER*1
= 'N': do not compute Z;
= 'I': Z is initialized to the unit matrix, and the
unitary matrix Z is returned;
= 'V': Z must contain a unitary matrix Z1 on entry,
and the product Z1*Z is returned.

N


N is INTEGER
The order of the matrices A and B. N >= 0.

ILO


ILO is INTEGER

IHI


IHI is 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 ZGGBAL; otherwise they
should be set to 1 and N respectively.
1 <= ILO <= IHI <= N, if N > 0; ILO=1 and IHI=0, if N=0.

A


A is COMPLEX*16 array, dimension (LDA, N)
On entry, the N-by-N general matrix to be reduced.
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.

LDA


LDA is INTEGER
The leading dimension of the array A. LDA >= max(1,N).

B


B is COMPLEX*16 array, dimension (LDB, N)
On entry, the N-by-N upper triangular matrix B.
On exit, the upper triangular matrix T = Q**H B Z. The
elements below the diagonal are set to zero.

LDB


LDB is INTEGER
The leading dimension of the array B. LDB >= max(1,N).

Q


Q is COMPLEX*16 array, dimension (LDQ, N)
On entry, if COMPQ = 'V', the unitary matrix Q1, typically
from the QR factorization of B.
On exit, if COMPQ='I', the unitary matrix Q, and if
COMPQ = 'V', the product Q1*Q.
Not referenced if COMPQ='N'.

LDQ


LDQ is INTEGER
The leading dimension of the array Q.
LDQ >= N if COMPQ='V' or 'I'; LDQ >= 1 otherwise.

Z


Z is COMPLEX*16 array, dimension (LDZ, N)
On entry, if COMPZ = 'V', the unitary matrix Z1.
On exit, if COMPZ='I', the unitary matrix Z, and if
COMPZ = 'V', the product Z1*Z.
Not referenced if COMPZ='N'.

LDZ


LDZ is INTEGER
The leading dimension of the array Z.
LDZ >= N if COMPZ='V' or 'I'; LDZ >= 1 otherwise.

WORK


WORK is COMPLEX*16 array, dimension (LWORK)
On exit, if INFO = 0, WORK(1) returns the optimal LWORK.

LWORK


LWORK is 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.

INFO


INFO is INTEGER
= 0: successful exit.
< 0: if INFO = -i, the i-th argument had an illegal value.

Author

Univ. of Tennessee

Univ. of California Berkeley

Univ. of Colorado Denver

NAG Ltd.

Further Details:


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

Author

Generated automatically by Doxygen for LAPACK from the source code.

Wed Feb 7 2024 11:30:40 Version 3.12.0