Developer Reference for Intel® oneAPI Math Kernel Library for C

ID 766684
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

?ggev3

Computes the generalized eigenvalues and the left and right generalized eigenvectors for a pair of matrices.

Syntax

lapack_int LAPACKE_sggev3 (int matrix_layout, char jobvl, char jobvr, lapack_int n, float * a, lapack_int lda, float * b, lapack_int ldb, float * alphar, float * alphai, float * beta, float * vl, lapack_int ldvl, float * vr, lapack_int ldvr);

lapack_int LAPACKE_dggev3 (int matrix_layout, char jobvl, char jobvr, lapack_int n, double * a, lapack_int lda, double * b, lapack_int ldb, double * alphar, double * alphai, double * beta, double * vl, lapack_int ldvl, double * vr, lapack_int ldvr);

lapack_int LAPACKE_cggev3 (int matrix_layout, char jobvl, char jobvr, lapack_int n, lapack_complex_float * a, lapack_int lda, lapack_complex_float * b, lapack_int ldb, lapack_complex_float * alpha, lapack_complex_float * beta, lapack_complex_float * vl, lapack_int ldvl, lapack_complex_float * vr, lapack_int ldvr);

lapack_int LAPACKE_zggev3 (int matrix_layout, char jobvl, char jobvr, lapack_int n, lapack_complex_double * a, lapack_int lda, lapack_complex_double * b, lapack_int ldb, lapack_complex_double * alpha, lapack_complex_double * beta, lapack_complex_double * vl, lapack_int ldvl, lapack_complex_double * vr, lapack_int ldvr);

Include Files

  • mkl.h

Description

For a pair of n-by-n real or complex nonsymmetric matrices (A, B), ?ggev3 computes the generalized eigenvalues, and optionally, the left and right generalized eigenvectors.

A generalized eigenvalue for a pair of matrices (A, B) is a scalar λ or a ratio alpha/beta = λ, such that A - λ*B is singular. It is usually represented as the pair (alpha,beta), as there is a reasonable interpretation for beta=0, and even for both being zero.

For real flavors:

The right eigenvector vj corresponding to the eigenvalue λj of (A, B) satisfies

A * vj = λj * B * vj.

The left eigenvector uj corresponding to the eigenvalue λj of (A, B) satisfies

ujH * A = λj * ujH * B

where ujH is the conjugate-transpose of uj.

For complex flavors:

The right generalized eigenvector vj corresponding to the generalized eigenvalue λj of (A, B) satisfies

A * vj = λj * B * vj.

The left generalized eigenvector uj corresponding to the generalized eigenvalues λj of (A, B) satisfies

ujH * A = λj * ujH * B

where ujH is the conjugate-transpose of uj.

Input Parameters

matrix_layout

Specifies whether matrix storage layout is row major (LAPACK_ROW_MAJOR) or column major (LAPACK_COL_MAJOR).

jobvl

= 'N': do not compute the left generalized eigenvectors;

= 'V': compute the left generalized eigenvectors.

jobvr

= 'N': do not compute the right generalized eigenvectors;

= 'V': compute the right generalized eigenvectors.

n

The order of the matrices A, B, VL, and VR.

n 0.

a

Array, size (lda*n).

On entry, the matrix A in the pair (A, B).

lda

The leading dimension of a.

lda max(1,n).

b

Array, size (ldb*n).

On entry, the matrix B in the pair (A, B).

ldb

The leading dimension of b.

ldb max(1,n).

ldvl

The leading dimension of the matrix VL.

ldvl 1, and if jobvl = 'V', ldvln.

ldvr

The leading dimension of the matrix VR.

ldvr 1, and if jobvr = 'V', ldvrn.

Output Parameters

a

On exit, a is overwritten.

b

On exit, b is overwritten.

alphar

Array, size (n).

alphai

Array, size (n).

alpha

Array, size (n).

beta

Array, size (n).

For real flavors:

On exit, (alphar[j] + alphai[j]*i)/beta[j], j=0,...,n - 1, are the generalized eigenvalues. If alphai[j - 1] is zero, then the j-th eigenvalue is real; if positive, then the j-th and (j+1)-st eigenvalues are a complex conjugate pair, with alphai[j] negative.

Note: the quotients alphar[j - 1]/beta[j - 1] and alphai[j - 1]/beta[j - 1] can easily over- or underflow, and beta(j) might even be zero. Thus, you should avoid computing the ratio alpha/beta by simply dividing alpha by beta. However, alphar and alphai are always less than and usually comparable with norm(A) in magnitude, and beta is always less than and usually comparable with norm(B).

For complex flavors:

On exit, alpha[j]/beta[j], j=0,...,n - 1, are the generalized eigenvalues.

Note: the quotients alpha[j - 1]/beta[j - 1] may easily over- or underflow, and beta(j) can even be zero. Thus, you should avoid computing the ratio alpha/beta by simply dividing alpha by beta. However, alpha is always less than and usually comparable with norm(A) in magnitude, and betais always less than and usually comparable with norm(B).

vl

Array, size (ldvl*n).

For real flavors:

If jobvl = 'V', the left eigenvectors uj are stored one after another in the columns of vl, in the same order as their eigenvalues. If the j-th eigenvalue is real, then uj = the j-th column of vl. If the j-th and (j+1)-st eigenvalues form a complex conjugate pair, then the real part of uj = the j-th column of vl and the imaginary part of vj = the (j + 1)-st column of vl.

Each eigenvector is scaled so the largest component has abs(real part)+abs(imag. part)=1.

Not referenced if jobvl = 'N'.

For complex flavors:

If jobvl = 'V', the left generalized eigenvectors uj are stored one after another in the columns of vl, in the same order as their eigenvalues.

Each eigenvector is scaled so the largest component has abs(real part) + abs(imag. part) = 1.

Not referenced if jobvl = 'N'.

vr

Array, size (ldvr*n).

For real flavors:

If jobvr = 'V', the right eigenvectors vj are stored one after another in the columns of vr, in the same order as their eigenvalues. If the j-th eigenvalue is real, then vj = the j-th column of vr. If the j-th and (j + 1)-st eigenvalues form a complex conjugate pair, then the real part of vj = the j-th column of vr and the imaginary part of vj = the (j + 1)-st column of vr.

Each eigenvector is scaled so the largest component has abs(real part)+abs(imag. part)=1.

Not referenced if jobvr = 'N'.

For complex flavors:

If jobvr = 'V', the right generalized eigenvectors vj are stored one after another in the columns of vr, in the same order as their eigenvalues. Each eigenvector is scaled so the largest component has abs(real part) + abs(imag. part) = 1.

Not referenced if jobvr = 'N'.

Return Values

This function returns a value info.

= 0: successful exit

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

=1,...,n:

for real flavors:

The QZ iteration failed. No eigenvectors have been calculated, but alphar[j], alphar[j] and beta[j] should be correct for j=info,...,n - 1.

for complex flavors:

The QZ iteration failed. No eigenvectors have been calculated, but alpha[j] and beta[j] should be correct for j=info,...,n - 1.

> n:

=n + 1: other than QZ iteration failed in ?hgeqz,

=n + 2: error return from ?tgevc.