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

?ggev

Computes the generalized eigenvalues, and the left and/or right generalized eigenvectors for a pair of nonsymmetric matrices.

Syntax

lapack_int LAPACKE_sggev( 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_dggev( 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_cggev( 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_zggev( 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

The ?ggev routine computes the generalized eigenvalues, and optionally, the left and/or right generalized eigenvectors for a pair of n-by-n real/complex nonsymmetric matrices (A,B).

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.

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

A*v(j) = λ(j)*B*v(j).

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

u(j)H*A = λ(j)*u(j)H*B

where u(j)H denotes the conjugate transpose of u(j).

The ?ggev routine replaces the deprecated ?gegv routine.

Input Parameters

matrix_layout

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

jobvl

Must be 'N' or 'V'.

If jobvl = 'N', the left generalized eigenvectors are not computed;

If jobvl = 'V', the left generalized eigenvectors are computed.

jobvr

Must be 'N' or 'V'.

If jobvr = 'N', the right generalized eigenvectors are not computed;

If jobvr = 'V', the right generalized eigenvectors are computed.

n

The order of the matrices A, B, vl, and vr (n 0).

a, b

Arrays:

a (size at least max(1, lda*n)) is an array containing the n-by-n matrix A (first of the pair of matrices).

b (size at least max(1, ldb*n)) is an array containing the n-by-n matrix B (second of the pair of matrices).

lda

The leading dimension of the array a. Must be at least max(1, n).

ldb

The leading dimension of the array b. Must be at least max(1, n).

ldvl, ldvr

The leading dimensions of the output matrices vl and vr, respectively.

Constraints:

ldvl 1. If jobvl = 'V', ldvl max(1, n).

ldvr 1. If jobvr = 'V', ldvr max(1, n).

Output Parameters

a, b

On exit, these arrays have been overwritten.

alphar, alphai

Arrays, size at least max(1, n) each. Contain values that form generalized eigenvalues in real flavors.

See beta.

alpha

Array, size at least max(1, n). Contain values that form generalized eigenvalues in complex flavors. See beta.

beta

Array, size at least max(1, n).

For real flavors:

On exit, (alphar[j] + alphai[j]*i)/beta[j], j=0,..., n - 1, are the generalized eigenvalues.

If alphai[j] 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+1] negative.

For complex flavors:

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

See also Application Notes below.

vl, vr

Arrays:

vl (size at least max(1, ldvl*n)). Contains the matrix of left generalized eigenvectors VL.

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(Re) + abs(Im) = 1.

If jobvl = 'N', vl is not referenced.

For real flavors:

If the j-th eigenvalue is real,then the k-th component of the j-th left eigenvector uj is stored in vl[(k - 1) + (j - 1)*ldvl] for column major layout and in vl[(k - 1)*ldvl + (j - 1)] for row major layout..

If the j-th and (j+1)-st eigenvalues form a complex conjugate pair, then for i = sqrt(-1), the k-th components of the j-th left eigenvector ujare vl[(k - 1) + (j - 1)*ldvl] + i*vl[(k - 1) + j*ldvl] for column major layout and vl[(k - 1)*ldvl + (j - 1)] + i*vl[(k - 1)*ldvl + j] for row major layout. Similarly, the k-th components of left eigenvector j+1 uj+1 are vl[(k - 1) + (j - 1)*ldvl] - i*vl[(k - 1) + j*ldvl] for column major layout and vl[(k - 1)*ldvl + (j - 1)] - i*vl[(k - 1)*ldvl + j] for row major layout..

For complex flavors:

The k-th component of the j-th left eigenvector uj is stored in vl[(k - 1) + (j - 1)*ldvl] for column major layout and in vl[(k - 1)*ldvl + (j - 1)] for row major layout.

vr (size at least max(1, ldvr*n)). Contains the matrix of right generalized eigenvectors VR.

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(Re) + abs(Im) = 1.

If jobvr = 'N', vr is not referenced.

For real flavors:

If the j-th eigenvalue is real, then The k-th component of the j-th right eigenvector vj is stored in vr[(k - 1) + (j - 1)*ldvr] for column major layout and in vr[(k - 1)*ldvr + (j - 1)] for row major layout..

If the j-th and (j+1)-st eigenvalues form a complex conjugate pair, then the k-th components of thej-th right eigenvector vj can be computed as vr[(k - 1) + (j - 1)*ldvr] + i*vr[(k - 1) + j*ldvr] for column major layout and vr[(k - 1)*ldvr + (j - 1)] + i*vr[(k - 1)*ldvr + j] for row major layout. Similarly, the k-th components of the right eigenvector j+1 v{j+1} can be computed as vr[(k - 1) + (j - 1)*ldvr] - i*vr[(k - 1) + j*ldvr] for column major layout and vr[(k - 1)*ldvr + (j - 1)] - i*vr[(k - 1)*ldvr + j] for row major layout..

For complex flavors:

The k-th component of the j-th right eigenvector vj is stored in vr[(k - 1) + (j - 1)*ldvr] for column major layout and in vr[(k - 1)*ldvr + (j - 1)] for row major layout.

Return Values

This function returns a value info.

If info=0, the execution is successful.

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

If info = i, and

in: the QZ iteration failed. No eigenvectors have been calculated, but alphar[j], alphai[j] (for real flavors), or alpha[j] (for complex flavors), and beta[j], j=info,..., n - 1 should be correct.

i > n: errors that usually indicate LAPACK problems:

i = n+1: other than QZ iteration failed in hgeqz;

i = n+2: error return from tgevc.

Application Notes

The quotients alphar[j]/beta[j] and alphai[j]/beta[j] may easily over- or underflow, and beta[j] may even be zero. Thus, you should avoid simply computing the ratio. However, alphar and alphai (for real flavors) or alpha (for complex flavors) will be always less than and usually comparable with norm(A) in magnitude, and beta always less than and usually comparable with norm(B).