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
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
|
*> \brief \b DTGEX2 swaps adjacent diagonal blocks in an upper (quasi) triangular matrix pair by an orthogonal equivalence transformation.
*
* =========== DOCUMENTATION ===========
*
* Online html documentation available at
* http://www.netlib.org/lapack/explore-html/
*
*> \htmlonly
*> Download DTGEX2 + dependencies
*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dtgex2.f">
*> [TGZ]</a>
*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dtgex2.f">
*> [ZIP]</a>
*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dtgex2.f">
*> [TXT]</a>
*> \endhtmlonly
*
* Definition:
* ===========
*
* SUBROUTINE DTGEX2( WANTQ, WANTZ, N, A, LDA, B, LDB, Q, LDQ, Z,
* LDZ, J1, N1, N2, WORK, LWORK, INFO )
*
* .. Scalar Arguments ..
* LOGICAL WANTQ, WANTZ
* INTEGER INFO, J1, LDA, LDB, LDQ, LDZ, LWORK, N, N1, N2
* ..
* .. Array Arguments ..
* DOUBLE PRECISION A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
* $ WORK( * ), Z( LDZ, * )
* ..
*
*
*> \par Purpose:
* =============
*>
*> \verbatim
*>
*> DTGEX2 swaps adjacent diagonal blocks (A11, B11) and (A22, B22)
*> of size 1-by-1 or 2-by-2 in an upper (quasi) triangular matrix pair
*> (A, B) by an orthogonal equivalence transformation.
*>
*> (A, B) must be in generalized real Schur canonical form (as returned
*> by DGGES), i.e. A is block upper triangular with 1-by-1 and 2-by-2
*> diagonal blocks. B is upper triangular.
*>
*> Optionally, the matrices Q and Z of generalized Schur vectors are
*> updated.
*>
*> Q(in) * A(in) * Z(in)**T = Q(out) * A(out) * Z(out)**T
*> Q(in) * B(in) * Z(in)**T = Q(out) * B(out) * Z(out)**T
*>
*> \endverbatim
*
* Arguments:
* ==========
*
*> \param[in] WANTQ
*> \verbatim
*> WANTQ is LOGICAL
*> .TRUE. : update the left transformation matrix Q;
*> .FALSE.: do not update Q.
*> \endverbatim
*>
*> \param[in] WANTZ
*> \verbatim
*> WANTZ is LOGICAL
*> .TRUE. : update the right transformation matrix Z;
*> .FALSE.: do not update Z.
*> \endverbatim
*>
*> \param[in] N
*> \verbatim
*> N is INTEGER
*> The order of the matrices A and B. N >= 0.
*> \endverbatim
*>
*> \param[in,out] A
*> \verbatim
*> A is DOUBLE PRECISION array, dimensions (LDA,N)
*> On entry, the matrix A in the pair (A, B).
*> On exit, the updated matrix A.
*> \endverbatim
*>
*> \param[in] LDA
*> \verbatim
*> LDA is INTEGER
*> The leading dimension of the array A. LDA >= max(1,N).
*> \endverbatim
*>
*> \param[in,out] B
*> \verbatim
*> B is DOUBLE PRECISION array, dimensions (LDB,N)
*> On entry, the matrix B in the pair (A, B).
*> On exit, the updated matrix B.
*> \endverbatim
*>
*> \param[in] LDB
*> \verbatim
*> LDB is INTEGER
*> The leading dimension of the array B. LDB >= max(1,N).
*> \endverbatim
*>
*> \param[in,out] Q
*> \verbatim
*> Q is DOUBLE PRECISION array, dimension (LDQ,N)
*> On entry, if WANTQ = .TRUE., the orthogonal matrix Q.
*> On exit, the updated matrix Q.
*> Not referenced if WANTQ = .FALSE..
*> \endverbatim
*>
*> \param[in] LDQ
*> \verbatim
*> LDQ is INTEGER
*> The leading dimension of the array Q. LDQ >= 1.
*> If WANTQ = .TRUE., LDQ >= N.
*> \endverbatim
*>
*> \param[in,out] Z
*> \verbatim
*> Z is DOUBLE PRECISION array, dimension (LDZ,N)
*> On entry, if WANTZ =.TRUE., the orthogonal matrix Z.
*> On exit, the updated matrix Z.
*> Not referenced if WANTZ = .FALSE..
*> \endverbatim
*>
*> \param[in] LDZ
*> \verbatim
*> LDZ is INTEGER
*> The leading dimension of the array Z. LDZ >= 1.
*> If WANTZ = .TRUE., LDZ >= N.
*> \endverbatim
*>
*> \param[in] J1
*> \verbatim
*> J1 is INTEGER
*> The index to the first block (A11, B11). 1 <= J1 <= N.
*> \endverbatim
*>
*> \param[in] N1
*> \verbatim
*> N1 is INTEGER
*> The order of the first block (A11, B11). N1 = 0, 1 or 2.
*> \endverbatim
*>
*> \param[in] N2
*> \verbatim
*> N2 is INTEGER
*> The order of the second block (A22, B22). N2 = 0, 1 or 2.
*> \endverbatim
*>
*> \param[out] WORK
*> \verbatim
*> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK)).
*> \endverbatim
*>
*> \param[in] LWORK
*> \verbatim
*> LWORK is INTEGER
*> The dimension of the array WORK.
*> LWORK >= MAX( 1, N*(N2+N1), (N2+N1)*(N2+N1)*2 )
*> \endverbatim
*>
*> \param[out] INFO
*> \verbatim
*> INFO is INTEGER
*> =0: Successful exit
*> >0: If INFO = 1, the transformed matrix (A, B) would be
*> too far from generalized Schur form; the blocks are
*> not swapped and (A, B) and (Q, Z) are unchanged.
*> The problem of swapping is too ill-conditioned.
*> <0: If INFO = -16: LWORK is too small. Appropriate value
*> for LWORK is returned in WORK(1).
*> \endverbatim
*
* Authors:
* ========
*
*> \author Univ. of Tennessee
*> \author Univ. of California Berkeley
*> \author Univ. of Colorado Denver
*> \author NAG Ltd.
*
*> \date September 2012
*
*> \ingroup doubleGEauxiliary
*
*> \par Further Details:
* =====================
*>
*> In the current code both weak and strong stability tests are
*> performed. The user can omit the strong stability test by changing
*> the internal logical parameter WANDS to .FALSE.. See ref. [2] for
*> details.
*
*> \par Contributors:
* ==================
*>
*> Bo Kagstrom and Peter Poromaa, Department of Computing Science,
*> Umea University, S-901 87 Umea, Sweden.
*
*> \par References:
* ================
*>
*> \verbatim
*>
*> [1] B. Kagstrom; A Direct Method for Reordering Eigenvalues in the
*> Generalized Real Schur Form of a Regular Matrix Pair (A, B), in
*> M.S. Moonen et al (eds), Linear Algebra for Large Scale and
*> Real-Time Applications, Kluwer Academic Publ. 1993, pp 195-218.
*>
*> [2] B. Kagstrom and P. Poromaa; Computing Eigenspaces with Specified
*> Eigenvalues of a Regular Matrix Pair (A, B) and Condition
*> Estimation: Theory, Algorithms and Software,
*> Report UMINF - 94.04, Department of Computing Science, Umea
*> University, S-901 87 Umea, Sweden, 1994. Also as LAPACK Working
*> Note 87. To appear in Numerical Algorithms, 1996.
*> \endverbatim
*>
* =====================================================================
SUBROUTINE DTGEX2( WANTQ, WANTZ, N, A, LDA, B, LDB, Q, LDQ, Z,
$ LDZ, J1, N1, N2, WORK, LWORK, INFO )
*
* -- LAPACK auxiliary routine (version 3.4.2) --
* -- LAPACK is a software package provided by Univ. of Tennessee, --
* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
* September 2012
*
* .. Scalar Arguments ..
LOGICAL WANTQ, WANTZ
INTEGER INFO, J1, LDA, LDB, LDQ, LDZ, LWORK, N, N1, N2
* ..
* .. Array Arguments ..
DOUBLE PRECISION A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
$ WORK( * ), Z( LDZ, * )
* ..
*
* =====================================================================
* Replaced various illegal calls to DCOPY by calls to DLASET, or by DO
* loops. Sven Hammarling, 1/5/02.
*
* .. Parameters ..
DOUBLE PRECISION ZERO, ONE
PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 )
DOUBLE PRECISION TWENTY
PARAMETER ( TWENTY = 2.0D+01 )
INTEGER LDST
PARAMETER ( LDST = 4 )
LOGICAL WANDS
PARAMETER ( WANDS = .TRUE. )
* ..
* .. Local Scalars ..
LOGICAL DTRONG, WEAK
INTEGER I, IDUM, LINFO, M
DOUBLE PRECISION BQRA21, BRQA21, DDUM, DNORM, DSCALE, DSUM, EPS,
$ F, G, SA, SB, SCALE, SMLNUM, SS, THRESH, WS
* ..
* .. Local Arrays ..
INTEGER IWORK( LDST )
DOUBLE PRECISION AI( 2 ), AR( 2 ), BE( 2 ), IR( LDST, LDST ),
$ IRCOP( LDST, LDST ), LI( LDST, LDST ),
$ LICOP( LDST, LDST ), S( LDST, LDST ),
$ SCPY( LDST, LDST ), T( LDST, LDST ),
$ TAUL( LDST ), TAUR( LDST ), TCPY( LDST, LDST )
* ..
* .. External Functions ..
DOUBLE PRECISION DLAMCH
EXTERNAL DLAMCH
* ..
* .. External Subroutines ..
EXTERNAL DGEMM, DGEQR2, DGERQ2, DLACPY, DLAGV2, DLARTG,
$ DLASET, DLASSQ, DORG2R, DORGR2, DORM2R, DORMR2,
$ DROT, DSCAL, DTGSY2
* ..
* .. Intrinsic Functions ..
INTRINSIC ABS, MAX, SQRT
* ..
* .. Executable Statements ..
*
INFO = 0
*
* Quick return if possible
*
IF( N.LE.1 .OR. N1.LE.0 .OR. N2.LE.0 )
$ RETURN
IF( N1.GT.N .OR. ( J1+N1 ).GT.N )
$ RETURN
M = N1 + N2
IF( LWORK.LT.MAX( 1, N*M, M*M*2 ) ) THEN
INFO = -16
WORK( 1 ) = MAX( 1, N*M, M*M*2 )
RETURN
END IF
*
WEAK = .FALSE.
DTRONG = .FALSE.
*
* Make a local copy of selected block
*
CALL DLASET( 'Full', LDST, LDST, ZERO, ZERO, LI, LDST )
CALL DLASET( 'Full', LDST, LDST, ZERO, ZERO, IR, LDST )
CALL DLACPY( 'Full', M, M, A( J1, J1 ), LDA, S, LDST )
CALL DLACPY( 'Full', M, M, B( J1, J1 ), LDB, T, LDST )
*
* Compute threshold for testing acceptance of swapping.
*
EPS = DLAMCH( 'P' )
SMLNUM = DLAMCH( 'S' ) / EPS
DSCALE = ZERO
DSUM = ONE
CALL DLACPY( 'Full', M, M, S, LDST, WORK, M )
CALL DLASSQ( M*M, WORK, 1, DSCALE, DSUM )
CALL DLACPY( 'Full', M, M, T, LDST, WORK, M )
CALL DLASSQ( M*M, WORK, 1, DSCALE, DSUM )
DNORM = DSCALE*SQRT( DSUM )
*
* THRES has been changed from
* THRESH = MAX( TEN*EPS*SA, SMLNUM )
* to
* THRESH = MAX( TWENTY*EPS*SA, SMLNUM )
* on 04/01/10.
* "Bug" reported by Ondra Kamenik, confirmed by Julie Langou, fixed by
* Jim Demmel and Guillaume Revy. See forum post 1783.
*
THRESH = MAX( TWENTY*EPS*DNORM, SMLNUM )
*
IF( M.EQ.2 ) THEN
*
* CASE 1: Swap 1-by-1 and 1-by-1 blocks.
*
* Compute orthogonal QL and RQ that swap 1-by-1 and 1-by-1 blocks
* using Givens rotations and perform the swap tentatively.
*
F = S( 2, 2 )*T( 1, 1 ) - T( 2, 2 )*S( 1, 1 )
G = S( 2, 2 )*T( 1, 2 ) - T( 2, 2 )*S( 1, 2 )
SB = ABS( T( 2, 2 ) )
SA = ABS( S( 2, 2 ) )
CALL DLARTG( F, G, IR( 1, 2 ), IR( 1, 1 ), DDUM )
IR( 2, 1 ) = -IR( 1, 2 )
IR( 2, 2 ) = IR( 1, 1 )
CALL DROT( 2, S( 1, 1 ), 1, S( 1, 2 ), 1, IR( 1, 1 ),
$ IR( 2, 1 ) )
CALL DROT( 2, T( 1, 1 ), 1, T( 1, 2 ), 1, IR( 1, 1 ),
$ IR( 2, 1 ) )
IF( SA.GE.SB ) THEN
CALL DLARTG( S( 1, 1 ), S( 2, 1 ), LI( 1, 1 ), LI( 2, 1 ),
$ DDUM )
ELSE
CALL DLARTG( T( 1, 1 ), T( 2, 1 ), LI( 1, 1 ), LI( 2, 1 ),
$ DDUM )
END IF
CALL DROT( 2, S( 1, 1 ), LDST, S( 2, 1 ), LDST, LI( 1, 1 ),
$ LI( 2, 1 ) )
CALL DROT( 2, T( 1, 1 ), LDST, T( 2, 1 ), LDST, LI( 1, 1 ),
$ LI( 2, 1 ) )
LI( 2, 2 ) = LI( 1, 1 )
LI( 1, 2 ) = -LI( 2, 1 )
*
* Weak stability test:
* |S21| + |T21| <= O(EPS * F-norm((S, T)))
*
WS = ABS( S( 2, 1 ) ) + ABS( T( 2, 1 ) )
WEAK = WS.LE.THRESH
IF( .NOT.WEAK )
$ GO TO 70
*
IF( WANDS ) THEN
*
* Strong stability test:
* F-norm((A-QL**T*S*QR, B-QL**T*T*QR)) <= O(EPS*F-norm((A,B)))
*
CALL DLACPY( 'Full', M, M, A( J1, J1 ), LDA, WORK( M*M+1 ),
$ M )
CALL DGEMM( 'N', 'N', M, M, M, ONE, LI, LDST, S, LDST, ZERO,
$ WORK, M )
CALL DGEMM( 'N', 'T', M, M, M, -ONE, WORK, M, IR, LDST, ONE,
$ WORK( M*M+1 ), M )
DSCALE = ZERO
DSUM = ONE
CALL DLASSQ( M*M, WORK( M*M+1 ), 1, DSCALE, DSUM )
*
CALL DLACPY( 'Full', M, M, B( J1, J1 ), LDB, WORK( M*M+1 ),
$ M )
CALL DGEMM( 'N', 'N', M, M, M, ONE, LI, LDST, T, LDST, ZERO,
$ WORK, M )
CALL DGEMM( 'N', 'T', M, M, M, -ONE, WORK, M, IR, LDST, ONE,
$ WORK( M*M+1 ), M )
CALL DLASSQ( M*M, WORK( M*M+1 ), 1, DSCALE, DSUM )
SS = DSCALE*SQRT( DSUM )
DTRONG = SS.LE.THRESH
IF( .NOT.DTRONG )
$ GO TO 70
END IF
*
* Update (A(J1:J1+M-1, M+J1:N), B(J1:J1+M-1, M+J1:N)) and
* (A(1:J1-1, J1:J1+M), B(1:J1-1, J1:J1+M)).
*
CALL DROT( J1+1, A( 1, J1 ), 1, A( 1, J1+1 ), 1, IR( 1, 1 ),
$ IR( 2, 1 ) )
CALL DROT( J1+1, B( 1, J1 ), 1, B( 1, J1+1 ), 1, IR( 1, 1 ),
$ IR( 2, 1 ) )
CALL DROT( N-J1+1, A( J1, J1 ), LDA, A( J1+1, J1 ), LDA,
$ LI( 1, 1 ), LI( 2, 1 ) )
CALL DROT( N-J1+1, B( J1, J1 ), LDB, B( J1+1, J1 ), LDB,
$ LI( 1, 1 ), LI( 2, 1 ) )
*
* Set N1-by-N2 (2,1) - blocks to ZERO.
*
A( J1+1, J1 ) = ZERO
B( J1+1, J1 ) = ZERO
*
* Accumulate transformations into Q and Z if requested.
*
IF( WANTZ )
$ CALL DROT( N, Z( 1, J1 ), 1, Z( 1, J1+1 ), 1, IR( 1, 1 ),
$ IR( 2, 1 ) )
IF( WANTQ )
$ CALL DROT( N, Q( 1, J1 ), 1, Q( 1, J1+1 ), 1, LI( 1, 1 ),
$ LI( 2, 1 ) )
*
* Exit with INFO = 0 if swap was successfully performed.
*
RETURN
*
ELSE
*
* CASE 2: Swap 1-by-1 and 2-by-2 blocks, or 2-by-2
* and 2-by-2 blocks.
*
* Solve the generalized Sylvester equation
* S11 * R - L * S22 = SCALE * S12
* T11 * R - L * T22 = SCALE * T12
* for R and L. Solutions in LI and IR.
*
CALL DLACPY( 'Full', N1, N2, T( 1, N1+1 ), LDST, LI, LDST )
CALL DLACPY( 'Full', N1, N2, S( 1, N1+1 ), LDST,
$ IR( N2+1, N1+1 ), LDST )
CALL DTGSY2( 'N', 0, N1, N2, S, LDST, S( N1+1, N1+1 ), LDST,
$ IR( N2+1, N1+1 ), LDST, T, LDST, T( N1+1, N1+1 ),
$ LDST, LI, LDST, SCALE, DSUM, DSCALE, IWORK, IDUM,
$ LINFO )
*
* Compute orthogonal matrix QL:
*
* QL**T * LI = [ TL ]
* [ 0 ]
* where
* LI = [ -L ]
* [ SCALE * identity(N2) ]
*
DO 10 I = 1, N2
CALL DSCAL( N1, -ONE, LI( 1, I ), 1 )
LI( N1+I, I ) = SCALE
10 CONTINUE
CALL DGEQR2( M, N2, LI, LDST, TAUL, WORK, LINFO )
IF( LINFO.NE.0 )
$ GO TO 70
CALL DORG2R( M, M, N2, LI, LDST, TAUL, WORK, LINFO )
IF( LINFO.NE.0 )
$ GO TO 70
*
* Compute orthogonal matrix RQ:
*
* IR * RQ**T = [ 0 TR],
*
* where IR = [ SCALE * identity(N1), R ]
*
DO 20 I = 1, N1
IR( N2+I, I ) = SCALE
20 CONTINUE
CALL DGERQ2( N1, M, IR( N2+1, 1 ), LDST, TAUR, WORK, LINFO )
IF( LINFO.NE.0 )
$ GO TO 70
CALL DORGR2( M, M, N1, IR, LDST, TAUR, WORK, LINFO )
IF( LINFO.NE.0 )
$ GO TO 70
*
* Perform the swapping tentatively:
*
CALL DGEMM( 'T', 'N', M, M, M, ONE, LI, LDST, S, LDST, ZERO,
$ WORK, M )
CALL DGEMM( 'N', 'T', M, M, M, ONE, WORK, M, IR, LDST, ZERO, S,
$ LDST )
CALL DGEMM( 'T', 'N', M, M, M, ONE, LI, LDST, T, LDST, ZERO,
$ WORK, M )
CALL DGEMM( 'N', 'T', M, M, M, ONE, WORK, M, IR, LDST, ZERO, T,
$ LDST )
CALL DLACPY( 'F', M, M, S, LDST, SCPY, LDST )
CALL DLACPY( 'F', M, M, T, LDST, TCPY, LDST )
CALL DLACPY( 'F', M, M, IR, LDST, IRCOP, LDST )
CALL DLACPY( 'F', M, M, LI, LDST, LICOP, LDST )
*
* Triangularize the B-part by an RQ factorization.
* Apply transformation (from left) to A-part, giving S.
*
CALL DGERQ2( M, M, T, LDST, TAUR, WORK, LINFO )
IF( LINFO.NE.0 )
$ GO TO 70
CALL DORMR2( 'R', 'T', M, M, M, T, LDST, TAUR, S, LDST, WORK,
$ LINFO )
IF( LINFO.NE.0 )
$ GO TO 70
CALL DORMR2( 'L', 'N', M, M, M, T, LDST, TAUR, IR, LDST, WORK,
$ LINFO )
IF( LINFO.NE.0 )
$ GO TO 70
*
* Compute F-norm(S21) in BRQA21. (T21 is 0.)
*
DSCALE = ZERO
DSUM = ONE
DO 30 I = 1, N2
CALL DLASSQ( N1, S( N2+1, I ), 1, DSCALE, DSUM )
30 CONTINUE
BRQA21 = DSCALE*SQRT( DSUM )
*
* Triangularize the B-part by a QR factorization.
* Apply transformation (from right) to A-part, giving S.
*
CALL DGEQR2( M, M, TCPY, LDST, TAUL, WORK, LINFO )
IF( LINFO.NE.0 )
$ GO TO 70
CALL DORM2R( 'L', 'T', M, M, M, TCPY, LDST, TAUL, SCPY, LDST,
$ WORK, INFO )
CALL DORM2R( 'R', 'N', M, M, M, TCPY, LDST, TAUL, LICOP, LDST,
$ WORK, INFO )
IF( LINFO.NE.0 )
$ GO TO 70
*
* Compute F-norm(S21) in BQRA21. (T21 is 0.)
*
DSCALE = ZERO
DSUM = ONE
DO 40 I = 1, N2
CALL DLASSQ( N1, SCPY( N2+1, I ), 1, DSCALE, DSUM )
40 CONTINUE
BQRA21 = DSCALE*SQRT( DSUM )
*
* Decide which method to use.
* Weak stability test:
* F-norm(S21) <= O(EPS * F-norm((S, T)))
*
IF( BQRA21.LE.BRQA21 .AND. BQRA21.LE.THRESH ) THEN
CALL DLACPY( 'F', M, M, SCPY, LDST, S, LDST )
CALL DLACPY( 'F', M, M, TCPY, LDST, T, LDST )
CALL DLACPY( 'F', M, M, IRCOP, LDST, IR, LDST )
CALL DLACPY( 'F', M, M, LICOP, LDST, LI, LDST )
ELSE IF( BRQA21.GE.THRESH ) THEN
GO TO 70
END IF
*
* Set lower triangle of B-part to zero
*
CALL DLASET( 'Lower', M-1, M-1, ZERO, ZERO, T(2,1), LDST )
*
IF( WANDS ) THEN
*
* Strong stability test:
* F-norm((A-QL*S*QR**T, B-QL*T*QR**T)) <= O(EPS*F-norm((A,B)))
*
CALL DLACPY( 'Full', M, M, A( J1, J1 ), LDA, WORK( M*M+1 ),
$ M )
CALL DGEMM( 'N', 'N', M, M, M, ONE, LI, LDST, S, LDST, ZERO,
$ WORK, M )
CALL DGEMM( 'N', 'N', M, M, M, -ONE, WORK, M, IR, LDST, ONE,
$ WORK( M*M+1 ), M )
DSCALE = ZERO
DSUM = ONE
CALL DLASSQ( M*M, WORK( M*M+1 ), 1, DSCALE, DSUM )
*
CALL DLACPY( 'Full', M, M, B( J1, J1 ), LDB, WORK( M*M+1 ),
$ M )
CALL DGEMM( 'N', 'N', M, M, M, ONE, LI, LDST, T, LDST, ZERO,
$ WORK, M )
CALL DGEMM( 'N', 'N', M, M, M, -ONE, WORK, M, IR, LDST, ONE,
$ WORK( M*M+1 ), M )
CALL DLASSQ( M*M, WORK( M*M+1 ), 1, DSCALE, DSUM )
SS = DSCALE*SQRT( DSUM )
DTRONG = ( SS.LE.THRESH )
IF( .NOT.DTRONG )
$ GO TO 70
*
END IF
*
* If the swap is accepted ("weakly" and "strongly"), apply the
* transformations and set N1-by-N2 (2,1)-block to zero.
*
CALL DLASET( 'Full', N1, N2, ZERO, ZERO, S(N2+1,1), LDST )
*
* copy back M-by-M diagonal block starting at index J1 of (A, B)
*
CALL DLACPY( 'F', M, M, S, LDST, A( J1, J1 ), LDA )
CALL DLACPY( 'F', M, M, T, LDST, B( J1, J1 ), LDB )
CALL DLASET( 'Full', LDST, LDST, ZERO, ZERO, T, LDST )
*
* Standardize existing 2-by-2 blocks.
*
DO 50 I = 1, M*M
WORK(I) = ZERO
50 CONTINUE
WORK( 1 ) = ONE
T( 1, 1 ) = ONE
IDUM = LWORK - M*M - 2
IF( N2.GT.1 ) THEN
CALL DLAGV2( A( J1, J1 ), LDA, B( J1, J1 ), LDB, AR, AI, BE,
$ WORK( 1 ), WORK( 2 ), T( 1, 1 ), T( 2, 1 ) )
WORK( M+1 ) = -WORK( 2 )
WORK( M+2 ) = WORK( 1 )
T( N2, N2 ) = T( 1, 1 )
T( 1, 2 ) = -T( 2, 1 )
END IF
WORK( M*M ) = ONE
T( M, M ) = ONE
*
IF( N1.GT.1 ) THEN
CALL DLAGV2( A( J1+N2, J1+N2 ), LDA, B( J1+N2, J1+N2 ), LDB,
$ TAUR, TAUL, WORK( M*M+1 ), WORK( N2*M+N2+1 ),
$ WORK( N2*M+N2+2 ), T( N2+1, N2+1 ),
$ T( M, M-1 ) )
WORK( M*M ) = WORK( N2*M+N2+1 )
WORK( M*M-1 ) = -WORK( N2*M+N2+2 )
T( M, M ) = T( N2+1, N2+1 )
T( M-1, M ) = -T( M, M-1 )
END IF
CALL DGEMM( 'T', 'N', N2, N1, N2, ONE, WORK, M, A( J1, J1+N2 ),
$ LDA, ZERO, WORK( M*M+1 ), N2 )
CALL DLACPY( 'Full', N2, N1, WORK( M*M+1 ), N2, A( J1, J1+N2 ),
$ LDA )
CALL DGEMM( 'T', 'N', N2, N1, N2, ONE, WORK, M, B( J1, J1+N2 ),
$ LDB, ZERO, WORK( M*M+1 ), N2 )
CALL DLACPY( 'Full', N2, N1, WORK( M*M+1 ), N2, B( J1, J1+N2 ),
$ LDB )
CALL DGEMM( 'N', 'N', M, M, M, ONE, LI, LDST, WORK, M, ZERO,
$ WORK( M*M+1 ), M )
CALL DLACPY( 'Full', M, M, WORK( M*M+1 ), M, LI, LDST )
CALL DGEMM( 'N', 'N', N2, N1, N1, ONE, A( J1, J1+N2 ), LDA,
$ T( N2+1, N2+1 ), LDST, ZERO, WORK, N2 )
CALL DLACPY( 'Full', N2, N1, WORK, N2, A( J1, J1+N2 ), LDA )
CALL DGEMM( 'N', 'N', N2, N1, N1, ONE, B( J1, J1+N2 ), LDB,
$ T( N2+1, N2+1 ), LDST, ZERO, WORK, N2 )
CALL DLACPY( 'Full', N2, N1, WORK, N2, B( J1, J1+N2 ), LDB )
CALL DGEMM( 'T', 'N', M, M, M, ONE, IR, LDST, T, LDST, ZERO,
$ WORK, M )
CALL DLACPY( 'Full', M, M, WORK, M, IR, LDST )
*
* Accumulate transformations into Q and Z if requested.
*
IF( WANTQ ) THEN
CALL DGEMM( 'N', 'N', N, M, M, ONE, Q( 1, J1 ), LDQ, LI,
$ LDST, ZERO, WORK, N )
CALL DLACPY( 'Full', N, M, WORK, N, Q( 1, J1 ), LDQ )
*
END IF
*
IF( WANTZ ) THEN
CALL DGEMM( 'N', 'N', N, M, M, ONE, Z( 1, J1 ), LDZ, IR,
$ LDST, ZERO, WORK, N )
CALL DLACPY( 'Full', N, M, WORK, N, Z( 1, J1 ), LDZ )
*
END IF
*
* Update (A(J1:J1+M-1, M+J1:N), B(J1:J1+M-1, M+J1:N)) and
* (A(1:J1-1, J1:J1+M), B(1:J1-1, J1:J1+M)).
*
I = J1 + M
IF( I.LE.N ) THEN
CALL DGEMM( 'T', 'N', M, N-I+1, M, ONE, LI, LDST,
$ A( J1, I ), LDA, ZERO, WORK, M )
CALL DLACPY( 'Full', M, N-I+1, WORK, M, A( J1, I ), LDA )
CALL DGEMM( 'T', 'N', M, N-I+1, M, ONE, LI, LDST,
$ B( J1, I ), LDA, ZERO, WORK, M )
CALL DLACPY( 'Full', M, N-I+1, WORK, M, B( J1, I ), LDB )
END IF
I = J1 - 1
IF( I.GT.0 ) THEN
CALL DGEMM( 'N', 'N', I, M, M, ONE, A( 1, J1 ), LDA, IR,
$ LDST, ZERO, WORK, I )
CALL DLACPY( 'Full', I, M, WORK, I, A( 1, J1 ), LDA )
CALL DGEMM( 'N', 'N', I, M, M, ONE, B( 1, J1 ), LDB, IR,
$ LDST, ZERO, WORK, I )
CALL DLACPY( 'Full', I, M, WORK, I, B( 1, J1 ), LDB )
END IF
*
* Exit with INFO = 0 if swap was successfully performed.
*
RETURN
*
END IF
*
* Exit with INFO = 1 if swap was rejected.
*
70 CONTINUE
*
INFO = 1
RETURN
*
* End of DTGEX2
*
END
|