1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347
|
<HTML>
<HEAD><TITLE>MB3LZP - SLICOT Library Routine Documentation</TITLE>
</HEAD>
<BODY>
<H2><A Name="MB3LZP">MB3LZP</A></H2>
<H3>
Eigenvalues and right deflating subspace of a complex skew-Hamiltonian/Hamiltonian pencil (applying transformations on panels of columns)
</H3>
<A HREF ="#Specification"><B>[Specification]</B></A>
<A HREF ="#Arguments"><B>[Arguments]</B></A>
<A HREF ="#Method"><B>[Method]</B></A>
<A HREF ="#References"><B>[References]</B></A>
<A HREF ="#Comments"><B>[Comments]</B></A>
<A HREF ="#Example"><B>[Example]</B></A>
<P>
<B><FONT SIZE="+1">Purpose</FONT></B>
<PRE>
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.
</PRE>
<A name="Specification"><B><FONT SIZE="+1">Specification</FONT></B></A>
<PRE>
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( * )
</PRE>
<A name="Arguments"><B><FONT SIZE="+1">Arguments</FONT></B></A>
<P>
<B>Mode Parameters</B>
<PRE>
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.
</PRE>
<B>Input/Output Parameters</B>
<PRE>
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.
</PRE>
<B>Workspace</B>
<PRE>
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 - 1, if COMPQ = 'C'.
</PRE>
<B>Error Indicator</B>
<PRE>
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.
</PRE>
<A name="Method"><B><FONT SIZE="+1">Method</FONT></B></A>
<PRE>
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.
</PRE>
<A name="References"><B><FONT SIZE="+1">References</FONT></B></A>
<PRE>
[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.
</PRE>
<A name="Numerical Aspects"><B><FONT SIZE="+1">Numerical Aspects</FONT></B></A>
<PRE> 3
The algorithm is numerically backward stable and needs O(N )
complex floating point operations.
</PRE>
<A name="Comments"><B><FONT SIZE="+1">Further Comments</FONT></B></A>
<PRE>
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.
</PRE>
<A name="Example"><B><FONT SIZE="+1">Example</FONT></B></A>
<P>
<B>Program Text</B>
<PRE>
None
</PRE>
<B>Program Data</B>
<PRE>
None
</PRE>
<B>Program Results</B>
<PRE>
None
</PRE>
<HR>
<p>
<A HREF=..\libindex.html><B>Return to index</B></A></BODY>
</HTML>
|