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.

dlaqz0.f 25 kB

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