1*> \brief \b SLAQR0 computes the eigenvalues of a Hessenberg matrix, and optionally the matrices from the Schur decomposition. 2* 3* =========== DOCUMENTATION =========== 4* 5* Online html documentation available at 6* http://www.netlib.org/lapack/explore-html/ 7* 8*> \htmlonly 9*> Download SLAQR0 + dependencies 10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/slaqr0.f"> 11*> [TGZ]</a> 12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/slaqr0.f"> 13*> [ZIP]</a> 14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/slaqr0.f"> 15*> [TXT]</a> 16*> \endhtmlonly 17* 18* Definition: 19* =========== 20* 21* SUBROUTINE SLAQR0( WANTT, WANTZ, N, ILO, IHI, H, LDH, WR, WI, 22* ILOZ, IHIZ, Z, LDZ, WORK, LWORK, INFO ) 23* 24* .. Scalar Arguments .. 25* INTEGER IHI, IHIZ, ILO, ILOZ, INFO, LDH, LDZ, LWORK, N 26* LOGICAL WANTT, WANTZ 27* .. 28* .. Array Arguments .. 29* REAL H( LDH, * ), WI( * ), WORK( * ), WR( * ), 30* $ Z( LDZ, * ) 31* .. 32* 33* 34*> \par Purpose: 35* ============= 36*> 37*> \verbatim 38*> 39*> SLAQR0 computes the eigenvalues of a Hessenberg matrix H 40*> and, optionally, the matrices T and Z from the Schur decomposition 41*> H = Z T Z**T, where T is an upper quasi-triangular matrix (the 42*> Schur form), and Z is the orthogonal matrix of Schur vectors. 43*> 44*> Optionally Z may be postmultiplied into an input orthogonal 45*> matrix Q so that this routine can give the Schur factorization 46*> of a matrix A which has been reduced to the Hessenberg form H 47*> by the orthogonal matrix Q: A = Q*H*Q**T = (QZ)*T*(QZ)**T. 48*> \endverbatim 49* 50* Arguments: 51* ========== 52* 53*> \param[in] WANTT 54*> \verbatim 55*> WANTT is LOGICAL 56*> = .TRUE. : the full Schur form T is required; 57*> = .FALSE.: only eigenvalues are required. 58*> \endverbatim 59*> 60*> \param[in] WANTZ 61*> \verbatim 62*> WANTZ is LOGICAL 63*> = .TRUE. : the matrix of Schur vectors Z is required; 64*> = .FALSE.: Schur vectors are not required. 65*> \endverbatim 66*> 67*> \param[in] N 68*> \verbatim 69*> N is INTEGER 70*> The order of the matrix H. N >= 0. 71*> \endverbatim 72*> 73*> \param[in] ILO 74*> \verbatim 75*> ILO is INTEGER 76*> \endverbatim 77*> 78*> \param[in] IHI 79*> \verbatim 80*> IHI is INTEGER 81*> It is assumed that H is already upper triangular in rows 82*> and columns 1:ILO-1 and IHI+1:N and, if ILO > 1, 83*> H(ILO,ILO-1) is zero. ILO and IHI are normally set by a 84*> previous call to SGEBAL, and then passed to SGEHRD when the 85*> matrix output by SGEBAL is reduced to Hessenberg form. 86*> Otherwise, ILO and IHI should be set to 1 and N, 87*> respectively. If N > 0, then 1 <= ILO <= IHI <= N. 88*> If N = 0, then ILO = 1 and IHI = 0. 89*> \endverbatim 90*> 91*> \param[in,out] H 92*> \verbatim 93*> H is REAL array, dimension (LDH,N) 94*> On entry, the upper Hessenberg matrix H. 95*> On exit, if INFO = 0 and WANTT is .TRUE., then H contains 96*> the upper quasi-triangular matrix T from the Schur 97*> decomposition (the Schur form); 2-by-2 diagonal blocks 98*> (corresponding to complex conjugate pairs of eigenvalues) 99*> are returned in standard form, with H(i,i) = H(i+1,i+1) 100*> and H(i+1,i)*H(i,i+1) < 0. If INFO = 0 and WANTT is 101*> .FALSE., then the contents of H are unspecified on exit. 102*> (The output value of H when INFO > 0 is given under the 103*> description of INFO below.) 104*> 105*> This subroutine may explicitly set H(i,j) = 0 for i > j and 106*> j = 1, 2, ... ILO-1 or j = IHI+1, IHI+2, ... N. 107*> \endverbatim 108*> 109*> \param[in] LDH 110*> \verbatim 111*> LDH is INTEGER 112*> The leading dimension of the array H. LDH >= max(1,N). 113*> \endverbatim 114*> 115*> \param[out] WR 116*> \verbatim 117*> WR is REAL array, dimension (IHI) 118*> \endverbatim 119*> 120*> \param[out] WI 121*> \verbatim 122*> WI is REAL array, dimension (IHI) 123*> The real and imaginary parts, respectively, of the computed 124*> eigenvalues of H(ILO:IHI,ILO:IHI) are stored in WR(ILO:IHI) 125*> and WI(ILO:IHI). If two eigenvalues are computed as a 126*> complex conjugate pair, they are stored in consecutive 127*> elements of WR and WI, say the i-th and (i+1)th, with 128*> WI(i) > 0 and WI(i+1) < 0. If WANTT is .TRUE., then 129*> the eigenvalues are stored in the same order as on the 130*> diagonal of the Schur form returned in H, with 131*> WR(i) = H(i,i) and, if H(i:i+1,i:i+1) is a 2-by-2 diagonal 132*> block, WI(i) = sqrt(-H(i+1,i)*H(i,i+1)) and 133*> WI(i+1) = -WI(i). 134*> \endverbatim 135*> 136*> \param[in] ILOZ 137*> \verbatim 138*> ILOZ is INTEGER 139*> \endverbatim 140*> 141*> \param[in] IHIZ 142*> \verbatim 143*> IHIZ is INTEGER 144*> Specify the rows of Z to which transformations must be 145*> applied if WANTZ is .TRUE.. 146*> 1 <= ILOZ <= ILO; IHI <= IHIZ <= N. 147*> \endverbatim 148*> 149*> \param[in,out] Z 150*> \verbatim 151*> Z is REAL array, dimension (LDZ,IHI) 152*> If WANTZ is .FALSE., then Z is not referenced. 153*> If WANTZ is .TRUE., then Z(ILO:IHI,ILOZ:IHIZ) is 154*> replaced by Z(ILO:IHI,ILOZ:IHIZ)*U where U is the 155*> orthogonal Schur factor of H(ILO:IHI,ILO:IHI). 156*> (The output value of Z when INFO > 0 is given under 157*> the description of INFO below.) 158*> \endverbatim 159*> 160*> \param[in] LDZ 161*> \verbatim 162*> LDZ is INTEGER 163*> The leading dimension of the array Z. if WANTZ is .TRUE. 164*> then LDZ >= MAX(1,IHIZ). Otherwise, LDZ >= 1. 165*> \endverbatim 166*> 167*> \param[out] WORK 168*> \verbatim 169*> WORK is REAL array, dimension LWORK 170*> On exit, if LWORK = -1, WORK(1) returns an estimate of 171*> the optimal value for LWORK. 172*> \endverbatim 173*> 174*> \param[in] LWORK 175*> \verbatim 176*> LWORK is INTEGER 177*> The dimension of the array WORK. LWORK >= max(1,N) 178*> is sufficient, but LWORK typically as large as 6*N may 179*> be required for optimal performance. A workspace query 180*> to determine the optimal workspace size is recommended. 181*> 182*> If LWORK = -1, then SLAQR0 does a workspace query. 183*> In this case, SLAQR0 checks the input parameters and 184*> estimates the optimal workspace size for the given 185*> values of N, ILO and IHI. The estimate is returned 186*> in WORK(1). No error message related to LWORK is 187*> issued by XERBLA. Neither H nor Z are accessed. 188*> \endverbatim 189*> 190*> \param[out] INFO 191*> \verbatim 192*> INFO is INTEGER 193*> = 0: successful exit 194*> > 0: if INFO = i, SLAQR0 failed to compute all of 195*> the eigenvalues. Elements 1:ilo-1 and i+1:n of WR 196*> and WI contain those eigenvalues which have been 197*> successfully computed. (Failures are rare.) 198*> 199*> If INFO > 0 and WANT is .FALSE., then on exit, 200*> the remaining unconverged eigenvalues are the eigen- 201*> values of the upper Hessenberg matrix rows and 202*> columns ILO through INFO of the final, output 203*> value of H. 204*> 205*> If INFO > 0 and WANTT is .TRUE., then on exit 206*> 207*> (*) (initial value of H)*U = U*(final value of H) 208*> 209*> where U is an orthogonal matrix. The final 210*> value of H is upper Hessenberg and quasi-triangular 211*> in rows and columns INFO+1 through IHI. 212*> 213*> If INFO > 0 and WANTZ is .TRUE., then on exit 214*> 215*> (final value of Z(ILO:IHI,ILOZ:IHIZ) 216*> = (initial value of Z(ILO:IHI,ILOZ:IHIZ)*U 217*> 218*> where U is the orthogonal matrix in (*) (regard- 219*> less of the value of WANTT.) 220*> 221*> If INFO > 0 and WANTZ is .FALSE., then Z is not 222*> accessed. 223*> \endverbatim 224* 225* Authors: 226* ======== 227* 228*> \author Univ. of Tennessee 229*> \author Univ. of California Berkeley 230*> \author Univ. of Colorado Denver 231*> \author NAG Ltd. 232* 233*> \ingroup realOTHERauxiliary 234* 235*> \par Contributors: 236* ================== 237*> 238*> Karen Braman and Ralph Byers, Department of Mathematics, 239*> University of Kansas, USA 240* 241*> \par References: 242* ================ 243*> 244*> K. Braman, R. Byers and R. Mathias, The Multi-Shift QR 245*> Algorithm Part I: Maintaining Well Focused Shifts, and Level 3 246*> Performance, SIAM Journal of Matrix Analysis, volume 23, pages 247*> 929--947, 2002. 248*> \n 249*> K. Braman, R. Byers and R. Mathias, The Multi-Shift QR 250*> Algorithm Part II: Aggressive Early Deflation, SIAM Journal 251*> of Matrix Analysis, volume 23, pages 948--973, 2002. 252*> 253* ===================================================================== 254 SUBROUTINE SLAQR0( WANTT, WANTZ, N, ILO, IHI, H, LDH, WR, WI, 255 $ ILOZ, IHIZ, Z, LDZ, WORK, LWORK, INFO ) 256* 257* -- LAPACK auxiliary routine -- 258* -- LAPACK is a software package provided by Univ. of Tennessee, -- 259* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- 260* 261* .. Scalar Arguments .. 262 INTEGER IHI, IHIZ, ILO, ILOZ, INFO, LDH, LDZ, LWORK, N 263 LOGICAL WANTT, WANTZ 264* .. 265* .. Array Arguments .. 266 REAL H( LDH, * ), WI( * ), WORK( * ), WR( * ), 267 $ Z( LDZ, * ) 268* .. 269* 270* ================================================================ 271* .. Parameters .. 272* 273* ==== Matrices of order NTINY or smaller must be processed by 274* . SLAHQR because of insufficient subdiagonal scratch space. 275* . (This is a hard limit.) ==== 276 INTEGER NTINY 277 PARAMETER ( NTINY = 15 ) 278* 279* ==== Exceptional deflation windows: try to cure rare 280* . slow convergence by varying the size of the 281* . deflation window after KEXNW iterations. ==== 282 INTEGER KEXNW 283 PARAMETER ( KEXNW = 5 ) 284* 285* ==== Exceptional shifts: try to cure rare slow convergence 286* . with ad-hoc exceptional shifts every KEXSH iterations. 287* . ==== 288 INTEGER KEXSH 289 PARAMETER ( KEXSH = 6 ) 290* 291* ==== The constants WILK1 and WILK2 are used to form the 292* . exceptional shifts. ==== 293 REAL WILK1, WILK2 294 PARAMETER ( WILK1 = 0.75e0, WILK2 = -0.4375e0 ) 295 REAL ZERO, ONE 296 PARAMETER ( ZERO = 0.0e0, ONE = 1.0e0 ) 297* .. 298* .. Local Scalars .. 299 REAL AA, BB, CC, CS, DD, SN, SS, SWAP 300 INTEGER I, INF, IT, ITMAX, K, KACC22, KBOT, KDU, KS, 301 $ KT, KTOP, KU, KV, KWH, KWTOP, KWV, LD, LS, 302 $ LWKOPT, NDEC, NDFL, NH, NHO, NIBBLE, NMIN, NS, 303 $ NSMAX, NSR, NVE, NW, NWMAX, NWR, NWUPBD 304 LOGICAL SORTED 305 CHARACTER JBCMPZ*2 306* .. 307* .. External Functions .. 308 INTEGER ILAENV 309 EXTERNAL ILAENV 310* .. 311* .. Local Arrays .. 312 REAL ZDUM( 1, 1 ) 313* .. 314* .. External Subroutines .. 315 EXTERNAL SLACPY, SLAHQR, SLANV2, SLAQR3, SLAQR4, SLAQR5 316* .. 317* .. Intrinsic Functions .. 318 INTRINSIC ABS, INT, MAX, MIN, MOD, REAL 319* .. 320* .. Executable Statements .. 321 INFO = 0 322* 323* ==== Quick return for N = 0: nothing to do. ==== 324* 325 IF( N.EQ.0 ) THEN 326 WORK( 1 ) = ONE 327 RETURN 328 END IF 329* 330 IF( N.LE.NTINY ) THEN 331* 332* ==== Tiny matrices must use SLAHQR. ==== 333* 334 LWKOPT = 1 335 IF( LWORK.NE.-1 ) 336 $ CALL SLAHQR( WANTT, WANTZ, N, ILO, IHI, H, LDH, WR, WI, 337 $ ILOZ, IHIZ, Z, LDZ, INFO ) 338 ELSE 339* 340* ==== Use small bulge multi-shift QR with aggressive early 341* . deflation on larger-than-tiny matrices. ==== 342* 343* ==== Hope for the best. ==== 344* 345 INFO = 0 346* 347* ==== Set up job flags for ILAENV. ==== 348* 349 IF( WANTT ) THEN 350 JBCMPZ( 1: 1 ) = 'S' 351 ELSE 352 JBCMPZ( 1: 1 ) = 'E' 353 END IF 354 IF( WANTZ ) THEN 355 JBCMPZ( 2: 2 ) = 'V' 356 ELSE 357 JBCMPZ( 2: 2 ) = 'N' 358 END IF 359* 360* ==== NWR = recommended deflation window size. At this 361* . point, N .GT. NTINY = 15, so there is enough 362* . subdiagonal workspace for NWR.GE.2 as required. 363* . (In fact, there is enough subdiagonal space for 364* . NWR.GE.4.) ==== 365* 366 NWR = ILAENV( 13, 'SLAQR0', JBCMPZ, N, ILO, IHI, LWORK ) 367 NWR = MAX( 2, NWR ) 368 NWR = MIN( IHI-ILO+1, ( N-1 ) / 3, NWR ) 369* 370* ==== NSR = recommended number of simultaneous shifts. 371* . At this point N .GT. NTINY = 15, so there is at 372* . enough subdiagonal workspace for NSR to be even 373* . and greater than or equal to two as required. ==== 374* 375 NSR = ILAENV( 15, 'SLAQR0', JBCMPZ, N, ILO, IHI, LWORK ) 376 NSR = MIN( NSR, ( N-3 ) / 6, IHI-ILO ) 377 NSR = MAX( 2, NSR-MOD( NSR, 2 ) ) 378* 379* ==== Estimate optimal workspace ==== 380* 381* ==== Workspace query call to SLAQR3 ==== 382* 383 CALL SLAQR3( WANTT, WANTZ, N, ILO, IHI, NWR+1, H, LDH, ILOZ, 384 $ IHIZ, Z, LDZ, LS, LD, WR, WI, H, LDH, N, H, LDH, 385 $ N, H, LDH, WORK, -1 ) 386* 387* ==== Optimal workspace = MAX(SLAQR5, SLAQR3) ==== 388* 389 LWKOPT = MAX( 3*NSR / 2, INT( WORK( 1 ) ) ) 390* 391* ==== Quick return in case of workspace query. ==== 392* 393 IF( LWORK.EQ.-1 ) THEN 394 WORK( 1 ) = REAL( LWKOPT ) 395 RETURN 396 END IF 397* 398* ==== SLAHQR/SLAQR0 crossover point ==== 399* 400 NMIN = ILAENV( 12, 'SLAQR0', JBCMPZ, N, ILO, IHI, LWORK ) 401 NMIN = MAX( NTINY, NMIN ) 402* 403* ==== Nibble crossover point ==== 404* 405 NIBBLE = ILAENV( 14, 'SLAQR0', JBCMPZ, N, ILO, IHI, LWORK ) 406 NIBBLE = MAX( 0, NIBBLE ) 407* 408* ==== Accumulate reflections during ttswp? Use block 409* . 2-by-2 structure during matrix-matrix multiply? ==== 410* 411 KACC22 = ILAENV( 16, 'SLAQR0', JBCMPZ, N, ILO, IHI, LWORK ) 412 KACC22 = MAX( 0, KACC22 ) 413 KACC22 = MIN( 2, KACC22 ) 414* 415* ==== NWMAX = the largest possible deflation window for 416* . which there is sufficient workspace. ==== 417* 418 NWMAX = MIN( ( N-1 ) / 3, LWORK / 2 ) 419 NW = NWMAX 420* 421* ==== NSMAX = the Largest number of simultaneous shifts 422* . for which there is sufficient workspace. ==== 423* 424 NSMAX = MIN( ( N-3 ) / 6, 2*LWORK / 3 ) 425 NSMAX = NSMAX - MOD( NSMAX, 2 ) 426* 427* ==== NDFL: an iteration count restarted at deflation. ==== 428* 429 NDFL = 1 430* 431* ==== ITMAX = iteration limit ==== 432* 433 ITMAX = MAX( 30, 2*KEXSH )*MAX( 10, ( IHI-ILO+1 ) ) 434* 435* ==== Last row and column in the active block ==== 436* 437 KBOT = IHI 438* 439* ==== Main Loop ==== 440* 441 DO 80 IT = 1, ITMAX 442* 443* ==== Done when KBOT falls below ILO ==== 444* 445 IF( KBOT.LT.ILO ) 446 $ GO TO 90 447* 448* ==== Locate active block ==== 449* 450 DO 10 K = KBOT, ILO + 1, -1 451 IF( H( K, K-1 ).EQ.ZERO ) 452 $ GO TO 20 453 10 CONTINUE 454 K = ILO 455 20 CONTINUE 456 KTOP = K 457* 458* ==== Select deflation window size: 459* . Typical Case: 460* . If possible and advisable, nibble the entire 461* . active block. If not, use size MIN(NWR,NWMAX) 462* . or MIN(NWR+1,NWMAX) depending upon which has 463* . the smaller corresponding subdiagonal entry 464* . (a heuristic). 465* . 466* . Exceptional Case: 467* . If there have been no deflations in KEXNW or 468* . more iterations, then vary the deflation window 469* . size. At first, because, larger windows are, 470* . in general, more powerful than smaller ones, 471* . rapidly increase the window to the maximum possible. 472* . Then, gradually reduce the window size. ==== 473* 474 NH = KBOT - KTOP + 1 475 NWUPBD = MIN( NH, NWMAX ) 476 IF( NDFL.LT.KEXNW ) THEN 477 NW = MIN( NWUPBD, NWR ) 478 ELSE 479 NW = MIN( NWUPBD, 2*NW ) 480 END IF 481 IF( NW.LT.NWMAX ) THEN 482 IF( NW.GE.NH-1 ) THEN 483 NW = NH 484 ELSE 485 KWTOP = KBOT - NW + 1 486 IF( ABS( H( KWTOP, KWTOP-1 ) ).GT. 487 $ ABS( H( KWTOP-1, KWTOP-2 ) ) )NW = NW + 1 488 END IF 489 END IF 490 IF( NDFL.LT.KEXNW ) THEN 491 NDEC = -1 492 ELSE IF( NDEC.GE.0 .OR. NW.GE.NWUPBD ) THEN 493 NDEC = NDEC + 1 494 IF( NW-NDEC.LT.2 ) 495 $ NDEC = 0 496 NW = NW - NDEC 497 END IF 498* 499* ==== Aggressive early deflation: 500* . split workspace under the subdiagonal into 501* . - an nw-by-nw work array V in the lower 502* . left-hand-corner, 503* . - an NW-by-at-least-NW-but-more-is-better 504* . (NW-by-NHO) horizontal work array along 505* . the bottom edge, 506* . - an at-least-NW-but-more-is-better (NHV-by-NW) 507* . vertical work array along the left-hand-edge. 508* . ==== 509* 510 KV = N - NW + 1 511 KT = NW + 1 512 NHO = ( N-NW-1 ) - KT + 1 513 KWV = NW + 2 514 NVE = ( N-NW ) - KWV + 1 515* 516* ==== Aggressive early deflation ==== 517* 518 CALL SLAQR3( WANTT, WANTZ, N, KTOP, KBOT, NW, H, LDH, ILOZ, 519 $ IHIZ, Z, LDZ, LS, LD, WR, WI, H( KV, 1 ), LDH, 520 $ NHO, H( KV, KT ), LDH, NVE, H( KWV, 1 ), LDH, 521 $ WORK, LWORK ) 522* 523* ==== Adjust KBOT accounting for new deflations. ==== 524* 525 KBOT = KBOT - LD 526* 527* ==== KS points to the shifts. ==== 528* 529 KS = KBOT - LS + 1 530* 531* ==== Skip an expensive QR sweep if there is a (partly 532* . heuristic) reason to expect that many eigenvalues 533* . will deflate without it. Here, the QR sweep is 534* . skipped if many eigenvalues have just been deflated 535* . or if the remaining active block is small. 536* 537 IF( ( LD.EQ.0 ) .OR. ( ( 100*LD.LE.NW*NIBBLE ) .AND. ( KBOT- 538 $ KTOP+1.GT.MIN( NMIN, NWMAX ) ) ) ) THEN 539* 540* ==== NS = nominal number of simultaneous shifts. 541* . This may be lowered (slightly) if SLAQR3 542* . did not provide that many shifts. ==== 543* 544 NS = MIN( NSMAX, NSR, MAX( 2, KBOT-KTOP ) ) 545 NS = NS - MOD( NS, 2 ) 546* 547* ==== If there have been no deflations 548* . in a multiple of KEXSH iterations, 549* . then try exceptional shifts. 550* . Otherwise use shifts provided by 551* . SLAQR3 above or from the eigenvalues 552* . of a trailing principal submatrix. ==== 553* 554 IF( MOD( NDFL, KEXSH ).EQ.0 ) THEN 555 KS = KBOT - NS + 1 556 DO 30 I = KBOT, MAX( KS+1, KTOP+2 ), -2 557 SS = ABS( H( I, I-1 ) ) + ABS( H( I-1, I-2 ) ) 558 AA = WILK1*SS + H( I, I ) 559 BB = SS 560 CC = WILK2*SS 561 DD = AA 562 CALL SLANV2( AA, BB, CC, DD, WR( I-1 ), WI( I-1 ), 563 $ WR( I ), WI( I ), CS, SN ) 564 30 CONTINUE 565 IF( KS.EQ.KTOP ) THEN 566 WR( KS+1 ) = H( KS+1, KS+1 ) 567 WI( KS+1 ) = ZERO 568 WR( KS ) = WR( KS+1 ) 569 WI( KS ) = WI( KS+1 ) 570 END IF 571 ELSE 572* 573* ==== Got NS/2 or fewer shifts? Use SLAQR4 or 574* . SLAHQR on a trailing principal submatrix to 575* . get more. (Since NS.LE.NSMAX.LE.(N-3)/6, 576* . there is enough space below the subdiagonal 577* . to fit an NS-by-NS scratch array.) ==== 578* 579 IF( KBOT-KS+1.LE.NS / 2 ) THEN 580 KS = KBOT - NS + 1 581 KT = N - NS + 1 582 CALL SLACPY( 'A', NS, NS, H( KS, KS ), LDH, 583 $ H( KT, 1 ), LDH ) 584 IF( NS.GT.NMIN ) THEN 585 CALL SLAQR4( .false., .false., NS, 1, NS, 586 $ H( KT, 1 ), LDH, WR( KS ), 587 $ WI( KS ), 1, 1, ZDUM, 1, WORK, 588 $ LWORK, INF ) 589 ELSE 590 CALL SLAHQR( .false., .false., NS, 1, NS, 591 $ H( KT, 1 ), LDH, WR( KS ), 592 $ WI( KS ), 1, 1, ZDUM, 1, INF ) 593 END IF 594 KS = KS + INF 595* 596* ==== In case of a rare QR failure use 597* . eigenvalues of the trailing 2-by-2 598* . principal submatrix. ==== 599* 600 IF( KS.GE.KBOT ) THEN 601 AA = H( KBOT-1, KBOT-1 ) 602 CC = H( KBOT, KBOT-1 ) 603 BB = H( KBOT-1, KBOT ) 604 DD = H( KBOT, KBOT ) 605 CALL SLANV2( AA, BB, CC, DD, WR( KBOT-1 ), 606 $ WI( KBOT-1 ), WR( KBOT ), 607 $ WI( KBOT ), CS, SN ) 608 KS = KBOT - 1 609 END IF 610 END IF 611* 612 IF( KBOT-KS+1.GT.NS ) THEN 613* 614* ==== Sort the shifts (Helps a little) 615* . Bubble sort keeps complex conjugate 616* . pairs together. ==== 617* 618 SORTED = .false. 619 DO 50 K = KBOT, KS + 1, -1 620 IF( SORTED ) 621 $ GO TO 60 622 SORTED = .true. 623 DO 40 I = KS, K - 1 624 IF( ABS( WR( I ) )+ABS( WI( I ) ).LT. 625 $ ABS( WR( I+1 ) )+ABS( WI( I+1 ) ) ) THEN 626 SORTED = .false. 627* 628 SWAP = WR( I ) 629 WR( I ) = WR( I+1 ) 630 WR( I+1 ) = SWAP 631* 632 SWAP = WI( I ) 633 WI( I ) = WI( I+1 ) 634 WI( I+1 ) = SWAP 635 END IF 636 40 CONTINUE 637 50 CONTINUE 638 60 CONTINUE 639 END IF 640* 641* ==== Shuffle shifts into pairs of real shifts 642* . and pairs of complex conjugate shifts 643* . assuming complex conjugate shifts are 644* . already adjacent to one another. (Yes, 645* . they are.) ==== 646* 647 DO 70 I = KBOT, KS + 2, -2 648 IF( WI( I ).NE.-WI( I-1 ) ) THEN 649* 650 SWAP = WR( I ) 651 WR( I ) = WR( I-1 ) 652 WR( I-1 ) = WR( I-2 ) 653 WR( I-2 ) = SWAP 654* 655 SWAP = WI( I ) 656 WI( I ) = WI( I-1 ) 657 WI( I-1 ) = WI( I-2 ) 658 WI( I-2 ) = SWAP 659 END IF 660 70 CONTINUE 661 END IF 662* 663* ==== If there are only two shifts and both are 664* . real, then use only one. ==== 665* 666 IF( KBOT-KS+1.EQ.2 ) THEN 667 IF( WI( KBOT ).EQ.ZERO ) THEN 668 IF( ABS( WR( KBOT )-H( KBOT, KBOT ) ).LT. 669 $ ABS( WR( KBOT-1 )-H( KBOT, KBOT ) ) ) THEN 670 WR( KBOT-1 ) = WR( KBOT ) 671 ELSE 672 WR( KBOT ) = WR( KBOT-1 ) 673 END IF 674 END IF 675 END IF 676* 677* ==== Use up to NS of the the smallest magnitude 678* . shifts. If there aren't NS shifts available, 679* . then use them all, possibly dropping one to 680* . make the number of shifts even. ==== 681* 682 NS = MIN( NS, KBOT-KS+1 ) 683 NS = NS - MOD( NS, 2 ) 684 KS = KBOT - NS + 1 685* 686* ==== Small-bulge multi-shift QR sweep: 687* . split workspace under the subdiagonal into 688* . - a KDU-by-KDU work array U in the lower 689* . left-hand-corner, 690* . - a KDU-by-at-least-KDU-but-more-is-better 691* . (KDU-by-NHo) horizontal work array WH along 692* . the bottom edge, 693* . - and an at-least-KDU-but-more-is-better-by-KDU 694* . (NVE-by-KDU) vertical work WV arrow along 695* . the left-hand-edge. ==== 696* 697 KDU = 2*NS 698 KU = N - KDU + 1 699 KWH = KDU + 1 700 NHO = ( N-KDU+1-4 ) - ( KDU+1 ) + 1 701 KWV = KDU + 4 702 NVE = N - KDU - KWV + 1 703* 704* ==== Small-bulge multi-shift QR sweep ==== 705* 706 CALL SLAQR5( WANTT, WANTZ, KACC22, N, KTOP, KBOT, NS, 707 $ WR( KS ), WI( KS ), H, LDH, ILOZ, IHIZ, Z, 708 $ LDZ, WORK, 3, H( KU, 1 ), LDH, NVE, 709 $ H( KWV, 1 ), LDH, NHO, H( KU, KWH ), LDH ) 710 END IF 711* 712* ==== Note progress (or the lack of it). ==== 713* 714 IF( LD.GT.0 ) THEN 715 NDFL = 1 716 ELSE 717 NDFL = NDFL + 1 718 END IF 719* 720* ==== End of main loop ==== 721 80 CONTINUE 722* 723* ==== Iteration limit exceeded. Set INFO to show where 724* . the problem occurred and exit. ==== 725* 726 INFO = KBOT 727 90 CONTINUE 728 END IF 729* 730* ==== Return the optimal value of LWORK. ==== 731* 732 WORK( 1 ) = REAL( LWKOPT ) 733* 734* ==== End of SLAQR0 ==== 735* 736 END 737