1DTGSEN(1)                LAPACK routine (version 3.2)                DTGSEN(1)
2
3
4

NAME

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

SYNOPSIS

13       SUBROUTINE DTGSEN( 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           DOUBLE         PRECISION PL, PR
22
23           LOGICAL        SELECT( * )
24
25           INTEGER        IWORK( * )
26
27           DOUBLE         PRECISION  A( LDA, * ), ALPHAI( * ), ALPHAR( * ), B(
28                          LDB, * ), BETA( * ), DIF( * ), Q( LDQ, * ), WORK(  *
29                          ), Z( LDZ, * )
30

PURPOSE

32       DTGSEN  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 DGGES), i.e. A is block upper triangular
40       with 1-by-1 and 2-by-2 diagonal blocks. B is upper triangular.
41       DTGSEN also computes the generalized eigenvalues
42                   w(j) = (ALPHAR(j) + i*ALPHAI(j))/BETA(j)
43       of the reordered matrix pair (A, B).
44       Optionally,  DTGSEN computes the estimates of reciprocal condition num‐
45       bers  for  eigenvalues  and  eigenspaces.  These  are   Difu[(A11,B11),
46       (A22,B22)]  and  Difl[(A11,B11),  (A22,B22)],  i.e.  the  separation(s)
47       between the matrix pairs (A11, B11) and (A22,B22)  that  correspond  to
48       the  selected  cluster  and the eigenvalues outside the cluster, resp.,
49       and norms of "projections" onto left and right eigenspaces w.r.t.   the
50       selected cluster in the (1,1)-block.
51

ARGUMENTS

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

FURTHER DETAILS

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