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.

zlaqz0.f 23 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692
  1. *> \brief \b ZLAQZ0
  2. *
  3. * =========== DOCUMENTATION ===========
  4. *
  5. * Online html documentation available at
  6. * http://www.netlib.org/lapack/explore-html/
  7. *
  8. *> \htmlonly
  9. *> Download ZLAQZ0 + dependencies
  10. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zlaqz0.f">
  11. *> [TGZ]</a>
  12. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zlaqz0.f">
  13. *> [ZIP]</a>
  14. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zlaqz0.f">
  15. *> [TXT]</a>
  16. *> \endhtmlonly
  17. *
  18. * Definition:
  19. * ===========
  20. *
  21. * SUBROUTINE ZLAQZ0( WANTS, WANTQ, WANTZ, N, ILO, IHI, A, LDA, B,
  22. * $ LDB, ALPHA, BETA, Q, LDQ, Z, LDZ, WORK, LWORK, RWORK, REC,
  23. * $ INFO )
  24. * IMPLICIT NONE
  25. *
  26. * Arguments
  27. * CHARACTER, INTENT( IN ) :: WANTS, WANTQ, WANTZ
  28. * INTEGER, INTENT( IN ) :: N, ILO, IHI, LDA, LDB, LDQ, LDZ, LWORK,
  29. * $ REC
  30. * INTEGER, INTENT( OUT ) :: INFO
  31. * COMPLEX*16, INTENT( INOUT ) :: A( LDA, * ), B( LDB, * ), Q( LDQ,
  32. * $ * ), Z( LDZ, * ), ALPHA( * ), BETA( * ), WORK( * )
  33. * DOUBLE PRECISION, INTENT( OUT ) :: RWORK( * )
  34. * ..
  35. *
  36. *
  37. *> \par Purpose:
  38. * =============
  39. *>
  40. *> \verbatim
  41. *>
  42. *> ZLAQZ0 computes the eigenvalues of a real matrix pair (H,T),
  43. *> where H is an upper Hessenberg matrix and T is upper triangular,
  44. *> using the double-shift QZ method.
  45. *> Matrix pairs of this type are produced by the reduction to
  46. *> generalized upper Hessenberg form of a real matrix pair (A,B):
  47. *>
  48. *> A = Q1*H*Z1**H, B = Q1*T*Z1**H,
  49. *>
  50. *> as computed by ZGGHRD.
  51. *>
  52. *> If JOB='S', then the Hessenberg-triangular pair (H,T) is
  53. *> also reduced to generalized Schur form,
  54. *>
  55. *> H = Q*S*Z**H, T = Q*P*Z**H,
  56. *>
  57. *> where Q and Z are unitary matrices, P and S are an upper triangular
  58. *> matrices.
  59. *>
  60. *> Optionally, the unitary matrix Q from the generalized Schur
  61. *> factorization may be postmultiplied into an input matrix Q1, and the
  62. *> unitary matrix Z may be postmultiplied into an input matrix Z1.
  63. *> If Q1 and Z1 are the unitary matrices from ZGGHRD that reduced
  64. *> the matrix pair (A,B) to generalized upper Hessenberg form, then the
  65. *> output matrices Q1*Q and Z1*Z are the unitary factors from the
  66. *> generalized Schur factorization of (A,B):
  67. *>
  68. *> A = (Q1*Q)*S*(Z1*Z)**H, B = (Q1*Q)*P*(Z1*Z)**H.
  69. *>
  70. *> To avoid overflow, eigenvalues of the matrix pair (H,T) (equivalently,
  71. *> of (A,B)) are computed as a pair of values (alpha,beta), where alpha is
  72. *> complex and beta real.
  73. *> If beta is nonzero, lambda = alpha / beta is an eigenvalue of the
  74. *> generalized nonsymmetric eigenvalue problem (GNEP)
  75. *> A*x = lambda*B*x
  76. *> and if alpha is nonzero, mu = beta / alpha is an eigenvalue of the
  77. *> alternate form of the GNEP
  78. *> mu*A*y = B*y.
  79. *> Eigenvalues can be read directly from the generalized Schur
  80. *> form:
  81. *> alpha = S(i,i), beta = P(i,i).
  82. *>
  83. *> Ref: C.B. Moler & G.W. Stewart, "An Algorithm for Generalized Matrix
  84. *> Eigenvalue Problems", SIAM J. Numer. Anal., 10(1973),
  85. *> pp. 241--256.
  86. *>
  87. *> Ref: B. Kagstrom, D. Kressner, "Multishift Variants of the QZ
  88. *> Algorithm with Aggressive Early Deflation", SIAM J. Numer.
  89. *> Anal., 29(2006), pp. 199--227.
  90. *>
  91. *> Ref: T. Steel, D. Camps, K. Meerbergen, R. Vandebril "A multishift,
  92. *> multipole rational QZ method with agressive early deflation"
  93. *> \endverbatim
  94. *
  95. * Arguments:
  96. * ==========
  97. *
  98. *> \param[in] WANTS
  99. *> \verbatim
  100. *> WANTS is CHARACTER*1
  101. *> = 'E': Compute eigenvalues only;
  102. *> = 'S': Compute eigenvalues and the Schur form.
  103. *> \endverbatim
  104. *>
  105. *> \param[in] WANTQ
  106. *> \verbatim
  107. *> WANTQ is CHARACTER*1
  108. *> = 'N': Left Schur vectors (Q) are not computed;
  109. *> = 'I': Q is initialized to the unit matrix and the matrix Q
  110. *> of left Schur vectors of (A,B) is returned;
  111. *> = 'V': Q must contain an unitary matrix Q1 on entry and
  112. *> the product Q1*Q is returned.
  113. *> \endverbatim
  114. *>
  115. *> \param[in] WANTZ
  116. *> \verbatim
  117. *> WANTZ is CHARACTER*1
  118. *> = 'N': Right Schur vectors (Z) are not computed;
  119. *> = 'I': Z is initialized to the unit matrix and the matrix Z
  120. *> of right Schur vectors of (A,B) is returned;
  121. *> = 'V': Z must contain an unitary matrix Z1 on entry and
  122. *> the product Z1*Z is returned.
  123. *> \endverbatim
  124. *>
  125. *> \param[in] N
  126. *> \verbatim
  127. *> N is INTEGER
  128. *> The order of the matrices A, B, Q, and Z. N >= 0.
  129. *> \endverbatim
  130. *>
  131. *> \param[in] ILO
  132. *> \verbatim
  133. *> ILO is INTEGER
  134. *> \endverbatim
  135. *>
  136. *> \param[in] IHI
  137. *> \verbatim
  138. *> IHI is INTEGER
  139. *> ILO and IHI mark the rows and columns of A which are in
  140. *> Hessenberg form. It is assumed that A is already upper
  141. *> triangular in rows and columns 1:ILO-1 and IHI+1:N.
  142. *> If N > 0, 1 <= ILO <= IHI <= N; if N = 0, ILO=1 and IHI=0.
  143. *> \endverbatim
  144. *>
  145. *> \param[in,out] A
  146. *> \verbatim
  147. *> A is COMPLEX*16 array, dimension (LDA, N)
  148. *> On entry, the N-by-N upper Hessenberg matrix A.
  149. *> On exit, if JOB = 'S', A contains the upper triangular
  150. *> matrix S from the generalized Schur factorization.
  151. *> If JOB = 'E', the diagonal blocks of A match those of S, but
  152. *> the rest of A is unspecified.
  153. *> \endverbatim
  154. *>
  155. *> \param[in] LDA
  156. *> \verbatim
  157. *> LDA is INTEGER
  158. *> The leading dimension of the array A. LDA >= max( 1, N ).
  159. *> \endverbatim
  160. *>
  161. *> \param[in,out] B
  162. *> \verbatim
  163. *> B is COMPLEX*16 array, dimension (LDB, N)
  164. *> On entry, the N-by-N upper triangular matrix B.
  165. *> On exit, if JOB = 'S', B contains the upper triangular
  166. *> matrix P from the generalized Schur factorization;
  167. *> If JOB = 'E', the diagonal blocks of B match those of P, but
  168. *> the rest of B is unspecified.
  169. *> \endverbatim
  170. *>
  171. *> \param[in] LDB
  172. *> \verbatim
  173. *> LDB is INTEGER
  174. *> The leading dimension of the array B. LDB >= max( 1, N ).
  175. *> \endverbatim
  176. *>
  177. *> \param[out] ALPHA
  178. *> \verbatim
  179. *> ALPHA is COMPLEX*16 array, dimension (N)
  180. *> Each scalar alpha defining an eigenvalue
  181. *> of GNEP.
  182. *> \endverbatim
  183. *>
  184. *> \param[out] BETA
  185. *> \verbatim
  186. *> BETA is COMPLEX*16 array, dimension (N)
  187. *> The scalars beta that define the eigenvalues of GNEP.
  188. *> Together, the quantities alpha = ALPHA(j) and
  189. *> beta = BETA(j) represent the j-th eigenvalue of the matrix
  190. *> pair (A,B), in one of the forms lambda = alpha/beta or
  191. *> mu = beta/alpha. Since either lambda or mu may overflow,
  192. *> they should not, in general, be computed.
  193. *> \endverbatim
  194. *>
  195. *> \param[in,out] Q
  196. *> \verbatim
  197. *> Q is COMPLEX*16 array, dimension (LDQ, N)
  198. *> On entry, if COMPQ = 'V', the unitary matrix Q1 used in
  199. *> the reduction of (A,B) to generalized Hessenberg form.
  200. *> On exit, if COMPQ = 'I', the unitary matrix of left Schur
  201. *> vectors of (A,B), and if COMPQ = 'V', the unitary matrix
  202. *> of left Schur vectors of (A,B).
  203. *> Not referenced if COMPQ = 'N'.
  204. *> \endverbatim
  205. *>
  206. *> \param[in] LDQ
  207. *> \verbatim
  208. *> LDQ is INTEGER
  209. *> The leading dimension of the array Q. LDQ >= 1.
  210. *> If COMPQ='V' or 'I', then LDQ >= N.
  211. *> \endverbatim
  212. *>
  213. *> \param[in,out] Z
  214. *> \verbatim
  215. *> Z is COMPLEX*16 array, dimension (LDZ, N)
  216. *> On entry, if COMPZ = 'V', the unitary matrix Z1 used in
  217. *> the reduction of (A,B) to generalized Hessenberg form.
  218. *> On exit, if COMPZ = 'I', the unitary matrix of
  219. *> right Schur vectors of (H,T), and if COMPZ = 'V', the
  220. *> unitary matrix of right Schur vectors of (A,B).
  221. *> Not referenced if COMPZ = 'N'.
  222. *> \endverbatim
  223. *>
  224. *> \param[in] LDZ
  225. *> \verbatim
  226. *> LDZ is INTEGER
  227. *> The leading dimension of the array Z. LDZ >= 1.
  228. *> If COMPZ='V' or 'I', then LDZ >= N.
  229. *> \endverbatim
  230. *>
  231. *> \param[out] WORK
  232. *> \verbatim
  233. *> WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
  234. *> On exit, if INFO >= 0, WORK(1) returns the optimal LWORK.
  235. *> \endverbatim
  236. *>
  237. *> \param[out] RWORK
  238. *> \verbatim
  239. *> RWORK is DOUBLE PRECISION array, dimension (N)
  240. *> \endverbatim
  241. *>
  242. *> \param[in] LWORK
  243. *> \verbatim
  244. *> LWORK is INTEGER
  245. *> The dimension of the array WORK. LWORK >= max(1,N).
  246. *>
  247. *> If LWORK = -1, then a workspace query is assumed; the routine
  248. *> only calculates the optimal size of the WORK array, returns
  249. *> this value as the first entry of the WORK array, and no error
  250. *> message related to LWORK is issued by XERBLA.
  251. *> \endverbatim
  252. *>
  253. *> \param[in] REC
  254. *> \verbatim
  255. *> REC is INTEGER
  256. *> REC indicates the current recursion level. Should be set
  257. *> to 0 on first call.
  258. *> \endverbatim
  259. *>
  260. *> \param[out] INFO
  261. *> \verbatim
  262. *> INFO is INTEGER
  263. *> = 0: successful exit
  264. *> < 0: if INFO = -i, the i-th argument had an illegal value
  265. *> = 1,...,N: the QZ iteration did not converge. (A,B) is not
  266. *> in Schur form, but ALPHA(i) and
  267. *> BETA(i), i=INFO+1,...,N should be correct.
  268. *> \endverbatim
  269. *
  270. * Authors:
  271. * ========
  272. *
  273. *> \author Thijs Steel, KU Leuven
  274. *
  275. *> \date May 2020
  276. *
  277. *> \ingroup complex16GEcomputational
  278. *>
  279. * =====================================================================
  280. RECURSIVE SUBROUTINE ZLAQZ0( WANTS, WANTQ, WANTZ, N, ILO, IHI, A,
  281. $ LDA, B, LDB, ALPHA, BETA, Q, LDQ, Z,
  282. $ LDZ, WORK, LWORK, RWORK, REC,
  283. $ INFO )
  284. IMPLICIT NONE
  285. * Arguments
  286. CHARACTER, INTENT( IN ) :: WANTS, WANTQ, WANTZ
  287. INTEGER, INTENT( IN ) :: N, ILO, IHI, LDA, LDB, LDQ, LDZ, LWORK,
  288. $ REC
  289. INTEGER, INTENT( OUT ) :: INFO
  290. COMPLEX*16, INTENT( INOUT ) :: A( LDA, * ), B( LDB, * ), Q( LDQ,
  291. $ * ), Z( LDZ, * ), ALPHA( * ), BETA( * ), WORK( * )
  292. DOUBLE PRECISION, INTENT( OUT ) :: RWORK( * )
  293. * Parameters
  294. COMPLEX*16 CZERO, CONE
  295. PARAMETER ( CZERO = ( 0.0D+0, 0.0D+0 ), CONE = ( 1.0D+0,
  296. $ 0.0D+0 ) )
  297. DOUBLE PRECISION :: ZERO, ONE, HALF
  298. PARAMETER( ZERO = 0.0D0, ONE = 1.0D0, HALF = 0.5D0 )
  299. * Local scalars
  300. DOUBLE PRECISION :: SMLNUM, ULP, SAFMIN, SAFMAX, C1, TEMPR
  301. COMPLEX*16 :: ESHIFT, S1, TEMP
  302. INTEGER :: ISTART, ISTOP, IITER, MAXIT, ISTART2, K, LD, NSHIFTS,
  303. $ NBLOCK, NW, NMIN, NIBBLE, N_UNDEFLATED, N_DEFLATED,
  304. $ NS, SWEEP_INFO, SHIFTPOS, LWORKREQ, K2, ISTARTM,
  305. $ ISTOPM, IWANTS, IWANTQ, IWANTZ, NORM_INFO, AED_INFO,
  306. $ NWR, NBR, NSR, ITEMP1, ITEMP2, RCOST
  307. LOGICAL :: ILSCHUR, ILQ, ILZ
  308. CHARACTER :: JBCMPZ*3
  309. * External Functions
  310. EXTERNAL :: XERBLA, ZHGEQZ, ZLAQZ2, ZLAQZ3, ZLASET, DLABAD,
  311. $ ZLARTG, ZROT
  312. DOUBLE PRECISION, EXTERNAL :: DLAMCH
  313. LOGICAL, EXTERNAL :: LSAME
  314. INTEGER, EXTERNAL :: ILAENV
  315. *
  316. * Decode wantS,wantQ,wantZ
  317. *
  318. IF( LSAME( WANTS, 'E' ) ) THEN
  319. ILSCHUR = .FALSE.
  320. IWANTS = 1
  321. ELSE IF( LSAME( WANTS, 'S' ) ) THEN
  322. ILSCHUR = .TRUE.
  323. IWANTS = 2
  324. ELSE
  325. IWANTS = 0
  326. END IF
  327. IF( LSAME( WANTQ, 'N' ) ) THEN
  328. ILQ = .FALSE.
  329. IWANTQ = 1
  330. ELSE IF( LSAME( WANTQ, 'V' ) ) THEN
  331. ILQ = .TRUE.
  332. IWANTQ = 2
  333. ELSE IF( LSAME( WANTQ, 'I' ) ) THEN
  334. ILQ = .TRUE.
  335. IWANTQ = 3
  336. ELSE
  337. IWANTQ = 0
  338. END IF
  339. IF( LSAME( WANTZ, 'N' ) ) THEN
  340. ILZ = .FALSE.
  341. IWANTZ = 1
  342. ELSE IF( LSAME( WANTZ, 'V' ) ) THEN
  343. ILZ = .TRUE.
  344. IWANTZ = 2
  345. ELSE IF( LSAME( WANTZ, 'I' ) ) THEN
  346. ILZ = .TRUE.
  347. IWANTZ = 3
  348. ELSE
  349. IWANTZ = 0
  350. END IF
  351. *
  352. * Check Argument Values
  353. *
  354. INFO = 0
  355. IF( IWANTS.EQ.0 ) THEN
  356. INFO = -1
  357. ELSE IF( IWANTQ.EQ.0 ) THEN
  358. INFO = -2
  359. ELSE IF( IWANTZ.EQ.0 ) THEN
  360. INFO = -3
  361. ELSE IF( N.LT.0 ) THEN
  362. INFO = -4
  363. ELSE IF( ILO.LT.1 ) THEN
  364. INFO = -5
  365. ELSE IF( IHI.GT.N .OR. IHI.LT.ILO-1 ) THEN
  366. INFO = -6
  367. ELSE IF( LDA.LT.N ) THEN
  368. INFO = -8
  369. ELSE IF( LDB.LT.N ) THEN
  370. INFO = -10
  371. ELSE IF( LDQ.LT.1 .OR. ( ILQ .AND. LDQ.LT.N ) ) THEN
  372. INFO = -15
  373. ELSE IF( LDZ.LT.1 .OR. ( ILZ .AND. LDZ.LT.N ) ) THEN
  374. INFO = -17
  375. END IF
  376. IF( INFO.NE.0 ) THEN
  377. CALL XERBLA( 'ZLAQZ0', -INFO )
  378. RETURN
  379. END IF
  380. *
  381. * Quick return if possible
  382. *
  383. IF( N.LE.0 ) THEN
  384. WORK( 1 ) = DBLE( 1 )
  385. RETURN
  386. END IF
  387. *
  388. * Get the parameters
  389. *
  390. JBCMPZ( 1:1 ) = WANTS
  391. JBCMPZ( 2:2 ) = WANTQ
  392. JBCMPZ( 3:3 ) = WANTZ
  393. NMIN = ILAENV( 12, 'ZLAQZ0', JBCMPZ, N, ILO, IHI, LWORK )
  394. NWR = ILAENV( 13, 'ZLAQZ0', JBCMPZ, N, ILO, IHI, LWORK )
  395. NWR = MAX( 2, NWR )
  396. NWR = MIN( IHI-ILO+1, ( N-1 ) / 3, NWR )
  397. NIBBLE = ILAENV( 14, 'ZLAQZ0', JBCMPZ, N, ILO, IHI, LWORK )
  398. NSR = ILAENV( 15, 'ZLAQZ0', JBCMPZ, N, ILO, IHI, LWORK )
  399. NSR = MIN( NSR, ( N+6 ) / 9, IHI-ILO )
  400. NSR = MAX( 2, NSR-MOD( NSR, 2 ) )
  401. RCOST = ILAENV( 17, 'ZLAQZ0', JBCMPZ, N, ILO, IHI, LWORK )
  402. ITEMP1 = INT( NSR/SQRT( 1+2*NSR/( DBLE( RCOST )/100*N ) ) )
  403. ITEMP1 = ( ( ITEMP1-1 )/4 )*4+4
  404. NBR = NSR+ITEMP1
  405. IF( N .LT. NMIN .OR. REC .GE. 2 ) THEN
  406. CALL ZHGEQZ( WANTS, WANTQ, WANTZ, N, ILO, IHI, A, LDA, B, LDB,
  407. $ ALPHA, BETA, Q, LDQ, Z, LDZ, WORK, LWORK, RWORK,
  408. $ INFO )
  409. RETURN
  410. END IF
  411. *
  412. * Find out required workspace
  413. *
  414. * Workspace query to ZLAQZ2
  415. NW = MAX( NWR, NMIN )
  416. CALL ZLAQZ2( ILSCHUR, ILQ, ILZ, N, ILO, IHI, NW, A, LDA, B, LDB,
  417. $ Q, LDQ, Z, LDZ, N_UNDEFLATED, N_DEFLATED, ALPHA,
  418. $ BETA, WORK, NW, WORK, NW, WORK, -1, RWORK, REC,
  419. $ AED_INFO )
  420. ITEMP1 = INT( WORK( 1 ) )
  421. * Workspace query to ZLAQZ3
  422. CALL ZLAQZ3( ILSCHUR, ILQ, ILZ, N, ILO, IHI, NSR, NBR, ALPHA,
  423. $ BETA, A, LDA, B, LDB, Q, LDQ, Z, LDZ, WORK, NBR,
  424. $ WORK, NBR, WORK, -1, SWEEP_INFO )
  425. ITEMP2 = INT( WORK( 1 ) )
  426. LWORKREQ = MAX( ITEMP1+2*NW**2, ITEMP2+2*NBR**2 )
  427. IF ( LWORK .EQ.-1 ) THEN
  428. WORK( 1 ) = DBLE( LWORKREQ )
  429. RETURN
  430. ELSE IF ( LWORK .LT. LWORKREQ ) THEN
  431. INFO = -19
  432. END IF
  433. IF( INFO.NE.0 ) THEN
  434. CALL XERBLA( 'ZLAQZ0', INFO )
  435. RETURN
  436. END IF
  437. *
  438. * Initialize Q and Z
  439. *
  440. IF( IWANTQ.EQ.3 ) CALL ZLASET( 'FULL', N, N, CZERO, CONE, Q,
  441. $ LDQ )
  442. IF( IWANTZ.EQ.3 ) CALL ZLASET( 'FULL', N, N, CZERO, CONE, Z,
  443. $ LDZ )
  444. * Get machine constants
  445. SAFMIN = DLAMCH( 'SAFE MINIMUM' )
  446. SAFMAX = ONE/SAFMIN
  447. CALL DLABAD( SAFMIN, SAFMAX )
  448. ULP = DLAMCH( 'PRECISION' )
  449. SMLNUM = SAFMIN*( DBLE( N )/ULP )
  450. ISTART = ILO
  451. ISTOP = IHI
  452. MAXIT = 30*( IHI-ILO+1 )
  453. LD = 0
  454. DO IITER = 1, MAXIT
  455. IF( IITER .GE. MAXIT ) THEN
  456. INFO = ISTOP+1
  457. GOTO 80
  458. END IF
  459. IF ( ISTART+1 .GE. ISTOP ) THEN
  460. ISTOP = ISTART
  461. EXIT
  462. END IF
  463. * Check deflations at the end
  464. IF ( ABS( A( ISTOP, ISTOP-1 ) ) .LE. MAX( SMLNUM,
  465. $ ULP*( ABS( A( ISTOP, ISTOP ) )+ABS( A( ISTOP-1,
  466. $ ISTOP-1 ) ) ) ) ) THEN
  467. A( ISTOP, ISTOP-1 ) = CZERO
  468. ISTOP = ISTOP-1
  469. LD = 0
  470. ESHIFT = CZERO
  471. END IF
  472. * Check deflations at the start
  473. IF ( ABS( A( ISTART+1, ISTART ) ) .LE. MAX( SMLNUM,
  474. $ ULP*( ABS( A( ISTART, ISTART ) )+ABS( A( ISTART+1,
  475. $ ISTART+1 ) ) ) ) ) THEN
  476. A( ISTART+1, ISTART ) = CZERO
  477. ISTART = ISTART+1
  478. LD = 0
  479. ESHIFT = CZERO
  480. END IF
  481. IF ( ISTART+1 .GE. ISTOP ) THEN
  482. EXIT
  483. END IF
  484. * Check interior deflations
  485. ISTART2 = ISTART
  486. DO K = ISTOP, ISTART+1, -1
  487. IF ( ABS( A( K, K-1 ) ) .LE. MAX( SMLNUM, ULP*( ABS( A( K,
  488. $ K ) )+ABS( A( K-1, K-1 ) ) ) ) ) THEN
  489. A( K, K-1 ) = CZERO
  490. ISTART2 = K
  491. EXIT
  492. END IF
  493. END DO
  494. * Get range to apply rotations to
  495. IF ( ILSCHUR ) THEN
  496. ISTARTM = 1
  497. ISTOPM = N
  498. ELSE
  499. ISTARTM = ISTART2
  500. ISTOPM = ISTOP
  501. END IF
  502. * Check infinite eigenvalues, this is done without blocking so might
  503. * slow down the method when many infinite eigenvalues are present
  504. K = ISTOP
  505. DO WHILE ( K.GE.ISTART2 )
  506. TEMPR = ZERO
  507. IF( K .LT. ISTOP ) THEN
  508. TEMPR = TEMPR+ABS( B( K, K+1 ) )
  509. END IF
  510. IF( K .GT. ISTART2 ) THEN
  511. TEMPR = TEMPR+ABS( B( K-1, K ) )
  512. END IF
  513. IF( ABS( B( K, K ) ) .LT. MAX( SMLNUM, ULP*TEMPR ) ) THEN
  514. * A diagonal element of B is negligable, move it
  515. * to the top and deflate it
  516. DO K2 = K, ISTART2+1, -1
  517. CALL ZLARTG( B( K2-1, K2 ), B( K2-1, K2-1 ), C1, S1,
  518. $ TEMP )
  519. B( K2-1, K2 ) = TEMP
  520. B( K2-1, K2-1 ) = CZERO
  521. CALL ZROT( K2-2-ISTARTM+1, B( ISTARTM, K2 ), 1,
  522. $ B( ISTARTM, K2-1 ), 1, C1, S1 )
  523. CALL ZROT( MIN( K2+1, ISTOP )-ISTARTM+1, A( ISTARTM,
  524. $ K2 ), 1, A( ISTARTM, K2-1 ), 1, C1, S1 )
  525. IF ( ILZ ) THEN
  526. CALL ZROT( N, Z( 1, K2 ), 1, Z( 1, K2-1 ), 1, C1,
  527. $ S1 )
  528. END IF
  529. IF( K2.LT.ISTOP ) THEN
  530. CALL ZLARTG( A( K2, K2-1 ), A( K2+1, K2-1 ), C1,
  531. $ S1, TEMP )
  532. A( K2, K2-1 ) = TEMP
  533. A( K2+1, K2-1 ) = CZERO
  534. CALL ZROT( ISTOPM-K2+1, A( K2, K2 ), LDA, A( K2+1,
  535. $ K2 ), LDA, C1, S1 )
  536. CALL ZROT( ISTOPM-K2+1, B( K2, K2 ), LDB, B( K2+1,
  537. $ K2 ), LDB, C1, S1 )
  538. IF( ILQ ) THEN
  539. CALL ZROT( N, Q( 1, K2 ), 1, Q( 1, K2+1 ), 1,
  540. $ C1, DCONJG( S1 ) )
  541. END IF
  542. END IF
  543. END DO
  544. IF( ISTART2.LT.ISTOP )THEN
  545. CALL ZLARTG( A( ISTART2, ISTART2 ), A( ISTART2+1,
  546. $ ISTART2 ), C1, S1, TEMP )
  547. A( ISTART2, ISTART2 ) = TEMP
  548. A( ISTART2+1, ISTART2 ) = CZERO
  549. CALL ZROT( ISTOPM-( ISTART2+1 )+1, A( ISTART2,
  550. $ ISTART2+1 ), LDA, A( ISTART2+1,
  551. $ ISTART2+1 ), LDA, C1, S1 )
  552. CALL ZROT( ISTOPM-( ISTART2+1 )+1, B( ISTART2,
  553. $ ISTART2+1 ), LDB, B( ISTART2+1,
  554. $ ISTART2+1 ), LDB, C1, S1 )
  555. IF( ILQ ) THEN
  556. CALL ZROT( N, Q( 1, ISTART2 ), 1, Q( 1,
  557. $ ISTART2+1 ), 1, C1, DCONJG( S1 ) )
  558. END IF
  559. END IF
  560. ISTART2 = ISTART2+1
  561. END IF
  562. K = K-1
  563. END DO
  564. * istart2 now points to the top of the bottom right
  565. * unreduced Hessenberg block
  566. IF ( ISTART2 .GE. ISTOP ) THEN
  567. ISTOP = ISTART2-1
  568. LD = 0
  569. ESHIFT = CZERO
  570. CYCLE
  571. END IF
  572. NW = NWR
  573. NSHIFTS = NSR
  574. NBLOCK = NBR
  575. IF ( ISTOP-ISTART2+1 .LT. NMIN ) THEN
  576. * Setting nw to the size of the subblock will make AED deflate
  577. * all the eigenvalues. This is slightly more efficient than just
  578. * using qz_small because the off diagonal part gets updated via BLAS.
  579. IF ( ISTOP-ISTART+1 .LT. NMIN ) THEN
  580. NW = ISTOP-ISTART+1
  581. ISTART2 = ISTART
  582. ELSE
  583. NW = ISTOP-ISTART2+1
  584. END IF
  585. END IF
  586. *
  587. * Time for AED
  588. *
  589. CALL ZLAQZ2( ILSCHUR, ILQ, ILZ, N, ISTART2, ISTOP, NW, A, LDA,
  590. $ B, LDB, Q, LDQ, Z, LDZ, N_UNDEFLATED, N_DEFLATED,
  591. $ ALPHA, BETA, WORK, NW, WORK( NW**2+1 ), NW,
  592. $ WORK( 2*NW**2+1 ), LWORK-2*NW**2, RWORK, REC,
  593. $ AED_INFO )
  594. IF ( N_DEFLATED > 0 ) THEN
  595. ISTOP = ISTOP-N_DEFLATED
  596. LD = 0
  597. ESHIFT = CZERO
  598. END IF
  599. IF ( 100*N_DEFLATED > NIBBLE*( N_DEFLATED+N_UNDEFLATED ) .OR.
  600. $ ISTOP-ISTART2+1 .LT. NMIN ) THEN
  601. * AED has uncovered many eigenvalues. Skip a QZ sweep and run
  602. * AED again.
  603. CYCLE
  604. END IF
  605. LD = LD+1
  606. NS = MIN( NSHIFTS, ISTOP-ISTART2 )
  607. NS = MIN( NS, N_UNDEFLATED )
  608. SHIFTPOS = ISTOP-N_DEFLATED-N_UNDEFLATED+1
  609. IF ( MOD( LD, 6 ) .EQ. 0 ) THEN
  610. *
  611. * Exceptional shift. Chosen for no particularly good reason.
  612. *
  613. IF( ( DBLE( MAXIT )*SAFMIN )*ABS( A( ISTOP,
  614. $ ISTOP-1 ) ).LT.ABS( A( ISTOP-1, ISTOP-1 ) ) ) THEN
  615. ESHIFT = A( ISTOP, ISTOP-1 )/B( ISTOP-1, ISTOP-1 )
  616. ELSE
  617. ESHIFT = ESHIFT+CONE/( SAFMIN*DBLE( MAXIT ) )
  618. END IF
  619. ALPHA( SHIFTPOS ) = CONE
  620. BETA( SHIFTPOS ) = ESHIFT
  621. NS = 1
  622. END IF
  623. *
  624. * Time for a QZ sweep
  625. *
  626. CALL ZLAQZ3( ILSCHUR, ILQ, ILZ, N, ISTART2, ISTOP, NS, NBLOCK,
  627. $ ALPHA( SHIFTPOS ), BETA( SHIFTPOS ), A, LDA, B,
  628. $ LDB, Q, LDQ, Z, LDZ, WORK, NBLOCK, WORK( NBLOCK**
  629. $ 2+1 ), NBLOCK, WORK( 2*NBLOCK**2+1 ),
  630. $ LWORK-2*NBLOCK**2, SWEEP_INFO )
  631. END DO
  632. *
  633. * Call ZHGEQZ to normalize the eigenvalue blocks and set the eigenvalues
  634. * If all the eigenvalues have been found, ZHGEQZ will not do any iterations
  635. * and only normalize the blocks. In case of a rare convergence failure,
  636. * the single shift might perform better.
  637. *
  638. 80 CALL ZHGEQZ( WANTS, WANTQ, WANTZ, N, ILO, IHI, A, LDA, B, LDB,
  639. $ ALPHA, BETA, Q, LDQ, Z, LDZ, WORK, LWORK, RWORK,
  640. $ NORM_INFO )
  641. INFO = NORM_INFO
  642. END SUBROUTINE