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

?tgsyl

Solves the generalized Sylvester equation.

Syntax

lapack_int LAPACKE_stgsyl( int matrix_layout, char trans, lapack_int ijob, lapack_int m, lapack_int n, const float* a, lapack_int lda, const float* b, lapack_int ldb, float* c, lapack_int ldc, const float* d, lapack_int ldd, const float* e, lapack_int lde, float* f, lapack_int ldf, float* scale, float* dif );

lapack_int LAPACKE_dtgsyl( int matrix_layout, char trans, lapack_int ijob, lapack_int m, lapack_int n, const double* a, lapack_int lda, const double* b, lapack_int ldb, double* c, lapack_int ldc, const double* d, lapack_int ldd, const double* e, lapack_int lde, double* f, lapack_int ldf, double* scale, double* dif );

lapack_int LAPACKE_ctgsyl( int matrix_layout, char trans, lapack_int ijob, lapack_int m, lapack_int n, const lapack_complex_float* a, lapack_int lda, const lapack_complex_float* b, lapack_int ldb, lapack_complex_float* c, lapack_int ldc, const lapack_complex_float* d, lapack_int ldd, const lapack_complex_float* e, lapack_int lde, lapack_complex_float* f, lapack_int ldf, float* scale, float* dif );

lapack_int LAPACKE_ztgsyl( int matrix_layout, char trans, lapack_int ijob, lapack_int m, lapack_int n, const lapack_complex_double* a, lapack_int lda, const lapack_complex_double* b, lapack_int ldb, lapack_complex_double* c, lapack_int ldc, const lapack_complex_double* d, lapack_int ldd, const lapack_complex_double* e, lapack_int lde, lapack_complex_double* f, lapack_int ldf, double* scale, double* dif );

Include Files

  • mkl.h

Description

The routine solves the generalized Sylvester equation:

A*R-L*B = scale*C

D*R-L*E = scale*F

where R and L are unknown m-by-n matrices, (A, D), (B, E) and (C, F) are given matrix pairs of size m-by-m, n-by-n and m-by-n, respectively, with real/complex entries. (A, D) and (B, E) must be in generalized real-Schur/Schur canonical form, that is, A, B are upper quasi-triangular/triangular and D, E are upper triangular.

The solution (R, L) overwrites (C, F). The factor scale, 0scale1, is an output scaling factor chosen to avoid overflow.

In matrix notation the above equation is equivalent to the following: solve Z*x = scale*b, where Z is defined as


Equation

Here Ik is the identity matrix of size k and XT is the transpose/conjugate-transpose of X. kron(X, Y) is the Kronecker product between the matrices X and Y.

If trans = 'T' (for real flavors), or trans = 'C' (for complex flavors), the routine ?tgsyl solves the transposed/conjugate-transposed system ZT*y = scale*b, which is equivalent to solve for R and L in

AT*R+DT*L = scale*C

R*BT+L*ET = scale*(-F)

This case (trans = 'T' for stgsyl/dtgsyl or trans = 'C' for ctgsyl/ztgsyl) is used to compute an one-norm-based estimate of Dif[(A, D), (B, E)], the separation between the matrix pairs (A,D) and (B,E).

If ijob ≥ 1, ?tgsyl computes a Frobenius norm-based estimate of Dif[(A, D), (B,E)]. That is, the reciprocal of a lower bound on the reciprocal of the smallest singular value of Z. This is a level 3 BLAS algorithm.

Input Parameters

matrix_layout

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

trans

Must be 'N', 'T', or 'C'.

If trans = 'N', solve the generalized Sylvester equation.

If trans = 'T', solve the 'transposed' system (for real flavors only).

If trans = 'C', solve the ' conjugate transposed' system (for complex flavors only).

ijob

Specifies what kind of functionality to be performed:

If ijob =0, solve the generalized Sylvester equation only;

If ijob =1, perform the functionality of ijob =0 and ijob =3;

If ijob =2, perform the functionality of ijob =0 and ijob =4;

If ijob =3, only an estimate of Dif[(A, D), (B, E)] is computed (look ahead strategy is used);

If ijob =4, only an estimate of Dif[(A, D), (B,E)] is computed (?gecon on sub-systems is used). If trans = 'T' or 'C', ijob is not referenced.

m

The order of the matrices A and D, and the row dimension of the matrices C, F, R and L.

n

The order of the matrices B and E, and the column dimension of the matrices C, F, R and L.

a, b, c, d, e, f

Arrays:

a (size max(1, lda*m)) contains the upper quasi-triangular (for real flavors) or upper triangular (for complex flavors) matrix A.

b (size max(1, ldb*n)) contains the upper quasi-triangular (for real flavors) or upper triangular (for complex flavors) matrix B.

c(size max(1, ldc*n) for column major layout and max(1, ldc*m) for row major layout) contains the right-hand-side of the first matrix equation in the generalized Sylvester equation (as defined by trans)

d (size max(1, ldd*m)) contains the upper triangular matrix D.

e (size max(1, lde*n)) contains the upper triangular matrix E.

f(size max(1, ldf*n) for column major layout and max(1, ldf*m) for row major layout) contains the right-hand-side of the second matrix equation in the generalized Sylvester equation (as defined by trans)

lda

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

ldb

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

ldc

The leading dimension of c; at least max(1, m) for column major layout and at least max(1, n) for row major layout .

ldd

The leading dimension of d; at least max(1, m).

lde

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

ldf

The leading dimension of f; at least max(1, m) for column major layout and at least max(1, n) for row major layout .

Output Parameters

c

If ijob=0, 1, or 2, overwritten by the solution R.

If ijob=3 or 4 and trans = 'N', c holds R, the solution achieved during the computation of the Dif-estimate.

f

If ijob=0, 1, or 2, overwritten by the solution L.

If ijob=3 or 4 and trans = 'N', f holds L, the solution achieved during the computation of the Dif-estimate.

dif

On exit, dif is the reciprocal of a lower bound of the reciprocal of the Dif-function, that is, dif is an upper bound of Dif[(A, D), (B, E)] = sigma_min(Z), where Z as defined in the description.

If ijob = 0, or trans = 'T' (for real flavors), or trans = 'C' (for complex flavors), dif is not touched.

scale

On exit, scale is the scaling factor in the generalized Sylvester equation.

If 0 < scale < 1, c and f hold the solutions R and L, respectively, to a slightly perturbed system but the input matrices A, B, D and E have not been changed.

If scale = 0, c and f hold the solutions R and L, respectively, to the homogeneous system with C = F = 0. Normally, scale = 1.

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 > 0, (A, D) and (B, E) have common or close eigenvalues.