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.

slaqz0.f 25 kB

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