Purpose
To compute the eigenvalues of a complex N-by-N skew-Hamiltonian/ Hamiltonian pencil aS - bH, with ( A D ) ( B F ) S = ( ) and H = ( ). (1) ( E A' ) ( G -B' ) The structured Schur form of the embedded real skew-Hamiltonian/ skew-Hamiltonian pencil aB_S - bB_T, defined as ( Re(A) -Im(A) | Re(D) -Im(D) ) ( | ) ( Im(A) Re(A) | Im(D) Re(D) ) ( | ) B_S = (-----------------+-----------------) , and ( | ) ( Re(E) -Im(E) | Re(A') Im(A') ) ( | ) ( Im(E) Re(E) | -Im(A') Re(A') ) (2) ( -Im(B) -Re(B) | -Im(F) -Re(F) ) ( | ) ( Re(B) -Im(B) | Re(F) -Im(F) ) ( | ) B_T = (-----------------+-----------------) , T = i*H, ( | ) ( -Im(G) -Re(G) | -Im(B') Re(B') ) ( | ) ( Re(G) -Im(G) | -Re(B') -Im(B') ) is determined and used to compute the eigenvalues. The notation M' denotes the conjugate transpose of the matrix M. Optionally, if COMPQ = 'C', an orthonormal basis of the right deflating subspace of the pencil aS - bH, corresponding to the eigenvalues with strictly negative real part, is computed. Namely, after transforming aB_S - bB_H by unitary matrices, we have ( BA BD ) ( BB BF ) B_Sout = ( ) and B_Hout = ( ), (3) ( 0 BA' ) ( 0 -BB' ) and the eigenvalues with strictly negative real part of the complex pencil aB_Sout - bB_Hout are moved to the top. The embedding doubles the multiplicities of the eigenvalues of the pencil aS - bH.Specification
SUBROUTINE MB3LZP( COMPQ, ORTH, N, A, LDA, DE, LDDE, B, LDB, FG, $ LDFG, NEIG, Q, LDQ, ALPHAR, ALPHAI, BETA, $ IWORK, DWORK, LDWORK, ZWORK, LZWORK, BWORK, $ INFO ) C .. Scalar Arguments .. CHARACTER COMPQ, ORTH INTEGER INFO, LDA, LDB, LDDE, LDFG, LDQ, LDWORK, $ LZWORK, N, NEIG C .. Array Arguments .. LOGICAL BWORK( * ) INTEGER IWORK( * ) DOUBLE PRECISION ALPHAI( * ), ALPHAR( * ), BETA( * ), DWORK( * ) COMPLEX*16 A( LDA, * ), B( LDB, * ), DE( LDDE, * ), $ FG( LDFG, * ), Q( LDQ, * ), ZWORK( * )Arguments
Mode Parameters
COMPQ CHARACTER*1 Specifies whether to compute the deflating subspace corresponding to the eigenvalues of aS - bH with strictly negative real part. = 'N': do not compute the deflating subspace; compute the eigenvalues only; = 'C': compute the deflating subspace and store it in the leading subarray of Q. ORTH CHARACTER*1 If COMPQ = 'C', specifies the technique for computing an orthonormal basis of the deflating subspace, as follows: = 'P': QR factorization with column pivoting; = 'S': singular value decomposition. If COMPQ = 'N', the ORTH value is not used.Input/Output Parameters
N (input) INTEGER The order of the pencil aS - bH. N >= 0, even. A (input/output) COMPLEX*16 array, dimension (LDA, N) On entry, the leading N/2-by-N/2 part of this array must contain the matrix A. On exit, if COMPQ = 'C', the leading N-by-N part of this array contains the upper triangular matrix BA in (3) (see also METHOD). The strictly lower triangular part is not zeroed; it is preserved in the leading N/2-by-N/2 part. If COMPQ = 'N', this array is unchanged on exit. LDA INTEGER The leading dimension of the array A. LDA >= MAX(1, N). DE (input/output) COMPLEX*16 array, dimension (LDDE, N) On entry, the leading N/2-by-N/2 lower triangular part of this array must contain the lower triangular part of the skew-Hermitian matrix E, and the N/2-by-N/2 upper triangular part of the submatrix in the columns 2 to N/2+1 of this array must contain the upper triangular part of the skew-Hermitian matrix D. On exit, if COMPQ = 'C', the leading N-by-N part of this array contains the skew-Hermitian matrix BD in (3) (see also METHOD). The strictly lower triangular part of the input matrix is preserved. If COMPQ = 'N', this array is unchanged on exit. LDDE INTEGER The leading dimension of the array DE. LDDE >= MAX(1, N). B (input/output) COMPLEX*16 array, dimension (LDB, N) On entry, the leading N/2-by-N/2 part of this array must contain the matrix B. On exit, if COMPQ = 'C', the leading N-by-N part of this array contains the upper triangular matrix BB in (3) (see also METHOD). The strictly lower triangular part is not zeroed; the elements below the first subdiagonal of the input matrix are preserved. If COMPQ = 'N', this array is unchanged on exit. LDB INTEGER The leading dimension of the array B. LDB >= MAX(1, N). FG (input/output) COMPLEX*16 array, dimension (LDFG, N) On entry, the leading N/2-by-N/2 lower triangular part of this array must contain the lower triangular part of the Hermitian matrix G, and the N/2-by-N/2 upper triangular part of the submatrix in the columns 2 to N/2+1 of this array must contain the upper triangular part of the Hermitian matrix F. On exit, if COMPQ = 'C', the leading N-by-N part of this array contains the Hermitian matrix BF in (3) (see also METHOD). The strictly lower triangular part of the input matrix is preserved. The diagonal elements might have tiny imaginary parts. If COMPQ = 'N', this array is unchanged on exit. LDFG INTEGER The leading dimension of the array FG. LDFG >= MAX(1, N). NEIG (output) INTEGER If COMPQ = 'C', the number of eigenvalues in aS - bH with strictly negative real part. Q (output) COMPLEX*16 array, dimension (LDQ, 2*N) On exit, if COMPQ = 'C', the leading N-by-NEIG part of this array contains an orthonormal basis of the right deflating subspace corresponding to the eigenvalues of the pencil aS - bH with strictly negative real part. The remaining entries are meaningless. If COMPQ = 'N', this array is not referenced. LDQ INTEGER The leading dimension of the array Q. LDQ >= 1, if COMPQ = 'N'; LDQ >= MAX(1, 2*N), if COMPQ = 'C'. ALPHAR (output) DOUBLE PRECISION array, dimension (N) The real parts of each scalar alpha defining an eigenvalue of the pencil aS - bH. ALPHAI (output) DOUBLE PRECISION array, dimension (N) The imaginary parts of each scalar alpha defining an eigenvalue of the pencil aS - bH. If ALPHAI(j) is zero, then the j-th eigenvalue is real. BETA (output) DOUBLE PRECISION array, dimension (N) The scalars beta that define the eigenvalues of the pencil aS - bH. Together, the quantities alpha = (ALPHAR(j),ALPHAI(j)) and beta = BETA(j) represent the j-th eigenvalue of the pencil aS - bH, in the form lambda = alpha/beta. Since lambda may overflow, the ratios should not, in general, be computed.Workspace
IWORK INTEGER array, dimension (N+1) DWORK DOUBLE PRECISION array, dimension (LDWORK) On exit, if INFO = 0, DWORK(1) returns the optimal LDWORK. On exit, if INFO = -20, DWORK(1) returns the minimum value of LDWORK. LDWORK INTEGER The dimension of the array DWORK. LDWORK >= MAX( 4*N*N + 2*N + MAX(3,N) ), if COMPQ = 'N'; LDWORK >= MAX( 1, 11*N*N + 2*N ), if COMPQ = 'C'. For good performance LDWORK should be generally larger. If LDWORK = -1, then a workspace query is assumed; the routine only calculates the optimal size of the DWORK array, returns this value as the first entry of the DWORK array, and no error message related to LDWORK is issued by XERBLA. ZWORK COMPLEX*16 array, dimension (LZWORK) On exit, if INFO = 0, ZWORK(1) returns the optimal LZWORK. On exit, if INFO = -22, ZWORK(1) returns the minimum value of LZWORK. LZWORK INTEGER The dimension of the array ZWORK. LZWORK >= 1, if COMPQ = 'N'; LZWORK >= 8*N + 4, if COMPQ = 'C'. For good performance LZWORK should be generally larger. If LZWORK = -1, then a workspace query is assumed; the routine only calculates the optimal size of the ZWORK array, returns this value as the first entry of the ZWORK array, and no error message related to LZWORK is issued by XERBLA. BWORK LOGICAL array, dimension (LBWORK) LBWORK >= 0, if COMPQ = 'N'; LBWORK >= N, if COMPQ = 'C'.Error Indicator
INFO INTEGER = 0: succesful exit; < 0: if INFO = -i, the i-th argument had an illegal value; = 1: QZ iteration failed in the SLICOT Library routine MB04FP (QZ iteration did not converge or computation of the shifts failed); = 2: QZ iteration failed in the LAPACK routine ZHGEQZ when trying to triangularize the 2-by-2 blocks; = 3: the singular value decomposition failed in the LAPACK routine ZGESVD (for ORTH = 'S'); = 4: warning: the pencil is numerically singular.Method
First, T = i*H is set. Then, the embeddings, B_S and B_T, of the matrices S and T, are determined and, subsequently, the SLICOT Library routine MB04FP is applied to compute the structured Schur form, i.e., the factorizations ~ ( S11 S12 ) B_S = J Q' J' B_S Q = ( ) and ( 0 S11' ) ~ ( T11 T12 ) ( 0 I ) B_T = J Q' J' B_T Q = ( ), with J = ( ), ( 0 T11' ) ( -I 0 ) where Q is real orthogonal, S11 is upper triangular, and T11 is upper quasi-triangular. Second, the SLICOT Library routine MB3JZP is applied, to compute a ~ unitary matrix Q, such that ~ ~ ~ ~ ~ ( S11 S12 ) J Q' J' B_S Q = ( ~ ) =: B_Sout, ( 0 S11' ) ~ ~ ~ ( H11 H12 ) J Q' J'(-i*B_T) Q = ( ) =: B_Hout, ( 0 -H11' ) ~ ~ ~ with S11, H11 upper triangular, and such that Spec_-(B_S, -i*B_T) is contained in the spectrum of the 2*NEIG-by-2*NEIG leading ~ principal subpencil aS11 - bH11. Finally, the right deflating subspace is computed. See also page 22 in [1] for more details.References
[1] Benner, P., Byers, R., Mehrmann, V. and Xu, H. Numerical Computation of Deflating Subspaces of Embedded Hamiltonian Pencils. Tech. Rep. SFB393/99-15, Technical University Chemnitz, Germany, June 1999.Numerical Aspects
3 The algorithm is numerically backward stable and needs O(N ) complex floating point operations.Further Comments
This routine does not perform any scaling of the matrices. Scaling might sometimes be useful, and it should be done externally. For large values of N, the routine applies the transformations on panels of columns. The user may specify in INFO the desired number of columns. If on entry INFO <= 0, then the routine estimates a suitable value of this number.Example
Program Text
NoneProgram Data
NoneProgram Results
None