LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dtgsen.f
Go to the documentation of this file.
1*> \brief \b DTGSEN
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download DTGSEN + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dtgsen.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dtgsen.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dtgsen.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE DTGSEN( IJOB, WANTQ, WANTZ, SELECT, N, A, LDA, B, LDB,
20* ALPHAR, ALPHAI, BETA, Q, LDQ, Z, LDZ, M, PL,
21* PR, DIF, WORK, LWORK, IWORK, LIWORK, INFO )
22*
23* .. Scalar Arguments ..
24* LOGICAL WANTQ, WANTZ
25* INTEGER IJOB, INFO, LDA, LDB, LDQ, LDZ, LIWORK, LWORK,
26* $ M, N
27* DOUBLE PRECISION PL, PR
28* ..
29* .. Array Arguments ..
30* LOGICAL SELECT( * )
31* INTEGER IWORK( * )
32* DOUBLE PRECISION A( LDA, * ), ALPHAI( * ), ALPHAR( * ),
33* $ B( LDB, * ), BETA( * ), DIF( * ), Q( LDQ, * ),
34* $ WORK( * ), Z( LDZ, * )
35* ..
36*
37*
38*> \par Purpose:
39* =============
40*>
41*> \verbatim
42*>
43*> DTGSEN reorders the generalized real Schur decomposition of a real
44*> matrix pair (A, B) (in terms of an orthonormal equivalence trans-
45*> formation Q**T * (A, B) * Z), so that a selected cluster of eigenvalues
46*> appears in the leading diagonal blocks of the upper quasi-triangular
47*> matrix A and the upper triangular B. The leading columns of Q and
48*> Z form orthonormal bases of the corresponding left and right eigen-
49*> spaces (deflating subspaces). (A, B) must be in generalized real
50*> Schur canonical form (as returned by DGGES), i.e. A is block upper
51*> triangular with 1-by-1 and 2-by-2 diagonal blocks. B is upper
52*> triangular.
53*>
54*> DTGSEN also computes the generalized eigenvalues
55*>
56*> w(j) = (ALPHAR(j) + i*ALPHAI(j))/BETA(j)
57*>
58*> of the reordered matrix pair (A, B).
59*>
60*> Optionally, DTGSEN computes the estimates of reciprocal condition
61*> numbers for eigenvalues and eigenspaces. These are Difu[(A11,B11),
62*> (A22,B22)] and Difl[(A11,B11), (A22,B22)], i.e. the separation(s)
63*> between the matrix pairs (A11, B11) and (A22,B22) that correspond to
64*> the selected cluster and the eigenvalues outside the cluster, resp.,
65*> and norms of "projections" onto left and right eigenspaces w.r.t.
66*> the selected cluster in the (1,1)-block.
67*> \endverbatim
68*
69* Arguments:
70* ==========
71*
72*> \param[in] IJOB
73*> \verbatim
74*> IJOB is INTEGER
75*> Specifies whether condition numbers are required for the
76*> cluster of eigenvalues (PL and PR) or the deflating subspaces
77*> (Difu and Difl):
78*> =0: Only reorder w.r.t. SELECT. No extras.
79*> =1: Reciprocal of norms of "projections" onto left and right
80*> eigenspaces w.r.t. the selected cluster (PL and PR).
81*> =2: Upper bounds on Difu and Difl. F-norm-based estimate
82*> (DIF(1:2)).
83*> =3: Estimate of Difu and Difl. 1-norm-based estimate
84*> (DIF(1:2)).
85*> About 5 times as expensive as IJOB = 2.
86*> =4: Compute PL, PR and DIF (i.e. 0, 1 and 2 above): Economic
87*> version to get it all.
88*> =5: Compute PL, PR and DIF (i.e. 0, 1 and 3 above)
89*> \endverbatim
90*>
91*> \param[in] WANTQ
92*> \verbatim
93*> WANTQ is LOGICAL
94*> .TRUE. : update the left transformation matrix Q;
95*> .FALSE.: do not update Q.
96*> \endverbatim
97*>
98*> \param[in] WANTZ
99*> \verbatim
100*> WANTZ is LOGICAL
101*> .TRUE. : update the right transformation matrix Z;
102*> .FALSE.: do not update Z.
103*> \endverbatim
104*>
105*> \param[in] SELECT
106*> \verbatim
107*> SELECT is LOGICAL array, dimension (N)
108*> SELECT specifies the eigenvalues in the selected cluster.
109*> To select a real eigenvalue w(j), SELECT(j) must be set to
110*> .TRUE.. To select a complex conjugate pair of eigenvalues
111*> w(j) and w(j+1), corresponding to a 2-by-2 diagonal block,
112*> either SELECT(j) or SELECT(j+1) or both must be set to
113*> .TRUE.; a complex conjugate pair of eigenvalues must be
114*> either both included in the cluster or both excluded.
115*> \endverbatim
116*>
117*> \param[in] N
118*> \verbatim
119*> N is INTEGER
120*> The order of the matrices A and B. N >= 0.
121*> \endverbatim
122*>
123*> \param[in,out] A
124*> \verbatim
125*> A is DOUBLE PRECISION array, dimension(LDA,N)
126*> On entry, the upper quasi-triangular matrix A, with (A, B) in
127*> generalized real Schur canonical form.
128*> On exit, A is overwritten by the reordered matrix A.
129*> \endverbatim
130*>
131*> \param[in] LDA
132*> \verbatim
133*> LDA is INTEGER
134*> The leading dimension of the array A. LDA >= max(1,N).
135*> \endverbatim
136*>
137*> \param[in,out] B
138*> \verbatim
139*> B is DOUBLE PRECISION array, dimension(LDB,N)
140*> On entry, the upper triangular matrix B, with (A, B) in
141*> generalized real Schur canonical form.
142*> On exit, B is overwritten by the reordered matrix B.
143*> \endverbatim
144*>
145*> \param[in] LDB
146*> \verbatim
147*> LDB is INTEGER
148*> The leading dimension of the array B. LDB >= max(1,N).
149*> \endverbatim
150*>
151*> \param[out] ALPHAR
152*> \verbatim
153*> ALPHAR is DOUBLE PRECISION array, dimension (N)
154*> \endverbatim
155*>
156*> \param[out] ALPHAI
157*> \verbatim
158*> ALPHAI is DOUBLE PRECISION array, dimension (N)
159*> \endverbatim
160*>
161*> \param[out] BETA
162*> \verbatim
163*> BETA is DOUBLE PRECISION array, dimension (N)
164*>
165*> On exit, (ALPHAR(j) + ALPHAI(j)*i)/BETA(j), j=1,...,N, will
166*> be the generalized eigenvalues. ALPHAR(j) + ALPHAI(j)*i
167*> and BETA(j),j=1,...,N are the diagonals of the complex Schur
168*> form (S,T) that would result if the 2-by-2 diagonal blocks of
169*> the real generalized Schur form of (A,B) were further reduced
170*> to triangular form using complex unitary transformations.
171*> If ALPHAI(j) is zero, then the j-th eigenvalue is real; if
172*> positive, then the j-th and (j+1)-st eigenvalues are a
173*> complex conjugate pair, with ALPHAI(j+1) negative.
174*> \endverbatim
175*>
176*> \param[in,out] Q
177*> \verbatim
178*> Q is DOUBLE PRECISION array, dimension (LDQ,N)
179*> On entry, if WANTQ = .TRUE., Q is an N-by-N matrix.
180*> On exit, Q has been postmultiplied by the left orthogonal
181*> transformation matrix which reorder (A, B); The leading M
182*> columns of Q form orthonormal bases for the specified pair of
183*> left eigenspaces (deflating subspaces).
184*> If WANTQ = .FALSE., Q is not referenced.
185*> \endverbatim
186*>
187*> \param[in] LDQ
188*> \verbatim
189*> LDQ is INTEGER
190*> The leading dimension of the array Q. LDQ >= 1;
191*> and if WANTQ = .TRUE., LDQ >= N.
192*> \endverbatim
193*>
194*> \param[in,out] Z
195*> \verbatim
196*> Z is DOUBLE PRECISION array, dimension (LDZ,N)
197*> On entry, if WANTZ = .TRUE., Z is an N-by-N matrix.
198*> On exit, Z has been postmultiplied by the left orthogonal
199*> transformation matrix which reorder (A, B); The leading M
200*> columns of Z form orthonormal bases for the specified pair of
201*> left eigenspaces (deflating subspaces).
202*> If WANTZ = .FALSE., Z is not referenced.
203*> \endverbatim
204*>
205*> \param[in] LDZ
206*> \verbatim
207*> LDZ is INTEGER
208*> The leading dimension of the array Z. LDZ >= 1;
209*> If WANTZ = .TRUE., LDZ >= N.
210*> \endverbatim
211*>
212*> \param[out] M
213*> \verbatim
214*> M is INTEGER
215*> The dimension of the specified pair of left and right eigen-
216*> spaces (deflating subspaces). 0 <= M <= N.
217*> \endverbatim
218*>
219*> \param[out] PL
220*> \verbatim
221*> PL is DOUBLE PRECISION
222*> \endverbatim
223*>
224*> \param[out] PR
225*> \verbatim
226*> PR is DOUBLE PRECISION
227*>
228*> If IJOB = 1, 4 or 5, PL, PR are lower bounds on the
229*> reciprocal of the norm of "projections" onto left and right
230*> eigenspaces with respect to the selected cluster.
231*> 0 < PL, PR <= 1.
232*> If M = 0 or M = N, PL = PR = 1.
233*> If IJOB = 0, 2 or 3, PL and PR are not referenced.
234*> \endverbatim
235*>
236*> \param[out] DIF
237*> \verbatim
238*> DIF is DOUBLE PRECISION array, dimension (2).
239*> If IJOB >= 2, DIF(1:2) store the estimates of Difu and Difl.
240*> If IJOB = 2 or 4, DIF(1:2) are F-norm-based upper bounds on
241*> Difu and Difl. If IJOB = 3 or 5, DIF(1:2) are 1-norm-based
242*> estimates of Difu and Difl.
243*> If M = 0 or N, DIF(1:2) = F-norm([A, B]).
244*> If IJOB = 0 or 1, DIF is not referenced.
245*> \endverbatim
246*>
247*> \param[out] WORK
248*> \verbatim
249*> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
250*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
251*> \endverbatim
252*>
253*> \param[in] LWORK
254*> \verbatim
255*> LWORK is INTEGER
256*> The dimension of the array WORK. LWORK >= 4*N+16.
257*> If IJOB = 1, 2 or 4, LWORK >= MAX(4*N+16, 2*M*(N-M)).
258*> If IJOB = 3 or 5, LWORK >= MAX(4*N+16, 4*M*(N-M)).
259*>
260*> If LWORK = -1, then a workspace query is assumed; the routine
261*> only calculates the optimal size of the WORK array, returns
262*> this value as the first entry of the WORK array, and no error
263*> message related to LWORK is issued by XERBLA.
264*> \endverbatim
265*>
266*> \param[out] IWORK
267*> \verbatim
268*> IWORK is INTEGER array, dimension (MAX(1,LIWORK))
269*> On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK.
270*> \endverbatim
271*>
272*> \param[in] LIWORK
273*> \verbatim
274*> LIWORK is INTEGER
275*> The dimension of the array IWORK. LIWORK >= 1.
276*> If IJOB = 1, 2 or 4, LIWORK >= N+6.
277*> If IJOB = 3 or 5, LIWORK >= MAX(2*M*(N-M), N+6).
278*>
279*> If LIWORK = -1, then a workspace query is assumed; the
280*> routine only calculates the optimal size of the IWORK array,
281*> returns this value as the first entry of the IWORK array, and
282*> no error message related to LIWORK is issued by XERBLA.
283*> \endverbatim
284*>
285*> \param[out] INFO
286*> \verbatim
287*> INFO is INTEGER
288*> =0: Successful exit.
289*> <0: If INFO = -i, the i-th argument had an illegal value.
290*> =1: Reordering of (A, B) failed because the transformed
291*> matrix pair (A, B) would be too far from generalized
292*> Schur form; the problem is very ill-conditioned.
293*> (A, B) may have been partially reordered.
294*> If requested, 0 is returned in DIF(*), PL and PR.
295*> \endverbatim
296*
297* Authors:
298* ========
299*
300*> \author Univ. of Tennessee
301*> \author Univ. of California Berkeley
302*> \author Univ. of Colorado Denver
303*> \author NAG Ltd.
304*
305*> \ingroup tgsen
306*
307*> \par Further Details:
308* =====================
309*>
310*> \verbatim
311*>
312*> DTGSEN first collects the selected eigenvalues by computing
313*> orthogonal U and W that move them to the top left corner of (A, B).
314*> In other words, the selected eigenvalues are the eigenvalues of
315*> (A11, B11) in:
316*>
317*> U**T*(A, B)*W = (A11 A12) (B11 B12) n1
318*> ( 0 A22),( 0 B22) n2
319*> n1 n2 n1 n2
320*>
321*> where N = n1+n2 and U**T means the transpose of U. The first n1 columns
322*> of U and W span the specified pair of left and right eigenspaces
323*> (deflating subspaces) of (A, B).
324*>
325*> If (A, B) has been obtained from the generalized real Schur
326*> decomposition of a matrix pair (C, D) = Q*(A, B)*Z**T, then the
327*> reordered generalized real Schur form of (C, D) is given by
328*>
329*> (C, D) = (Q*U)*(U**T*(A, B)*W)*(Z*W)**T,
330*>
331*> and the first n1 columns of Q*U and Z*W span the corresponding
332*> deflating subspaces of (C, D) (Q and Z store Q*U and Z*W, resp.).
333*>
334*> Note that if the selected eigenvalue is sufficiently ill-conditioned,
335*> then its value may differ significantly from its value before
336*> reordering.
337*>
338*> The reciprocal condition numbers of the left and right eigenspaces
339*> spanned by the first n1 columns of U and W (or Q*U and Z*W) may
340*> be returned in DIF(1:2), corresponding to Difu and Difl, resp.
341*>
342*> The Difu and Difl are defined as:
343*>
344*> Difu[(A11, B11), (A22, B22)] = sigma-min( Zu )
345*> and
346*> Difl[(A11, B11), (A22, B22)] = Difu[(A22, B22), (A11, B11)],
347*>
348*> where sigma-min(Zu) is the smallest singular value of the
349*> (2*n1*n2)-by-(2*n1*n2) matrix
350*>
351*> Zu = [ kron(In2, A11) -kron(A22**T, In1) ]
352*> [ kron(In2, B11) -kron(B22**T, In1) ].
353*>
354*> Here, Inx is the identity matrix of size nx and A22**T is the
355*> transpose of A22. kron(X, Y) is the Kronecker product between
356*> the matrices X and Y.
357*>
358*> When DIF(2) is small, small changes in (A, B) can cause large changes
359*> in the deflating subspace. An approximate (asymptotic) bound on the
360*> maximum angular error in the computed deflating subspaces is
361*>
362*> EPS * norm((A, B)) / DIF(2),
363*>
364*> where EPS is the machine precision.
365*>
366*> The reciprocal norm of the projectors on the left and right
367*> eigenspaces associated with (A11, B11) may be returned in PL and PR.
368*> They are computed as follows. First we compute L and R so that
369*> P*(A, B)*Q is block diagonal, where
370*>
371*> P = ( I -L ) n1 Q = ( I R ) n1
372*> ( 0 I ) n2 and ( 0 I ) n2
373*> n1 n2 n1 n2
374*>
375*> and (L, R) is the solution to the generalized Sylvester equation
376*>
377*> A11*R - L*A22 = -A12
378*> B11*R - L*B22 = -B12
379*>
380*> Then PL = (F-norm(L)**2+1)**(-1/2) and PR = (F-norm(R)**2+1)**(-1/2).
381*> An approximate (asymptotic) bound on the average absolute error of
382*> the selected eigenvalues is
383*>
384*> EPS * norm((A, B)) / PL.
385*>
386*> There are also global error bounds which valid for perturbations up
387*> to a certain restriction: A lower bound (x) on the smallest
388*> F-norm(E,F) for which an eigenvalue of (A11, B11) may move and
389*> coalesce with an eigenvalue of (A22, B22) under perturbation (E,F),
390*> (i.e. (A + E, B + F), is
391*>
392*> x = min(Difu,Difl)/((1/(PL*PL)+1/(PR*PR))**(1/2)+2*max(1/PL,1/PR)).
393*>
394*> An approximate bound on x can be computed from DIF(1:2), PL and PR.
395*>
396*> If y = ( F-norm(E,F) / x) <= 1, the angles between the perturbed
397*> (L', R') and unperturbed (L, R) left and right deflating subspaces
398*> associated with the selected cluster in the (1,1)-blocks can be
399*> bounded as
400*>
401*> max-angle(L, L') <= arctan( y * PL / (1 - y * (1 - PL * PL)**(1/2))
402*> max-angle(R, R') <= arctan( y * PR / (1 - y * (1 - PR * PR)**(1/2))
403*>
404*> See LAPACK User's Guide section 4.11 or the following references
405*> for more information.
406*>
407*> Note that if the default method for computing the Frobenius-norm-
408*> based estimate DIF is not wanted (see DLATDF), then the parameter
409*> IDIFJB (see below) should be changed from 3 to 4 (routine DLATDF
410*> (IJOB = 2 will be used)). See DTGSYL for more details.
411*> \endverbatim
412*
413*> \par Contributors:
414* ==================
415*>
416*> Bo Kagstrom and Peter Poromaa, Department of Computing Science,
417*> Umea University, S-901 87 Umea, Sweden.
418*
419*> \par References:
420* ================
421*>
422*> \verbatim
423*>
424*> [1] B. Kagstrom; A Direct Method for Reordering Eigenvalues in the
425*> Generalized Real Schur Form of a Regular Matrix Pair (A, B), in
426*> M.S. Moonen et al (eds), Linear Algebra for Large Scale and
427*> Real-Time Applications, Kluwer Academic Publ. 1993, pp 195-218.
428*>
429*> [2] B. Kagstrom and P. Poromaa; Computing Eigenspaces with Specified
430*> Eigenvalues of a Regular Matrix Pair (A, B) and Condition
431*> Estimation: Theory, Algorithms and Software,
432*> Report UMINF - 94.04, Department of Computing Science, Umea
433*> University, S-901 87 Umea, Sweden, 1994. Also as LAPACK Working
434*> Note 87. To appear in Numerical Algorithms, 1996.
435*>
436*> [3] B. Kagstrom and P. Poromaa, LAPACK-Style Algorithms and Software
437*> for Solving the Generalized Sylvester Equation and Estimating the
438*> Separation between Regular Matrix Pairs, Report UMINF - 93.23,
439*> Department of Computing Science, Umea University, S-901 87 Umea,
440*> Sweden, December 1993, Revised April 1994, Also as LAPACK Working
441*> Note 75. To appear in ACM Trans. on Math. Software, Vol 22, No 1,
442*> 1996.
443*> \endverbatim
444*>
445* =====================================================================
446 SUBROUTINE dtgsen( IJOB, WANTQ, WANTZ, SELECT, N, A, LDA, B,
447 $ LDB,
448 $ ALPHAR, ALPHAI, BETA, Q, LDQ, Z, LDZ, M, PL,
449 $ PR, DIF, WORK, LWORK, IWORK, LIWORK, INFO )
450*
451* -- LAPACK computational routine --
452* -- LAPACK is a software package provided by Univ. of Tennessee, --
453* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
454*
455* .. Scalar Arguments ..
456 LOGICAL WANTQ, WANTZ
457 INTEGER IJOB, INFO, LDA, LDB, LDQ, LDZ, LIWORK, LWORK,
458 $ M, N
459 DOUBLE PRECISION PL, PR
460* ..
461* .. Array Arguments ..
462 LOGICAL SELECT( * )
463 INTEGER IWORK( * )
464 DOUBLE PRECISION A( LDA, * ), ALPHAI( * ), ALPHAR( * ),
465 $ b( ldb, * ), beta( * ), dif( * ), q( ldq, * ),
466 $ work( * ), z( ldz, * )
467* ..
468*
469* =====================================================================
470*
471* .. Parameters ..
472 INTEGER IDIFJB
473 PARAMETER ( IDIFJB = 3 )
474 double precision zero, one
475 parameter( zero = 0.0d+0, one = 1.0d+0 )
476* ..
477* .. Local Scalars ..
478 LOGICAL LQUERY, PAIR, SWAP, WANTD, WANTD1, WANTD2,
479 $ WANTP
480 INTEGER I, IERR, IJB, K, KASE, KK, KS, LIWMIN, LWMIN,
481 $ mn2, n1, n2
482 DOUBLE PRECISION DSCALE, DSUM, EPS, RDSCAL, SMLNUM
483* ..
484* .. Local Arrays ..
485 INTEGER ISAVE( 3 )
486* ..
487* .. External Subroutines ..
488 EXTERNAL dlacn2, dlacpy, dlag2, dlassq, dtgexc,
489 $ dtgsyl,
490 $ xerbla
491* ..
492* .. External Functions ..
493 DOUBLE PRECISION DLAMCH
494 EXTERNAL DLAMCH
495* ..
496* .. Intrinsic Functions ..
497 INTRINSIC max, sign, sqrt
498* ..
499* .. Executable Statements ..
500*
501* Decode and test the input parameters
502*
503 info = 0
504 lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
505*
506 IF( ijob.LT.0 .OR. ijob.GT.5 ) THEN
507 info = -1
508 ELSE IF( n.LT.0 ) THEN
509 info = -5
510 ELSE IF( lda.LT.max( 1, n ) ) THEN
511 info = -7
512 ELSE IF( ldb.LT.max( 1, n ) ) THEN
513 info = -9
514 ELSE IF( ldq.LT.1 .OR. ( wantq .AND. ldq.LT.n ) ) THEN
515 info = -14
516 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
517 info = -16
518 END IF
519*
520 IF( info.NE.0 ) THEN
521 CALL xerbla( 'DTGSEN', -info )
522 RETURN
523 END IF
524*
525* Get machine constants
526*
527 eps = dlamch( 'P' )
528 smlnum = dlamch( 'S' ) / eps
529 ierr = 0
530*
531 wantp = ijob.EQ.1 .OR. ijob.GE.4
532 wantd1 = ijob.EQ.2 .OR. ijob.EQ.4
533 wantd2 = ijob.EQ.3 .OR. ijob.EQ.5
534 wantd = wantd1 .OR. wantd2
535*
536* Set M to the dimension of the specified pair of deflating
537* subspaces.
538*
539 m = 0
540 pair = .false.
541 IF( .NOT.lquery .OR. ijob.NE.0 ) THEN
542 DO 10 k = 1, n
543 IF( pair ) THEN
544 pair = .false.
545 ELSE
546 IF( k.LT.n ) THEN
547 IF( a( k+1, k ).EQ.zero ) THEN
548 IF( SELECT( k ) )
549 $ m = m + 1
550 ELSE
551 pair = .true.
552 IF( SELECT( k ) .OR. SELECT( k+1 ) )
553 $ m = m + 2
554 END IF
555 ELSE
556 IF( SELECT( n ) )
557 $ m = m + 1
558 END IF
559 END IF
560 10 CONTINUE
561 END IF
562*
563 IF( ijob.EQ.1 .OR. ijob.EQ.2 .OR. ijob.EQ.4 ) THEN
564 lwmin = max( 1, 4*n+16, 2*m*( n-m ) )
565 liwmin = max( 1, n+6 )
566 ELSE IF( ijob.EQ.3 .OR. ijob.EQ.5 ) THEN
567 lwmin = max( 1, 4*n+16, 4*m*( n-m ) )
568 liwmin = max( 1, 2*m*( n-m ), n+6 )
569 ELSE
570 lwmin = max( 1, 4*n+16 )
571 liwmin = 1
572 END IF
573*
574 work( 1 ) = lwmin
575 iwork( 1 ) = liwmin
576*
577 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
578 info = -22
579 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
580 info = -24
581 END IF
582*
583 IF( info.NE.0 ) THEN
584 CALL xerbla( 'DTGSEN', -info )
585 RETURN
586 ELSE IF( lquery ) THEN
587 RETURN
588 END IF
589*
590* Quick return if possible.
591*
592 IF( m.EQ.n .OR. m.EQ.0 ) THEN
593 IF( wantp ) THEN
594 pl = one
595 pr = one
596 END IF
597 IF( wantd ) THEN
598 dscale = zero
599 dsum = one
600 DO 20 i = 1, n
601 CALL dlassq( n, a( 1, i ), 1, dscale, dsum )
602 CALL dlassq( n, b( 1, i ), 1, dscale, dsum )
603 20 CONTINUE
604 dif( 1 ) = dscale*sqrt( dsum )
605 dif( 2 ) = dif( 1 )
606 END IF
607 GO TO 60
608 END IF
609*
610* Collect the selected blocks at the top-left corner of (A, B).
611*
612 ks = 0
613 pair = .false.
614 DO 30 k = 1, n
615 IF( pair ) THEN
616 pair = .false.
617 ELSE
618*
619 swap = SELECT( k )
620 IF( k.LT.n ) THEN
621 IF( a( k+1, k ).NE.zero ) THEN
622 pair = .true.
623 swap = swap .OR. SELECT( k+1 )
624 END IF
625 END IF
626*
627 IF( swap ) THEN
628 ks = ks + 1
629*
630* Swap the K-th block to position KS.
631* Perform the reordering of diagonal blocks in (A, B)
632* by orthogonal transformation matrices and update
633* Q and Z accordingly (if requested):
634*
635 kk = k
636 IF( k.NE.ks )
637 $ CALL dtgexc( wantq, wantz, n, a, lda, b, ldb, q,
638 $ ldq,
639 $ z, ldz, kk, ks, work, lwork, ierr )
640*
641 IF( ierr.GT.0 ) THEN
642*
643* Swap is rejected: exit.
644*
645 info = 1
646 IF( wantp ) THEN
647 pl = zero
648 pr = zero
649 END IF
650 IF( wantd ) THEN
651 dif( 1 ) = zero
652 dif( 2 ) = zero
653 END IF
654 GO TO 60
655 END IF
656*
657 IF( pair )
658 $ ks = ks + 1
659 END IF
660 END IF
661 30 CONTINUE
662 IF( wantp ) THEN
663*
664* Solve generalized Sylvester equation for R and L
665* and compute PL and PR.
666*
667 n1 = m
668 n2 = n - m
669 i = n1 + 1
670 ijb = 0
671 CALL dlacpy( 'Full', n1, n2, a( 1, i ), lda, work, n1 )
672 CALL dlacpy( 'Full', n1, n2, b( 1, i ), ldb,
673 $ work( n1*n2+1 ),
674 $ n1 )
675 CALL dtgsyl( 'N', ijb, n1, n2, a, lda, a( i, i ), lda, work,
676 $ n1, b, ldb, b( i, i ), ldb, work( n1*n2+1 ), n1,
677 $ dscale, dif( 1 ), work( n1*n2*2+1 ),
678 $ lwork-2*n1*n2, iwork, ierr )
679*
680* Estimate the reciprocal of norms of "projections" onto left
681* and right eigenspaces.
682*
683 rdscal = zero
684 dsum = one
685 CALL dlassq( n1*n2, work, 1, rdscal, dsum )
686 pl = rdscal*sqrt( dsum )
687 IF( pl.EQ.zero ) THEN
688 pl = one
689 ELSE
690 pl = dscale / ( sqrt( dscale*dscale / pl+pl )*sqrt( pl ) )
691 END IF
692 rdscal = zero
693 dsum = one
694 CALL dlassq( n1*n2, work( n1*n2+1 ), 1, rdscal, dsum )
695 pr = rdscal*sqrt( dsum )
696 IF( pr.EQ.zero ) THEN
697 pr = one
698 ELSE
699 pr = dscale / ( sqrt( dscale*dscale / pr+pr )*sqrt( pr ) )
700 END IF
701 END IF
702*
703 IF( wantd ) THEN
704*
705* Compute estimates of Difu and Difl.
706*
707 IF( wantd1 ) THEN
708 n1 = m
709 n2 = n - m
710 i = n1 + 1
711 ijb = idifjb
712*
713* Frobenius norm-based Difu-estimate.
714*
715 CALL dtgsyl( 'N', ijb, n1, n2, a, lda, a( i, i ), lda,
716 $ work,
717 $ n1, b, ldb, b( i, i ), ldb, work( n1*n2+1 ),
718 $ n1, dscale, dif( 1 ), work( 2*n1*n2+1 ),
719 $ lwork-2*n1*n2, iwork, ierr )
720*
721* Frobenius norm-based Difl-estimate.
722*
723 CALL dtgsyl( 'N', ijb, n2, n1, a( i, i ), lda, a, lda,
724 $ work,
725 $ n2, b( i, i ), ldb, b, ldb, work( n1*n2+1 ),
726 $ n2, dscale, dif( 2 ), work( 2*n1*n2+1 ),
727 $ lwork-2*n1*n2, iwork, ierr )
728 ELSE
729*
730*
731* Compute 1-norm-based estimates of Difu and Difl using
732* reversed communication with DLACN2. In each step a
733* generalized Sylvester equation or a transposed variant
734* is solved.
735*
736 kase = 0
737 n1 = m
738 n2 = n - m
739 i = n1 + 1
740 ijb = 0
741 mn2 = 2*n1*n2
742*
743* 1-norm-based estimate of Difu.
744*
745 40 CONTINUE
746 CALL dlacn2( mn2, work( mn2+1 ), work, iwork, dif( 1 ),
747 $ kase, isave )
748 IF( kase.NE.0 ) THEN
749 IF( kase.EQ.1 ) THEN
750*
751* Solve generalized Sylvester equation.
752*
753 CALL dtgsyl( 'N', ijb, n1, n2, a, lda, a( i, i ),
754 $ lda,
755 $ work, n1, b, ldb, b( i, i ), ldb,
756 $ work( n1*n2+1 ), n1, dscale, dif( 1 ),
757 $ work( 2*n1*n2+1 ), lwork-2*n1*n2, iwork,
758 $ ierr )
759 ELSE
760*
761* Solve the transposed variant.
762*
763 CALL dtgsyl( 'T', ijb, n1, n2, a, lda, a( i, i ),
764 $ lda,
765 $ work, n1, b, ldb, b( i, i ), ldb,
766 $ work( n1*n2+1 ), n1, dscale, dif( 1 ),
767 $ work( 2*n1*n2+1 ), lwork-2*n1*n2, iwork,
768 $ ierr )
769 END IF
770 GO TO 40
771 END IF
772 dif( 1 ) = dscale / dif( 1 )
773*
774* 1-norm-based estimate of Difl.
775*
776 50 CONTINUE
777 CALL dlacn2( mn2, work( mn2+1 ), work, iwork, dif( 2 ),
778 $ kase, isave )
779 IF( kase.NE.0 ) THEN
780 IF( kase.EQ.1 ) THEN
781*
782* Solve generalized Sylvester equation.
783*
784 CALL dtgsyl( 'N', ijb, n2, n1, a( i, i ), lda, a,
785 $ lda,
786 $ work, n2, b( i, i ), ldb, b, ldb,
787 $ work( n1*n2+1 ), n2, dscale, dif( 2 ),
788 $ work( 2*n1*n2+1 ), lwork-2*n1*n2, iwork,
789 $ ierr )
790 ELSE
791*
792* Solve the transposed variant.
793*
794 CALL dtgsyl( 'T', ijb, n2, n1, a( i, i ), lda, a,
795 $ lda,
796 $ work, n2, b( i, i ), ldb, b, ldb,
797 $ work( n1*n2+1 ), n2, dscale, dif( 2 ),
798 $ work( 2*n1*n2+1 ), lwork-2*n1*n2, iwork,
799 $ ierr )
800 END IF
801 GO TO 50
802 END IF
803 dif( 2 ) = dscale / dif( 2 )
804*
805 END IF
806 END IF
807*
808 60 CONTINUE
809*
810* Compute generalized eigenvalues of reordered pair (A, B) and
811* normalize the generalized Schur form.
812*
813 pair = .false.
814 DO 80 k = 1, n
815 IF( pair ) THEN
816 pair = .false.
817 ELSE
818*
819 IF( k.LT.n ) THEN
820 IF( a( k+1, k ).NE.zero ) THEN
821 pair = .true.
822 END IF
823 END IF
824*
825 IF( pair ) THEN
826*
827* Compute the eigenvalue(s) at position K.
828*
829 work( 1 ) = a( k, k )
830 work( 2 ) = a( k+1, k )
831 work( 3 ) = a( k, k+1 )
832 work( 4 ) = a( k+1, k+1 )
833 work( 5 ) = b( k, k )
834 work( 6 ) = b( k+1, k )
835 work( 7 ) = b( k, k+1 )
836 work( 8 ) = b( k+1, k+1 )
837 CALL dlag2( work, 2, work( 5 ), 2, smlnum*eps,
838 $ beta( k ),
839 $ beta( k+1 ), alphar( k ), alphar( k+1 ),
840 $ alphai( k ) )
841 alphai( k+1 ) = -alphai( k )
842*
843 ELSE
844*
845 IF( sign( one, b( k, k ) ).LT.zero ) THEN
846*
847* If B(K,K) is negative, make it positive
848*
849 DO 70 i = 1, n
850 a( k, i ) = -a( k, i )
851 b( k, i ) = -b( k, i )
852 IF( wantq ) q( i, k ) = -q( i, k )
853 70 CONTINUE
854 END IF
855*
856 alphar( k ) = a( k, k )
857 alphai( k ) = zero
858 beta( k ) = b( k, k )
859*
860 END IF
861 END IF
862 80 CONTINUE
863*
864 work( 1 ) = lwmin
865 iwork( 1 ) = liwmin
866*
867 RETURN
868*
869* End of DTGSEN
870*
871 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dlacn2(n, v, x, isgn, est, kase, isave)
DLACN2 estimates the 1-norm of a square matrix, using reverse communication for evaluating matrix-vec...
Definition dlacn2.f:134
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:101
subroutine dlag2(a, lda, b, ldb, safmin, scale1, scale2, wr1, wr2, wi)
DLAG2 computes the eigenvalues of a 2-by-2 generalized eigenvalue problem, with scaling as necessary ...
Definition dlag2.f:154
subroutine dlassq(n, x, incx, scale, sumsq)
DLASSQ updates a sum of squares represented in scaled form.
Definition dlassq.f90:122
subroutine dtgexc(wantq, wantz, n, a, lda, b, ldb, q, ldq, z, ldz, ifst, ilst, work, lwork, info)
DTGEXC
Definition dtgexc.f:218
subroutine dtgsen(ijob, wantq, wantz, select, n, a, lda, b, ldb, alphar, alphai, beta, q, ldq, z, ldz, m, pl, pr, dif, work, lwork, iwork, liwork, info)
DTGSEN
Definition dtgsen.f:450
subroutine dtgsyl(trans, ijob, m, n, a, lda, b, ldb, c, ldc, d, ldd, e, lde, f, ldf, scale, dif, work, lwork, iwork, info)
DTGSYL
Definition dtgsyl.f:298