1STGSEN(1)               LAPACK routine (version 3.1.1)               STGSEN(1)
2
3
4

NAME

6       STGSEN - the generalized real Schur decomposition of a real matrix pair
7       (A, B) (in terms of an orthonormal equivalence trans-  formation  Q'  *
8       (A,  B)  * Z), so that a selected cluster of eigenvalues appears in the
9       leading diagonal blocks of the upper quasi-triangular matrix A and  the
10       upper triangular B
11

SYNOPSIS

13       SUBROUTINE STGSEN( IJOB,  WANTQ,  WANTZ,  SELECT,  N,  A,  LDA, B, LDB,
14                          ALPHAR, ALPHAI, BETA, Q, LDQ, Z,  LDZ,  M,  PL,  PR,
15                          DIF, WORK, LWORK, IWORK, LIWORK, INFO )
16
17           LOGICAL        WANTQ, WANTZ
18
19           INTEGER        IJOB, INFO, LDA, LDB, LDQ, LDZ, LIWORK, LWORK, M, N
20
21           REAL           PL, PR
22
23           LOGICAL        SELECT( * )
24
25           INTEGER        IWORK( * )
26
27           REAL           A(  LDA, * ), ALPHAI( * ), ALPHAR( * ), B( LDB, * ),
28                          BETA( * ), DIF( * ), Q( LDQ, * ), WORK( * ), Z( LDZ,
29                          * )
30

PURPOSE

32       STGSEN  reorders  the  generalized  real  Schur decomposition of a real
33       matrix pair (A, B) (in terms of an orthonormal equivalence trans-  for‐
34       mation  Q'  *  (A,  B)  * Z), so that a selected cluster of eigenvalues
35       appears in the leading diagonal blocks of  the  upper  quasi-triangular
36       matrix  A  and  the  upper triangular B. The leading columns of Q and Z
37       form orthonormal bases of the corresponding left and right eigen-  spa‐
38       ces  (deflating  subspaces).  (A,  B) must be in generalized real Schur
39       canonical form (as returned by SGGES), i.e. A is block upper triangular
40       with 1-by-1 and 2-by-2 diagonal blocks. B is upper triangular.
41
42       STGSEN also computes the generalized eigenvalues
43
44                   w(j) = (ALPHAR(j) + i*ALPHAI(j))/BETA(j)
45
46       of the reordered matrix pair (A, B).
47
48       Optionally,  STGSEN computes the estimates of reciprocal condition num‐
49       bers  for  eigenvalues  and  eigenspaces.  These  are   Difu[(A11,B11),
50       (A22,B22)]  and  Difl[(A11,B11),  (A22,B22)],  i.e.  the  separation(s)
51       between the matrix pairs (A11, B11) and (A22,B22)  that  correspond  to
52       the  selected  cluster  and the eigenvalues outside the cluster, resp.,
53       and norms of "projections" onto left and right eigenspaces w.r.t.   the
54       selected cluster in the (1,1)-block.
55
56

ARGUMENTS

58       IJOB    (input) INTEGER
59               Specifies  whether condition numbers are required for the clus‐
60               ter of eigenvalues (PL and PR) or the deflating subspaces (Difu
61               and Difl):
62               =0: Only reorder w.r.t. SELECT. No extras.
63               =1:  Reciprocal  of  norms of "projections" onto left and right
64               eigenspaces w.r.t. the selected cluster (PL and PR).  =2: Upper
65               bounds on Difu and Difl. F-norm-based estimate
66               (DIF(1:2)).
67               =3: Estimate of Difu and Difl. 1-norm-based estimate
68               (DIF(1:2)).   About 5 times as expensive as IJOB = 2.  =4: Com‐
69               pute PL, PR and DIF (i.e. 0, 1 and 2 above):  Economic  version
70               to  get  it  all.   =5: Compute PL, PR and DIF (i.e. 0, 1 and 3
71               above)
72
73       WANTQ   (input) LOGICAL
74
75       WANTZ   (input) LOGICAL
76
77       SELECT  (input) LOGICAL array, dimension (N)
78               SELECT specifies the eigenvalues in the selected  cluster.   To
79               select  a  real  eigenvalue w(j), SELECT(j) must be set to w(j)
80               and w(j+1), corresponding to a 2-by-2  diagonal  block,  either
81               SELECT(j)  or  SELECT(j+1)  or  both must be set to either both
82               included in the cluster or both excluded.
83
84       N       (input) INTEGER
85               The order of the matrices A and B. N >= 0.
86
87       A       (input/output) REAL array, dimension(LDA,N)
88               On entry, the upper quasi-triangular matrix A, with (A,  B)  in
89               generalized real Schur canonical form.  On exit, A is overwrit‐
90               ten by the reordered matrix A.
91
92       LDA     (input) INTEGER
93               The leading dimension of the array A. LDA >= max(1,N).
94
95       B       (input/output) REAL array, dimension(LDB,N)
96               On entry, the upper triangular matrix B, with (A, B) in  gener‐
97               alized real Schur canonical form.  On exit, B is overwritten by
98               the reordered matrix B.
99
100       LDB     (input) INTEGER
101               The leading dimension of the array B. LDB >= max(1,N).
102
103       ALPHAR  (output) REAL array, dimension (N)
104               ALPHAI  (output) REAL array,  dimension  (N)  BETA     (output)
105               REAL    array,    dimension   (N)   On   exit,   (ALPHAR(j)   +
106               ALPHAI(j)*i)/BETA(j), j=1,...,N, will be the generalized eigen‐
107               values.  ALPHAR(j) + ALPHAI(j)*i and BETA(j),j=1,...,N  are the
108               diagonals of the complex Schur form (S,T) that would result  if
109               the  2-by-2  diagonal blocks of the real generalized Schur form
110               of (A,B) were further reduced to triangular form using  complex
111               unitary  transformations.   If ALPHAI(j) is zero, then the j-th
112               eigenvalue is real; if positive, then the j-th and (j+1)-st ei‐
113               genvalues  are a complex conjugate pair, with ALPHAI(j+1) nega‐
114               tive.
115
116       Q       (input/output) REAL array, dimension (LDQ,N)
117               On entry, if WANTQ = .TRUE., Q is an N-by-N matrix.  On exit, Q
118               has  been  postmultiplied by the left orthogonal transformation
119               matrix which reorder (A, B); The leading M columns  of  Q  form
120               orthonormal  bases  for  the specified pair of left eigenspaces
121               (deflating subspaces).  If WANTQ = .FALSE.,  Q  is  not  refer‐
122               enced.
123
124       LDQ     (input) INTEGER
125               The leading dimension of the array Q.  LDQ >= 1; and if WANTQ =
126               .TRUE., LDQ >= N.
127
128       Z       (input/output) REAL array, dimension (LDZ,N)
129               On entry, if WANTZ = .TRUE., Z is an N-by-N matrix.  On exit, Z
130               has  been  postmultiplied by the left orthogonal transformation
131               matrix which reorder (A, B); The leading M columns  of  Z  form
132               orthonormal  bases  for  the specified pair of left eigenspaces
133               (deflating subspaces).  If WANTZ = .FALSE.,  Z  is  not  refer‐
134               enced.
135
136       LDZ     (input) INTEGER
137               The  leading  dimension  of  the  array Z. LDZ >= 1; If WANTZ =
138               .TRUE., LDZ >= N.
139
140       M       (output) INTEGER
141               The dimension of the specified pair of left  and  right  eigen-
142               spaces (deflating subspaces). 0 <= M <= N.
143
144       PL      (output) REAL
145               PR       (output)  REAL  If  IJOB = 1, 4 or 5, PL, PR are lower
146               bounds on the reciprocal of the norm of "projections" onto left
147               and  right eigenspaces with respect to the selected cluster.  0
148               < PL, PR <= 1.  If M = 0 or M = N, PL = PR  = 1.  If IJOB =  0,
149               2 or 3, PL and PR are not referenced.
150
151       DIF     (output) REAL array, dimension (2).
152               If IJOB >= 2, DIF(1:2) store the estimates of Difu and Difl.
153               If IJOB = 2 or 4, DIF(1:2) are F-norm-based upper bounds on
154               Difu  and  Difl.  If  IJOB  = 3 or 5, DIF(1:2) are 1-norm-based
155               estimates of Difu and Difl.  If M =  0  or  N,  DIF(1:2)  =  F-
156               norm([A, B]).  If IJOB = 0 or 1, DIF is not referenced.
157
158       WORK    (workspace/output) REAL array, dimension (MAX(1,LWORK))
159               On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
160
161       LWORK   (input) INTEGER
162               The  dimension  of the array WORK. LWORK >=  4*N+16.  If IJOB =
163               1, 2 or 4, LWORK >= MAX(4*N+16, 2*M*(N-M)).  If IJOB = 3 or  5,
164               LWORK >= MAX(4*N+16, 4*M*(N-M)).
165
166               If  LWORK  = -1, then a workspace query is assumed; the routine
167               only calculates the optimal size of  the  WORK  array,  returns
168               this  value  as the first entry of the WORK array, and no error
169               message related to LWORK is issued by XERBLA.
170
171       IWORK   (workspace/output) INTEGER array, dimension (MAX(1,LIWORK))
172               IF IJOB = 0, IWORK is not referenced.  Otherwise, on  exit,  if
173               INFO = 0, IWORK(1) returns the optimal LIWORK.
174
175       LIWORK  (input) INTEGER
176               The  dimension of the array IWORK. LIWORK >= 1.  If IJOB = 1, 2
177               or 4, LIWORK >=  N+6.  If IJOB = 3 or 5, LIWORK >=  MAX(2*M*(N-
178               M), N+6).
179
180               If  LIWORK = -1, then a workspace query is assumed; the routine
181               only calculates the optimal size of the  IWORK  array,  returns
182               this  value as the first entry of the IWORK array, and no error
183               message related to LIWORK is issued by XERBLA.
184
185       INFO    (output) INTEGER
186               =0: Successful exit.
187               <0: If INFO = -i, the i-th argument had an illegal value.
188               =1: Reordering of (A, B) failed because the transformed  matrix
189               pair  (A,  B) would be too far from generalized Schur form; the
190               problem is very ill-conditioned.  (A, B)  may  have  been  par‐
191               tially  reordered.   If  requested, 0 is returned in DIF(*), PL
192               and PR.
193

FURTHER DETAILS

195       STGSEN first collects the selected eigenvalues by computing  orthogonal
196       U  and  W  that  move  them to the top left corner of (A, B).  In other
197       words, the selected eigenvalues are the eigenvalues of (A11, B11) in:
198
199                     U'*(A, B)*W = (A11 A12) (B11 B12) n1
200                                   ( 0  A22),( 0  B22) n2
201                                     n1  n2    n1  n2
202
203       where N = n1+n2 and U' means the transpose of U. The first  n1  columns
204       of  U  and  W  span  the  specified  pair of left and right eigenspaces
205       (deflating subspaces) of (A, B).
206
207       If (A, B) has been obtained from the generalized real Schur  decomposi‐
208       tion of a matrix pair (C, D) = Q*(A, B)*Z', then the reordered general‐
209       ized real Schur form of (C, D) is given by
210
211                (C, D) = (Q*U)*(U'*(A, B)*W)*(Z*W)',
212
213       and the first n1 columns of Q*U and Z*W span the corresponding  deflat‐
214       ing subspaces of (C, D) (Q and Z store Q*U and Z*W, resp.).
215
216       Note  that  if the selected eigenvalue is sufficiently ill-conditioned,
217       then its value may differ significantly from its value before  reorder‐
218       ing.
219
220       The  reciprocal  condition  numbers  of  the left and right eigenspaces
221       spanned by the first n1 columns of U and W (or  Q*U  and  Z*W)  may  be
222       returned in DIF(1:2), corresponding to Difu and Difl, resp.
223
224       The Difu and Difl are defined as:
225
226            Difu[(A11, B11), (A22, B22)] = sigma-min( Zu )
227       and
228            Difl[(A11, B11), (A22, B22)] = Difu[(A22, B22), (A11, B11)],
229
230       where   sigma-min(Zu)   is   the   smallest   singular   value  of  the
231       (2*n1*n2)-by-(2*n1*n2) matrix
232
233            Zu = [ kron(In2, A11)  -kron(A22', In1) ]
234                 [ kron(In2, B11)  -kron(B22', In1) ].
235
236       Here, Inx is the identity matrix of size nx and A22' is  the  transpose
237       of  A22. kron(X, Y) is the Kronecker product between the matrices X and
238       Y.
239
240       When DIF(2) is small, small changes in (A, B) can cause  large  changes
241       in  the  deflating  subspace.  An approximate (asymptotic) bound on the
242       maximum angular error in the computed deflating subspaces is
243
244            EPS * norm((A, B)) / DIF(2),
245
246       where EPS is the machine precision.
247
248       The reciprocal norm of the projectors on the left and right eigenspaces
249       associated with (A11, B11) may be returned in PL and PR.  They are com‐
250       puted as follows. First we compute L and R so that P*(A, B)*Q is  block
251       diagonal, where
252
253            P = ( I -L ) n1           Q = ( I R ) n1
254                ( 0  I ) n2    and        ( 0 I ) n2
255                  n1 n2                    n1 n2
256
257       and (L, R) is the solution to the generalized Sylvester equation
258
259            A11*R - L*A22 = -A12
260            B11*R - L*B22 = -B12
261
262       Then  PL  = (F-norm(L)**2+1)**(-1/2) and PR = (F-norm(R)**2+1)**(-1/2).
263       An approximate (asymptotic) bound on the average absolute error of  the
264       selected eigenvalues is
265
266            EPS * norm((A, B)) / PL.
267
268       There  are also global error bounds which valid for perturbations up to
269       a certain restriction:  A lower bound (x) on the  smallest  F-norm(E,F)
270       for which an eigenvalue of (A11, B11) may move and coalesce with an ei‐
271       genvalue of (A22, B22) under perturbation (E,F), (i.e. (A + E, B +  F),
272       is
273
274        x = min(Difu,Difl)/((1/(PL*PL)+1/(PR*PR))**(1/2)+2*max(1/PL,1/PR)).
275
276       An approximate bound on x can be computed from DIF(1:2), PL and PR.
277
278       If  y  = ( F-norm(E,F) / x) <= 1, the angles between the perturbed (L',
279       R') and unperturbed (L, R) left and right deflating  subspaces  associ‐
280       ated with the selected cluster in the (1,1)-blocks can be bounded as
281
282        max-angle(L, L') <= arctan( y * PL / (1 - y * (1 - PL * PL)**(1/2))
283        max-angle(R, R') <= arctan( y * PR / (1 - y * (1 - PR * PR)**(1/2))
284
285       See  LAPACK  User's  Guide section 4.11 or the following references for
286       more information.
287
288       Note that if the default method for computing the Frobenius-norm- based
289       estimate DIF is not wanted (see SLATDF), then the parameter IDIFJB (see
290       below) should be changed from 3 to 4 (routine SLATDF (IJOB = 2 will  be
291       used)). See STGSYL for more details.
292
293       Based on contributions by
294          Bo Kagstrom and Peter Poromaa, Department of Computing Science,
295          Umea University, S-901 87 Umea, Sweden.
296
297       References
298       ==========
299
300       [1] B. Kagstrom; A Direct Method for Reordering Eigenvalues in the
301           Generalized Real Schur Form of a Regular Matrix Pair (A, B), in
302           M.S. Moonen et al (eds), Linear Algebra for Large Scale and
303           Real-Time Applications, Kluwer Academic Publ. 1993, pp 195-218.
304
305       [2] B. Kagstrom and P. Poromaa; Computing Eigenspaces with Specified
306           Eigenvalues of a Regular Matrix Pair (A, B) and Condition
307           Estimation: Theory, Algorithms and Software,
308           Report UMINF - 94.04, Department of Computing Science, Umea
309           University, S-901 87 Umea, Sweden, 1994. Also as LAPACK Working
310           Note 87. To appear in Numerical Algorithms, 1996.
311
312       [3] B. Kagstrom and P. Poromaa, LAPACK-Style Algorithms and Software
313           for Solving the Generalized Sylvester Equation and Estimating the
314           Separation between Regular Matrix Pairs, Report UMINF - 93.23,
315           Department of Computing Science, Umea University, S-901 87 Umea,
316           Sweden, December 1993, Revised April 1994, Also as LAPACK Working
317           Note 75. To appear in ACM Trans. on Math. Software, Vol 22, No 1,
318           1996.
319
320
321
322
323 LAPACK routine (version 3.1.1)  February 2007                       STGSEN(1)
Impressum