You can not select more than 25 topics Topics must start with a chinese character,a letter or number, can include dashes ('-') and can be up to 35 characters long.

dtgex2.f 25 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711
  1. *> \brief \b DTGEX2 swaps adjacent diagonal blocks in an upper (quasi) triangular matrix pair by an orthogonal equivalence transformation.
  2. *
  3. * =========== DOCUMENTATION ===========
  4. *
  5. * Online html documentation available at
  6. * http://www.netlib.org/lapack/explore-html/
  7. *
  8. *> \htmlonly
  9. *> Download DTGEX2 + dependencies
  10. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dtgex2.f">
  11. *> [TGZ]</a>
  12. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dtgex2.f">
  13. *> [ZIP]</a>
  14. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dtgex2.f">
  15. *> [TXT]</a>
  16. *> \endhtmlonly
  17. *
  18. * Definition:
  19. * ===========
  20. *
  21. * SUBROUTINE DTGEX2( WANTQ, WANTZ, N, A, LDA, B, LDB, Q, LDQ, Z,
  22. * LDZ, J1, N1, N2, WORK, LWORK, INFO )
  23. *
  24. * .. Scalar Arguments ..
  25. * LOGICAL WANTQ, WANTZ
  26. * INTEGER INFO, J1, LDA, LDB, LDQ, LDZ, LWORK, N, N1, N2
  27. * ..
  28. * .. Array Arguments ..
  29. * DOUBLE PRECISION A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
  30. * $ WORK( * ), Z( LDZ, * )
  31. * ..
  32. *
  33. *
  34. *> \par Purpose:
  35. * =============
  36. *>
  37. *> \verbatim
  38. *>
  39. *> DTGEX2 swaps adjacent diagonal blocks (A11, B11) and (A22, B22)
  40. *> of size 1-by-1 or 2-by-2 in an upper (quasi) triangular matrix pair
  41. *> (A, B) by an orthogonal equivalence transformation.
  42. *>
  43. *> (A, B) must be in generalized real Schur canonical form (as returned
  44. *> by DGGES), i.e. A is block upper triangular with 1-by-1 and 2-by-2
  45. *> diagonal blocks. B is upper triangular.
  46. *>
  47. *> Optionally, the matrices Q and Z of generalized Schur vectors are
  48. *> updated.
  49. *>
  50. *> Q(in) * A(in) * Z(in)**T = Q(out) * A(out) * Z(out)**T
  51. *> Q(in) * B(in) * Z(in)**T = Q(out) * B(out) * Z(out)**T
  52. *>
  53. *> \endverbatim
  54. *
  55. * Arguments:
  56. * ==========
  57. *
  58. *> \param[in] WANTQ
  59. *> \verbatim
  60. *> WANTQ is LOGICAL
  61. *> .TRUE. : update the left transformation matrix Q;
  62. *> .FALSE.: do not update Q.
  63. *> \endverbatim
  64. *>
  65. *> \param[in] WANTZ
  66. *> \verbatim
  67. *> WANTZ is LOGICAL
  68. *> .TRUE. : update the right transformation matrix Z;
  69. *> .FALSE.: do not update Z.
  70. *> \endverbatim
  71. *>
  72. *> \param[in] N
  73. *> \verbatim
  74. *> N is INTEGER
  75. *> The order of the matrices A and B. N >= 0.
  76. *> \endverbatim
  77. *>
  78. *> \param[in,out] A
  79. *> \verbatim
  80. *> A is DOUBLE PRECISION array, dimensions (LDA,N)
  81. *> On entry, the matrix A in the pair (A, B).
  82. *> On exit, the updated matrix A.
  83. *> \endverbatim
  84. *>
  85. *> \param[in] LDA
  86. *> \verbatim
  87. *> LDA is INTEGER
  88. *> The leading dimension of the array A. LDA >= max(1,N).
  89. *> \endverbatim
  90. *>
  91. *> \param[in,out] B
  92. *> \verbatim
  93. *> B is DOUBLE PRECISION array, dimensions (LDB,N)
  94. *> On entry, the matrix B in the pair (A, B).
  95. *> On exit, the updated matrix B.
  96. *> \endverbatim
  97. *>
  98. *> \param[in] LDB
  99. *> \verbatim
  100. *> LDB is INTEGER
  101. *> The leading dimension of the array B. LDB >= max(1,N).
  102. *> \endverbatim
  103. *>
  104. *> \param[in,out] Q
  105. *> \verbatim
  106. *> Q is DOUBLE PRECISION array, dimension (LDQ,N)
  107. *> On entry, if WANTQ = .TRUE., the orthogonal matrix Q.
  108. *> On exit, the updated matrix Q.
  109. *> Not referenced if WANTQ = .FALSE..
  110. *> \endverbatim
  111. *>
  112. *> \param[in] LDQ
  113. *> \verbatim
  114. *> LDQ is INTEGER
  115. *> The leading dimension of the array Q. LDQ >= 1.
  116. *> If WANTQ = .TRUE., LDQ >= N.
  117. *> \endverbatim
  118. *>
  119. *> \param[in,out] Z
  120. *> \verbatim
  121. *> Z is DOUBLE PRECISION array, dimension (LDZ,N)
  122. *> On entry, if WANTZ =.TRUE., the orthogonal matrix Z.
  123. *> On exit, the updated matrix Z.
  124. *> Not referenced if WANTZ = .FALSE..
  125. *> \endverbatim
  126. *>
  127. *> \param[in] LDZ
  128. *> \verbatim
  129. *> LDZ is INTEGER
  130. *> The leading dimension of the array Z. LDZ >= 1.
  131. *> If WANTZ = .TRUE., LDZ >= N.
  132. *> \endverbatim
  133. *>
  134. *> \param[in] J1
  135. *> \verbatim
  136. *> J1 is INTEGER
  137. *> The index to the first block (A11, B11). 1 <= J1 <= N.
  138. *> \endverbatim
  139. *>
  140. *> \param[in] N1
  141. *> \verbatim
  142. *> N1 is INTEGER
  143. *> The order of the first block (A11, B11). N1 = 0, 1 or 2.
  144. *> \endverbatim
  145. *>
  146. *> \param[in] N2
  147. *> \verbatim
  148. *> N2 is INTEGER
  149. *> The order of the second block (A22, B22). N2 = 0, 1 or 2.
  150. *> \endverbatim
  151. *>
  152. *> \param[out] WORK
  153. *> \verbatim
  154. *> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK)).
  155. *> \endverbatim
  156. *>
  157. *> \param[in] LWORK
  158. *> \verbatim
  159. *> LWORK is INTEGER
  160. *> The dimension of the array WORK.
  161. *> LWORK >= MAX( 1, N*(N2+N1), (N2+N1)*(N2+N1)*2 )
  162. *> \endverbatim
  163. *>
  164. *> \param[out] INFO
  165. *> \verbatim
  166. *> INFO is INTEGER
  167. *> =0: Successful exit
  168. *> >0: If INFO = 1, the transformed matrix (A, B) would be
  169. *> too far from generalized Schur form; the blocks are
  170. *> not swapped and (A, B) and (Q, Z) are unchanged.
  171. *> The problem of swapping is too ill-conditioned.
  172. *> <0: If INFO = -16: LWORK is too small. Appropriate value
  173. *> for LWORK is returned in WORK(1).
  174. *> \endverbatim
  175. *
  176. * Authors:
  177. * ========
  178. *
  179. *> \author Univ. of Tennessee
  180. *> \author Univ. of California Berkeley
  181. *> \author Univ. of Colorado Denver
  182. *> \author NAG Ltd.
  183. *
  184. *> \ingroup doubleGEauxiliary
  185. *
  186. *> \par Further Details:
  187. * =====================
  188. *>
  189. *> In the current code both weak and strong stability tests are
  190. *> performed. The user can omit the strong stability test by changing
  191. *> the internal logical parameter WANDS to .FALSE.. See ref. [2] for
  192. *> details.
  193. *
  194. *> \par Contributors:
  195. * ==================
  196. *>
  197. *> Bo Kagstrom and Peter Poromaa, Department of Computing Science,
  198. *> Umea University, S-901 87 Umea, Sweden.
  199. *
  200. *> \par References:
  201. * ================
  202. *>
  203. *> \verbatim
  204. *>
  205. *> [1] B. Kagstrom; A Direct Method for Reordering Eigenvalues in the
  206. *> Generalized Real Schur Form of a Regular Matrix Pair (A, B), in
  207. *> M.S. Moonen et al (eds), Linear Algebra for Large Scale and
  208. *> Real-Time Applications, Kluwer Academic Publ. 1993, pp 195-218.
  209. *>
  210. *> [2] B. Kagstrom and P. Poromaa; Computing Eigenspaces with Specified
  211. *> Eigenvalues of a Regular Matrix Pair (A, B) and Condition
  212. *> Estimation: Theory, Algorithms and Software,
  213. *> Report UMINF - 94.04, Department of Computing Science, Umea
  214. *> University, S-901 87 Umea, Sweden, 1994. Also as LAPACK Working
  215. *> Note 87. To appear in Numerical Algorithms, 1996.
  216. *> \endverbatim
  217. *>
  218. * =====================================================================
  219. SUBROUTINE DTGEX2( WANTQ, WANTZ, N, A, LDA, B, LDB, Q, LDQ, Z,
  220. $ LDZ, J1, N1, N2, WORK, LWORK, INFO )
  221. *
  222. * -- LAPACK auxiliary routine --
  223. * -- LAPACK is a software package provided by Univ. of Tennessee, --
  224. * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  225. *
  226. * .. Scalar Arguments ..
  227. LOGICAL WANTQ, WANTZ
  228. INTEGER INFO, J1, LDA, LDB, LDQ, LDZ, LWORK, N, N1, N2
  229. * ..
  230. * .. Array Arguments ..
  231. DOUBLE PRECISION A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
  232. $ WORK( * ), Z( LDZ, * )
  233. * ..
  234. *
  235. * =====================================================================
  236. * Replaced various illegal calls to DCOPY by calls to DLASET, or by DO
  237. * loops. Sven Hammarling, 1/5/02.
  238. *
  239. * .. Parameters ..
  240. DOUBLE PRECISION ZERO, ONE
  241. PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 )
  242. DOUBLE PRECISION TWENTY
  243. PARAMETER ( TWENTY = 2.0D+01 )
  244. INTEGER LDST
  245. PARAMETER ( LDST = 4 )
  246. LOGICAL WANDS
  247. PARAMETER ( WANDS = .TRUE. )
  248. * ..
  249. * .. Local Scalars ..
  250. LOGICAL STRONG, WEAK
  251. INTEGER I, IDUM, LINFO, M
  252. DOUBLE PRECISION BQRA21, BRQA21, DDUM, DNORMA, DNORMB, DSCALE,
  253. $ DSUM, EPS, F, G, SA, SB, SCALE, SMLNUM,
  254. $ THRESHA, THRESHB
  255. * ..
  256. * .. Local Arrays ..
  257. INTEGER IWORK( LDST + 2 )
  258. DOUBLE PRECISION AI( 2 ), AR( 2 ), BE( 2 ), IR( LDST, LDST ),
  259. $ IRCOP( LDST, LDST ), LI( LDST, LDST ),
  260. $ LICOP( LDST, LDST ), S( LDST, LDST ),
  261. $ SCPY( LDST, LDST ), T( LDST, LDST ),
  262. $ TAUL( LDST ), TAUR( LDST ), TCPY( LDST, LDST )
  263. * ..
  264. * .. External Functions ..
  265. DOUBLE PRECISION DLAMCH
  266. EXTERNAL DLAMCH
  267. * ..
  268. * .. External Subroutines ..
  269. EXTERNAL DGEMM, DGEQR2, DGERQ2, DLACPY, DLAGV2, DLARTG,
  270. $ DLASET, DLASSQ, DORG2R, DORGR2, DORM2R, DORMR2,
  271. $ DROT, DSCAL, DTGSY2
  272. * ..
  273. * .. Intrinsic Functions ..
  274. INTRINSIC ABS, MAX, SQRT
  275. * ..
  276. * .. Executable Statements ..
  277. *
  278. INFO = 0
  279. *
  280. * Quick return if possible
  281. *
  282. IF( N.LE.1 .OR. N1.LE.0 .OR. N2.LE.0 )
  283. $ RETURN
  284. IF( N1.GT.N .OR. ( J1+N1 ).GT.N )
  285. $ RETURN
  286. M = N1 + N2
  287. IF( LWORK.LT.MAX( 1, N*M, M*M*2 ) ) THEN
  288. INFO = -16
  289. WORK( 1 ) = MAX( 1, N*M, M*M*2 )
  290. RETURN
  291. END IF
  292. *
  293. WEAK = .FALSE.
  294. STRONG = .FALSE.
  295. *
  296. * Make a local copy of selected block
  297. *
  298. CALL DLASET( 'Full', LDST, LDST, ZERO, ZERO, LI, LDST )
  299. CALL DLASET( 'Full', LDST, LDST, ZERO, ZERO, IR, LDST )
  300. CALL DLACPY( 'Full', M, M, A( J1, J1 ), LDA, S, LDST )
  301. CALL DLACPY( 'Full', M, M, B( J1, J1 ), LDB, T, LDST )
  302. *
  303. * Compute threshold for testing acceptance of swapping.
  304. *
  305. EPS = DLAMCH( 'P' )
  306. SMLNUM = DLAMCH( 'S' ) / EPS
  307. DSCALE = ZERO
  308. DSUM = ONE
  309. CALL DLACPY( 'Full', M, M, S, LDST, WORK, M )
  310. CALL DLASSQ( M*M, WORK, 1, DSCALE, DSUM )
  311. DNORMA = DSCALE*SQRT( DSUM )
  312. DSCALE = ZERO
  313. DSUM = ONE
  314. CALL DLACPY( 'Full', M, M, T, LDST, WORK, M )
  315. CALL DLASSQ( M*M, WORK, 1, DSCALE, DSUM )
  316. DNORMB = DSCALE*SQRT( DSUM )
  317. *
  318. * THRES has been changed from
  319. * THRESH = MAX( TEN*EPS*SA, SMLNUM )
  320. * to
  321. * THRESH = MAX( TWENTY*EPS*SA, SMLNUM )
  322. * on 04/01/10.
  323. * "Bug" reported by Ondra Kamenik, confirmed by Julie Langou, fixed by
  324. * Jim Demmel and Guillaume Revy. See forum post 1783.
  325. *
  326. THRESHA = MAX( TWENTY*EPS*DNORMA, SMLNUM )
  327. THRESHB = MAX( TWENTY*EPS*DNORMB, SMLNUM )
  328. *
  329. IF( M.EQ.2 ) THEN
  330. *
  331. * CASE 1: Swap 1-by-1 and 1-by-1 blocks.
  332. *
  333. * Compute orthogonal QL and RQ that swap 1-by-1 and 1-by-1 blocks
  334. * using Givens rotations and perform the swap tentatively.
  335. *
  336. F = S( 2, 2 )*T( 1, 1 ) - T( 2, 2 )*S( 1, 1 )
  337. G = S( 2, 2 )*T( 1, 2 ) - T( 2, 2 )*S( 1, 2 )
  338. SA = ABS( S( 2, 2 ) ) * ABS( T( 1, 1 ) )
  339. SB = ABS( S( 1, 1 ) ) * ABS( T( 2, 2 ) )
  340. CALL DLARTG( F, G, IR( 1, 2 ), IR( 1, 1 ), DDUM )
  341. IR( 2, 1 ) = -IR( 1, 2 )
  342. IR( 2, 2 ) = IR( 1, 1 )
  343. CALL DROT( 2, S( 1, 1 ), 1, S( 1, 2 ), 1, IR( 1, 1 ),
  344. $ IR( 2, 1 ) )
  345. CALL DROT( 2, T( 1, 1 ), 1, T( 1, 2 ), 1, IR( 1, 1 ),
  346. $ IR( 2, 1 ) )
  347. IF( SA.GE.SB ) THEN
  348. CALL DLARTG( S( 1, 1 ), S( 2, 1 ), LI( 1, 1 ), LI( 2, 1 ),
  349. $ DDUM )
  350. ELSE
  351. CALL DLARTG( T( 1, 1 ), T( 2, 1 ), LI( 1, 1 ), LI( 2, 1 ),
  352. $ DDUM )
  353. END IF
  354. CALL DROT( 2, S( 1, 1 ), LDST, S( 2, 1 ), LDST, LI( 1, 1 ),
  355. $ LI( 2, 1 ) )
  356. CALL DROT( 2, T( 1, 1 ), LDST, T( 2, 1 ), LDST, LI( 1, 1 ),
  357. $ LI( 2, 1 ) )
  358. LI( 2, 2 ) = LI( 1, 1 )
  359. LI( 1, 2 ) = -LI( 2, 1 )
  360. *
  361. * Weak stability test: |S21| <= O(EPS F-norm((A)))
  362. * and |T21| <= O(EPS F-norm((B)))
  363. *
  364. WEAK = ABS( S( 2, 1 ) ) .LE. THRESHA .AND.
  365. $ ABS( T( 2, 1 ) ) .LE. THRESHB
  366. IF( .NOT.WEAK )
  367. $ GO TO 70
  368. *
  369. IF( WANDS ) THEN
  370. *
  371. * Strong stability test:
  372. * F-norm((A-QL**H*S*QR)) <= O(EPS*F-norm((A)))
  373. * and
  374. * F-norm((B-QL**H*T*QR)) <= O(EPS*F-norm((B)))
  375. *
  376. CALL DLACPY( 'Full', M, M, A( J1, J1 ), LDA, WORK( M*M+1 ),
  377. $ M )
  378. CALL DGEMM( 'N', 'N', M, M, M, ONE, LI, LDST, S, LDST, ZERO,
  379. $ WORK, M )
  380. CALL DGEMM( 'N', 'T', M, M, M, -ONE, WORK, M, IR, LDST, ONE,
  381. $ WORK( M*M+1 ), M )
  382. DSCALE = ZERO
  383. DSUM = ONE
  384. CALL DLASSQ( M*M, WORK( M*M+1 ), 1, DSCALE, DSUM )
  385. SA = DSCALE*SQRT( DSUM )
  386. *
  387. CALL DLACPY( 'Full', M, M, B( J1, J1 ), LDB, WORK( M*M+1 ),
  388. $ M )
  389. CALL DGEMM( 'N', 'N', M, M, M, ONE, LI, LDST, T, LDST, ZERO,
  390. $ WORK, M )
  391. CALL DGEMM( 'N', 'T', M, M, M, -ONE, WORK, M, IR, LDST, ONE,
  392. $ WORK( M*M+1 ), M )
  393. DSCALE = ZERO
  394. DSUM = ONE
  395. CALL DLASSQ( M*M, WORK( M*M+1 ), 1, DSCALE, DSUM )
  396. SB = DSCALE*SQRT( DSUM )
  397. STRONG = SA.LE.THRESHA .AND. SB.LE.THRESHB
  398. IF( .NOT.STRONG )
  399. $ GO TO 70
  400. END IF
  401. *
  402. * Update (A(J1:J1+M-1, M+J1:N), B(J1:J1+M-1, M+J1:N)) and
  403. * (A(1:J1-1, J1:J1+M), B(1:J1-1, J1:J1+M)).
  404. *
  405. CALL DROT( J1+1, A( 1, J1 ), 1, A( 1, J1+1 ), 1, IR( 1, 1 ),
  406. $ IR( 2, 1 ) )
  407. CALL DROT( J1+1, B( 1, J1 ), 1, B( 1, J1+1 ), 1, IR( 1, 1 ),
  408. $ IR( 2, 1 ) )
  409. CALL DROT( N-J1+1, A( J1, J1 ), LDA, A( J1+1, J1 ), LDA,
  410. $ LI( 1, 1 ), LI( 2, 1 ) )
  411. CALL DROT( N-J1+1, B( J1, J1 ), LDB, B( J1+1, J1 ), LDB,
  412. $ LI( 1, 1 ), LI( 2, 1 ) )
  413. *
  414. * Set N1-by-N2 (2,1) - blocks to ZERO.
  415. *
  416. A( J1+1, J1 ) = ZERO
  417. B( J1+1, J1 ) = ZERO
  418. *
  419. * Accumulate transformations into Q and Z if requested.
  420. *
  421. IF( WANTZ )
  422. $ CALL DROT( N, Z( 1, J1 ), 1, Z( 1, J1+1 ), 1, IR( 1, 1 ),
  423. $ IR( 2, 1 ) )
  424. IF( WANTQ )
  425. $ CALL DROT( N, Q( 1, J1 ), 1, Q( 1, J1+1 ), 1, LI( 1, 1 ),
  426. $ LI( 2, 1 ) )
  427. *
  428. * Exit with INFO = 0 if swap was successfully performed.
  429. *
  430. RETURN
  431. *
  432. ELSE
  433. *
  434. * CASE 2: Swap 1-by-1 and 2-by-2 blocks, or 2-by-2
  435. * and 2-by-2 blocks.
  436. *
  437. * Solve the generalized Sylvester equation
  438. * S11 * R - L * S22 = SCALE * S12
  439. * T11 * R - L * T22 = SCALE * T12
  440. * for R and L. Solutions in LI and IR.
  441. *
  442. CALL DLACPY( 'Full', N1, N2, T( 1, N1+1 ), LDST, LI, LDST )
  443. CALL DLACPY( 'Full', N1, N2, S( 1, N1+1 ), LDST,
  444. $ IR( N2+1, N1+1 ), LDST )
  445. CALL DTGSY2( 'N', 0, N1, N2, S, LDST, S( N1+1, N1+1 ), LDST,
  446. $ IR( N2+1, N1+1 ), LDST, T, LDST, T( N1+1, N1+1 ),
  447. $ LDST, LI, LDST, SCALE, DSUM, DSCALE, IWORK, IDUM,
  448. $ LINFO )
  449. IF( LINFO.NE.0 )
  450. $ GO TO 70
  451. *
  452. * Compute orthogonal matrix QL:
  453. *
  454. * QL**T * LI = [ TL ]
  455. * [ 0 ]
  456. * where
  457. * LI = [ -L ]
  458. * [ SCALE * identity(N2) ]
  459. *
  460. DO 10 I = 1, N2
  461. CALL DSCAL( N1, -ONE, LI( 1, I ), 1 )
  462. LI( N1+I, I ) = SCALE
  463. 10 CONTINUE
  464. CALL DGEQR2( M, N2, LI, LDST, TAUL, WORK, LINFO )
  465. IF( LINFO.NE.0 )
  466. $ GO TO 70
  467. CALL DORG2R( M, M, N2, LI, LDST, TAUL, WORK, LINFO )
  468. IF( LINFO.NE.0 )
  469. $ GO TO 70
  470. *
  471. * Compute orthogonal matrix RQ:
  472. *
  473. * IR * RQ**T = [ 0 TR],
  474. *
  475. * where IR = [ SCALE * identity(N1), R ]
  476. *
  477. DO 20 I = 1, N1
  478. IR( N2+I, I ) = SCALE
  479. 20 CONTINUE
  480. CALL DGERQ2( N1, M, IR( N2+1, 1 ), LDST, TAUR, WORK, LINFO )
  481. IF( LINFO.NE.0 )
  482. $ GO TO 70
  483. CALL DORGR2( M, M, N1, IR, LDST, TAUR, WORK, LINFO )
  484. IF( LINFO.NE.0 )
  485. $ GO TO 70
  486. *
  487. * Perform the swapping tentatively:
  488. *
  489. CALL DGEMM( 'T', 'N', M, M, M, ONE, LI, LDST, S, LDST, ZERO,
  490. $ WORK, M )
  491. CALL DGEMM( 'N', 'T', M, M, M, ONE, WORK, M, IR, LDST, ZERO, S,
  492. $ LDST )
  493. CALL DGEMM( 'T', 'N', M, M, M, ONE, LI, LDST, T, LDST, ZERO,
  494. $ WORK, M )
  495. CALL DGEMM( 'N', 'T', M, M, M, ONE, WORK, M, IR, LDST, ZERO, T,
  496. $ LDST )
  497. CALL DLACPY( 'F', M, M, S, LDST, SCPY, LDST )
  498. CALL DLACPY( 'F', M, M, T, LDST, TCPY, LDST )
  499. CALL DLACPY( 'F', M, M, IR, LDST, IRCOP, LDST )
  500. CALL DLACPY( 'F', M, M, LI, LDST, LICOP, LDST )
  501. *
  502. * Triangularize the B-part by an RQ factorization.
  503. * Apply transformation (from left) to A-part, giving S.
  504. *
  505. CALL DGERQ2( M, M, T, LDST, TAUR, WORK, LINFO )
  506. IF( LINFO.NE.0 )
  507. $ GO TO 70
  508. CALL DORMR2( 'R', 'T', M, M, M, T, LDST, TAUR, S, LDST, WORK,
  509. $ LINFO )
  510. IF( LINFO.NE.0 )
  511. $ GO TO 70
  512. CALL DORMR2( 'L', 'N', M, M, M, T, LDST, TAUR, IR, LDST, WORK,
  513. $ LINFO )
  514. IF( LINFO.NE.0 )
  515. $ GO TO 70
  516. *
  517. * Compute F-norm(S21) in BRQA21. (T21 is 0.)
  518. *
  519. DSCALE = ZERO
  520. DSUM = ONE
  521. DO 30 I = 1, N2
  522. CALL DLASSQ( N1, S( N2+1, I ), 1, DSCALE, DSUM )
  523. 30 CONTINUE
  524. BRQA21 = DSCALE*SQRT( DSUM )
  525. *
  526. * Triangularize the B-part by a QR factorization.
  527. * Apply transformation (from right) to A-part, giving S.
  528. *
  529. CALL DGEQR2( M, M, TCPY, LDST, TAUL, WORK, LINFO )
  530. IF( LINFO.NE.0 )
  531. $ GO TO 70
  532. CALL DORM2R( 'L', 'T', M, M, M, TCPY, LDST, TAUL, SCPY, LDST,
  533. $ WORK, INFO )
  534. CALL DORM2R( 'R', 'N', M, M, M, TCPY, LDST, TAUL, LICOP, LDST,
  535. $ WORK, INFO )
  536. IF( LINFO.NE.0 )
  537. $ GO TO 70
  538. *
  539. * Compute F-norm(S21) in BQRA21. (T21 is 0.)
  540. *
  541. DSCALE = ZERO
  542. DSUM = ONE
  543. DO 40 I = 1, N2
  544. CALL DLASSQ( N1, SCPY( N2+1, I ), 1, DSCALE, DSUM )
  545. 40 CONTINUE
  546. BQRA21 = DSCALE*SQRT( DSUM )
  547. *
  548. * Decide which method to use.
  549. * Weak stability test:
  550. * F-norm(S21) <= O(EPS * F-norm((S)))
  551. *
  552. IF( BQRA21.LE.BRQA21 .AND. BQRA21.LE.THRESHA ) THEN
  553. CALL DLACPY( 'F', M, M, SCPY, LDST, S, LDST )
  554. CALL DLACPY( 'F', M, M, TCPY, LDST, T, LDST )
  555. CALL DLACPY( 'F', M, M, IRCOP, LDST, IR, LDST )
  556. CALL DLACPY( 'F', M, M, LICOP, LDST, LI, LDST )
  557. ELSE IF( BRQA21.GE.THRESHA ) THEN
  558. GO TO 70
  559. END IF
  560. *
  561. * Set lower triangle of B-part to zero
  562. *
  563. CALL DLASET( 'Lower', M-1, M-1, ZERO, ZERO, T(2,1), LDST )
  564. *
  565. IF( WANDS ) THEN
  566. *
  567. * Strong stability test:
  568. * F-norm((A-QL**H*S*QR)) <= O(EPS*F-norm((A)))
  569. * and
  570. * F-norm((B-QL**H*T*QR)) <= O(EPS*F-norm((B)))
  571. *
  572. CALL DLACPY( 'Full', M, M, A( J1, J1 ), LDA, WORK( M*M+1 ),
  573. $ M )
  574. CALL DGEMM( 'N', 'N', M, M, M, ONE, LI, LDST, S, LDST, ZERO,
  575. $ WORK, M )
  576. CALL DGEMM( 'N', 'N', M, M, M, -ONE, WORK, M, IR, LDST, ONE,
  577. $ WORK( M*M+1 ), M )
  578. DSCALE = ZERO
  579. DSUM = ONE
  580. CALL DLASSQ( M*M, WORK( M*M+1 ), 1, DSCALE, DSUM )
  581. SA = DSCALE*SQRT( DSUM )
  582. *
  583. CALL DLACPY( 'Full', M, M, B( J1, J1 ), LDB, WORK( M*M+1 ),
  584. $ M )
  585. CALL DGEMM( 'N', 'N', M, M, M, ONE, LI, LDST, T, LDST, ZERO,
  586. $ WORK, M )
  587. CALL DGEMM( 'N', 'N', M, M, M, -ONE, WORK, M, IR, LDST, ONE,
  588. $ WORK( M*M+1 ), M )
  589. DSCALE = ZERO
  590. DSUM = ONE
  591. CALL DLASSQ( M*M, WORK( M*M+1 ), 1, DSCALE, DSUM )
  592. SB = DSCALE*SQRT( DSUM )
  593. STRONG = SA.LE.THRESHA .AND. SB.LE.THRESHB
  594. IF( .NOT.STRONG )
  595. $ GO TO 70
  596. *
  597. END IF
  598. *
  599. * If the swap is accepted ("weakly" and "strongly"), apply the
  600. * transformations and set N1-by-N2 (2,1)-block to zero.
  601. *
  602. CALL DLASET( 'Full', N1, N2, ZERO, ZERO, S(N2+1,1), LDST )
  603. *
  604. * copy back M-by-M diagonal block starting at index J1 of (A, B)
  605. *
  606. CALL DLACPY( 'F', M, M, S, LDST, A( J1, J1 ), LDA )
  607. CALL DLACPY( 'F', M, M, T, LDST, B( J1, J1 ), LDB )
  608. CALL DLASET( 'Full', LDST, LDST, ZERO, ZERO, T, LDST )
  609. *
  610. * Standardize existing 2-by-2 blocks.
  611. *
  612. CALL DLASET( 'Full', M, M, ZERO, ZERO, WORK, M )
  613. WORK( 1 ) = ONE
  614. T( 1, 1 ) = ONE
  615. IDUM = LWORK - M*M - 2
  616. IF( N2.GT.1 ) THEN
  617. CALL DLAGV2( A( J1, J1 ), LDA, B( J1, J1 ), LDB, AR, AI, BE,
  618. $ WORK( 1 ), WORK( 2 ), T( 1, 1 ), T( 2, 1 ) )
  619. WORK( M+1 ) = -WORK( 2 )
  620. WORK( M+2 ) = WORK( 1 )
  621. T( N2, N2 ) = T( 1, 1 )
  622. T( 1, 2 ) = -T( 2, 1 )
  623. END IF
  624. WORK( M*M ) = ONE
  625. T( M, M ) = ONE
  626. *
  627. IF( N1.GT.1 ) THEN
  628. CALL DLAGV2( A( J1+N2, J1+N2 ), LDA, B( J1+N2, J1+N2 ), LDB,
  629. $ TAUR, TAUL, WORK( M*M+1 ), WORK( N2*M+N2+1 ),
  630. $ WORK( N2*M+N2+2 ), T( N2+1, N2+1 ),
  631. $ T( M, M-1 ) )
  632. WORK( M*M ) = WORK( N2*M+N2+1 )
  633. WORK( M*M-1 ) = -WORK( N2*M+N2+2 )
  634. T( M, M ) = T( N2+1, N2+1 )
  635. T( M-1, M ) = -T( M, M-1 )
  636. END IF
  637. CALL DGEMM( 'T', 'N', N2, N1, N2, ONE, WORK, M, A( J1, J1+N2 ),
  638. $ LDA, ZERO, WORK( M*M+1 ), N2 )
  639. CALL DLACPY( 'Full', N2, N1, WORK( M*M+1 ), N2, A( J1, J1+N2 ),
  640. $ LDA )
  641. CALL DGEMM( 'T', 'N', N2, N1, N2, ONE, WORK, M, B( J1, J1+N2 ),
  642. $ LDB, ZERO, WORK( M*M+1 ), N2 )
  643. CALL DLACPY( 'Full', N2, N1, WORK( M*M+1 ), N2, B( J1, J1+N2 ),
  644. $ LDB )
  645. CALL DGEMM( 'N', 'N', M, M, M, ONE, LI, LDST, WORK, M, ZERO,
  646. $ WORK( M*M+1 ), M )
  647. CALL DLACPY( 'Full', M, M, WORK( M*M+1 ), M, LI, LDST )
  648. CALL DGEMM( 'N', 'N', N2, N1, N1, ONE, A( J1, J1+N2 ), LDA,
  649. $ T( N2+1, N2+1 ), LDST, ZERO, WORK, N2 )
  650. CALL DLACPY( 'Full', N2, N1, WORK, N2, A( J1, J1+N2 ), LDA )
  651. CALL DGEMM( 'N', 'N', N2, N1, N1, ONE, B( J1, J1+N2 ), LDB,
  652. $ T( N2+1, N2+1 ), LDST, ZERO, WORK, N2 )
  653. CALL DLACPY( 'Full', N2, N1, WORK, N2, B( J1, J1+N2 ), LDB )
  654. CALL DGEMM( 'T', 'N', M, M, M, ONE, IR, LDST, T, LDST, ZERO,
  655. $ WORK, M )
  656. CALL DLACPY( 'Full', M, M, WORK, M, IR, LDST )
  657. *
  658. * Accumulate transformations into Q and Z if requested.
  659. *
  660. IF( WANTQ ) THEN
  661. CALL DGEMM( 'N', 'N', N, M, M, ONE, Q( 1, J1 ), LDQ, LI,
  662. $ LDST, ZERO, WORK, N )
  663. CALL DLACPY( 'Full', N, M, WORK, N, Q( 1, J1 ), LDQ )
  664. *
  665. END IF
  666. *
  667. IF( WANTZ ) THEN
  668. CALL DGEMM( 'N', 'N', N, M, M, ONE, Z( 1, J1 ), LDZ, IR,
  669. $ LDST, ZERO, WORK, N )
  670. CALL DLACPY( 'Full', N, M, WORK, N, Z( 1, J1 ), LDZ )
  671. *
  672. END IF
  673. *
  674. * Update (A(J1:J1+M-1, M+J1:N), B(J1:J1+M-1, M+J1:N)) and
  675. * (A(1:J1-1, J1:J1+M), B(1:J1-1, J1:J1+M)).
  676. *
  677. I = J1 + M
  678. IF( I.LE.N ) THEN
  679. CALL DGEMM( 'T', 'N', M, N-I+1, M, ONE, LI, LDST,
  680. $ A( J1, I ), LDA, ZERO, WORK, M )
  681. CALL DLACPY( 'Full', M, N-I+1, WORK, M, A( J1, I ), LDA )
  682. CALL DGEMM( 'T', 'N', M, N-I+1, M, ONE, LI, LDST,
  683. $ B( J1, I ), LDB, ZERO, WORK, M )
  684. CALL DLACPY( 'Full', M, N-I+1, WORK, M, B( J1, I ), LDB )
  685. END IF
  686. I = J1 - 1
  687. IF( I.GT.0 ) THEN
  688. CALL DGEMM( 'N', 'N', I, M, M, ONE, A( 1, J1 ), LDA, IR,
  689. $ LDST, ZERO, WORK, I )
  690. CALL DLACPY( 'Full', I, M, WORK, I, A( 1, J1 ), LDA )
  691. CALL DGEMM( 'N', 'N', I, M, M, ONE, B( 1, J1 ), LDB, IR,
  692. $ LDST, ZERO, WORK, I )
  693. CALL DLACPY( 'Full', I, M, WORK, I, B( 1, J1 ), LDB )
  694. END IF
  695. *
  696. * Exit with INFO = 0 if swap was successfully performed.
  697. *
  698. RETURN
  699. *
  700. END IF
  701. *
  702. * Exit with INFO = 1 if swap was rejected.
  703. *
  704. 70 CONTINUE
  705. *
  706. INFO = 1
  707. RETURN
  708. *
  709. * End of DTGEX2
  710. *
  711. END