1*> \brief \b DHGEQZ 2* 3* =========== DOCUMENTATION =========== 4* 5* Online html documentation available at 6* http://www.netlib.org/lapack/explore-html/ 7* 8*> \htmlonly 9*> Download DHGEQZ + dependencies 10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dhgeqz.f"> 11*> [TGZ]</a> 12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dhgeqz.f"> 13*> [ZIP]</a> 14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dhgeqz.f"> 15*> [TXT]</a> 16*> \endhtmlonly 17* 18* Definition: 19* =========== 20* 21* SUBROUTINE DHGEQZ( JOB, COMPQ, COMPZ, N, ILO, IHI, H, LDH, T, LDT, 22* ALPHAR, ALPHAI, BETA, Q, LDQ, Z, LDZ, WORK, 23* LWORK, INFO ) 24* 25* .. Scalar Arguments .. 26* CHARACTER COMPQ, COMPZ, JOB 27* INTEGER IHI, ILO, INFO, LDH, LDQ, LDT, LDZ, LWORK, N 28* .. 29* .. Array Arguments .. 30* DOUBLE PRECISION ALPHAI( * ), ALPHAR( * ), BETA( * ), 31* $ H( LDH, * ), Q( LDQ, * ), T( LDT, * ), 32* $ WORK( * ), Z( LDZ, * ) 33* .. 34* 35* 36*> \par Purpose: 37* ============= 38*> 39*> \verbatim 40*> 41*> DHGEQZ computes the eigenvalues of a real matrix pair (H,T), 42*> where H is an upper Hessenberg matrix and T is upper triangular, 43*> using the double-shift QZ method. 44*> Matrix pairs of this type are produced by the reduction to 45*> generalized upper Hessenberg form of a real matrix pair (A,B): 46*> 47*> A = Q1*H*Z1**T, B = Q1*T*Z1**T, 48*> 49*> as computed by DGGHRD. 50*> 51*> If JOB='S', then the Hessenberg-triangular pair (H,T) is 52*> also reduced to generalized Schur form, 53*> 54*> H = Q*S*Z**T, T = Q*P*Z**T, 55*> 56*> where Q and Z are orthogonal matrices, P is an upper triangular 57*> matrix, and S is a quasi-triangular matrix with 1-by-1 and 2-by-2 58*> diagonal blocks. 59*> 60*> The 1-by-1 blocks correspond to real eigenvalues of the matrix pair 61*> (H,T) and the 2-by-2 blocks correspond to complex conjugate pairs of 62*> eigenvalues. 63*> 64*> Additionally, the 2-by-2 upper triangular diagonal blocks of P 65*> corresponding to 2-by-2 blocks of S are reduced to positive diagonal 66*> form, i.e., if S(j+1,j) is non-zero, then P(j+1,j) = P(j,j+1) = 0, 67*> P(j,j) > 0, and P(j+1,j+1) > 0. 68*> 69*> Optionally, the orthogonal matrix Q from the generalized Schur 70*> factorization may be postmultiplied into an input matrix Q1, and the 71*> orthogonal matrix Z may be postmultiplied into an input matrix Z1. 72*> If Q1 and Z1 are the orthogonal matrices from DGGHRD that reduced 73*> the matrix pair (A,B) to generalized upper Hessenberg form, then the 74*> output matrices Q1*Q and Z1*Z are the orthogonal factors from the 75*> generalized Schur factorization of (A,B): 76*> 77*> A = (Q1*Q)*S*(Z1*Z)**T, B = (Q1*Q)*P*(Z1*Z)**T. 78*> 79*> To avoid overflow, eigenvalues of the matrix pair (H,T) (equivalently, 80*> of (A,B)) are computed as a pair of values (alpha,beta), where alpha is 81*> complex and beta real. 82*> If beta is nonzero, lambda = alpha / beta is an eigenvalue of the 83*> generalized nonsymmetric eigenvalue problem (GNEP) 84*> A*x = lambda*B*x 85*> and if alpha is nonzero, mu = beta / alpha is an eigenvalue of the 86*> alternate form of the GNEP 87*> mu*A*y = B*y. 88*> Real eigenvalues can be read directly from the generalized Schur 89*> form: 90*> alpha = S(i,i), beta = P(i,i). 91*> 92*> Ref: C.B. Moler & G.W. Stewart, "An Algorithm for Generalized Matrix 93*> Eigenvalue Problems", SIAM J. Numer. Anal., 10(1973), 94*> pp. 241--256. 95*> \endverbatim 96* 97* Arguments: 98* ========== 99* 100*> \param[in] JOB 101*> \verbatim 102*> JOB is CHARACTER*1 103*> = 'E': Compute eigenvalues only; 104*> = 'S': Compute eigenvalues and the Schur form. 105*> \endverbatim 106*> 107*> \param[in] COMPQ 108*> \verbatim 109*> COMPQ is CHARACTER*1 110*> = 'N': Left Schur vectors (Q) are not computed; 111*> = 'I': Q is initialized to the unit matrix and the matrix Q 112*> of left Schur vectors of (H,T) is returned; 113*> = 'V': Q must contain an orthogonal matrix Q1 on entry and 114*> the product Q1*Q is returned. 115*> \endverbatim 116*> 117*> \param[in] COMPZ 118*> \verbatim 119*> COMPZ is CHARACTER*1 120*> = 'N': Right Schur vectors (Z) are not computed; 121*> = 'I': Z is initialized to the unit matrix and the matrix Z 122*> of right Schur vectors of (H,T) is returned; 123*> = 'V': Z must contain an orthogonal matrix Z1 on entry and 124*> the product Z1*Z is returned. 125*> \endverbatim 126*> 127*> \param[in] N 128*> \verbatim 129*> N is INTEGER 130*> The order of the matrices H, T, Q, and Z. N >= 0. 131*> \endverbatim 132*> 133*> \param[in] ILO 134*> \verbatim 135*> ILO is INTEGER 136*> \endverbatim 137*> 138*> \param[in] IHI 139*> \verbatim 140*> IHI is INTEGER 141*> ILO and IHI mark the rows and columns of H which are in 142*> Hessenberg form. It is assumed that A is already upper 143*> triangular in rows and columns 1:ILO-1 and IHI+1:N. 144*> If N > 0, 1 <= ILO <= IHI <= N; if N = 0, ILO=1 and IHI=0. 145*> \endverbatim 146*> 147*> \param[in,out] H 148*> \verbatim 149*> H is DOUBLE PRECISION array, dimension (LDH, N) 150*> On entry, the N-by-N upper Hessenberg matrix H. 151*> On exit, if JOB = 'S', H contains the upper quasi-triangular 152*> matrix S from the generalized Schur factorization. 153*> If JOB = 'E', the diagonal blocks of H match those of S, but 154*> the rest of H is unspecified. 155*> \endverbatim 156*> 157*> \param[in] LDH 158*> \verbatim 159*> LDH is INTEGER 160*> The leading dimension of the array H. LDH >= max( 1, N ). 161*> \endverbatim 162*> 163*> \param[in,out] T 164*> \verbatim 165*> T is DOUBLE PRECISION array, dimension (LDT, N) 166*> On entry, the N-by-N upper triangular matrix T. 167*> On exit, if JOB = 'S', T contains the upper triangular 168*> matrix P from the generalized Schur factorization; 169*> 2-by-2 diagonal blocks of P corresponding to 2-by-2 blocks of S 170*> are reduced to positive diagonal form, i.e., if H(j+1,j) is 171*> non-zero, then T(j+1,j) = T(j,j+1) = 0, T(j,j) > 0, and 172*> T(j+1,j+1) > 0. 173*> If JOB = 'E', the diagonal blocks of T match those of P, but 174*> the rest of T is unspecified. 175*> \endverbatim 176*> 177*> \param[in] LDT 178*> \verbatim 179*> LDT is INTEGER 180*> The leading dimension of the array T. LDT >= max( 1, N ). 181*> \endverbatim 182*> 183*> \param[out] ALPHAR 184*> \verbatim 185*> ALPHAR is DOUBLE PRECISION array, dimension (N) 186*> The real parts of each scalar alpha defining an eigenvalue 187*> of GNEP. 188*> \endverbatim 189*> 190*> \param[out] ALPHAI 191*> \verbatim 192*> ALPHAI is DOUBLE PRECISION array, dimension (N) 193*> The imaginary parts of each scalar alpha defining an 194*> eigenvalue of GNEP. 195*> If ALPHAI(j) is zero, then the j-th eigenvalue is real; if 196*> positive, then the j-th and (j+1)-st eigenvalues are a 197*> complex conjugate pair, with ALPHAI(j+1) = -ALPHAI(j). 198*> \endverbatim 199*> 200*> \param[out] BETA 201*> \verbatim 202*> BETA is DOUBLE PRECISION array, dimension (N) 203*> The scalars beta that define the eigenvalues of GNEP. 204*> Together, the quantities alpha = (ALPHAR(j),ALPHAI(j)) and 205*> beta = BETA(j) represent the j-th eigenvalue of the matrix 206*> pair (A,B), in one of the forms lambda = alpha/beta or 207*> mu = beta/alpha. Since either lambda or mu may overflow, 208*> they should not, in general, be computed. 209*> \endverbatim 210*> 211*> \param[in,out] Q 212*> \verbatim 213*> Q is DOUBLE PRECISION array, dimension (LDQ, N) 214*> On entry, if COMPQ = 'V', the orthogonal matrix Q1 used in 215*> the reduction of (A,B) to generalized Hessenberg form. 216*> On exit, if COMPQ = 'I', the orthogonal matrix of left Schur 217*> vectors of (H,T), and if COMPQ = 'V', the orthogonal matrix 218*> of left Schur vectors of (A,B). 219*> Not referenced if COMPQ = 'N'. 220*> \endverbatim 221*> 222*> \param[in] LDQ 223*> \verbatim 224*> LDQ is INTEGER 225*> The leading dimension of the array Q. LDQ >= 1. 226*> If COMPQ='V' or 'I', then LDQ >= N. 227*> \endverbatim 228*> 229*> \param[in,out] Z 230*> \verbatim 231*> Z is DOUBLE PRECISION array, dimension (LDZ, N) 232*> On entry, if COMPZ = 'V', the orthogonal matrix Z1 used in 233*> the reduction of (A,B) to generalized Hessenberg form. 234*> On exit, if COMPZ = 'I', the orthogonal matrix of 235*> right Schur vectors of (H,T), and if COMPZ = 'V', the 236*> orthogonal matrix of right Schur vectors of (A,B). 237*> Not referenced if COMPZ = 'N'. 238*> \endverbatim 239*> 240*> \param[in] LDZ 241*> \verbatim 242*> LDZ is INTEGER 243*> The leading dimension of the array Z. LDZ >= 1. 244*> If COMPZ='V' or 'I', then LDZ >= N. 245*> \endverbatim 246*> 247*> \param[out] WORK 248*> \verbatim 249*> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK)) 250*> On exit, if INFO >= 0, WORK(1) returns the optimal LWORK. 251*> \endverbatim 252*> 253*> \param[in] LWORK 254*> \verbatim 255*> LWORK is INTEGER 256*> The dimension of the array WORK. LWORK >= max(1,N). 257*> 258*> If LWORK = -1, then a workspace query is assumed; the routine 259*> only calculates the optimal size of the WORK array, returns 260*> this value as the first entry of the WORK array, and no error 261*> message related to LWORK is issued by XERBLA. 262*> \endverbatim 263*> 264*> \param[out] INFO 265*> \verbatim 266*> INFO is INTEGER 267*> = 0: successful exit 268*> < 0: if INFO = -i, the i-th argument had an illegal value 269*> = 1,...,N: the QZ iteration did not converge. (H,T) is not 270*> in Schur form, but ALPHAR(i), ALPHAI(i), and 271*> BETA(i), i=INFO+1,...,N should be correct. 272*> = N+1,...,2*N: the shift calculation failed. (H,T) is not 273*> in Schur form, but ALPHAR(i), ALPHAI(i), and 274*> BETA(i), i=INFO-N+1,...,N should be correct. 275*> \endverbatim 276* 277* Authors: 278* ======== 279* 280*> \author Univ. of Tennessee 281*> \author Univ. of California Berkeley 282*> \author Univ. of Colorado Denver 283*> \author NAG Ltd. 284* 285*> \ingroup doubleGEcomputational 286* 287*> \par Further Details: 288* ===================== 289*> 290*> \verbatim 291*> 292*> Iteration counters: 293*> 294*> JITER -- counts iterations. 295*> IITER -- counts iterations run since ILAST was last 296*> changed. This is therefore reset only when a 1-by-1 or 297*> 2-by-2 block deflates off the bottom. 298*> \endverbatim 299*> 300* ===================================================================== 301 SUBROUTINE DHGEQZ( JOB, COMPQ, COMPZ, N, ILO, IHI, H, LDH, T, LDT, 302 $ ALPHAR, ALPHAI, BETA, Q, LDQ, Z, LDZ, WORK, 303 $ LWORK, INFO ) 304* 305* -- LAPACK computational routine -- 306* -- LAPACK is a software package provided by Univ. of Tennessee, -- 307* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- 308* 309* .. Scalar Arguments .. 310 CHARACTER COMPQ, COMPZ, JOB 311 INTEGER IHI, ILO, INFO, LDH, LDQ, LDT, LDZ, LWORK, N 312* .. 313* .. Array Arguments .. 314 DOUBLE PRECISION ALPHAI( * ), ALPHAR( * ), BETA( * ), 315 $ H( LDH, * ), Q( LDQ, * ), T( LDT, * ), 316 $ WORK( * ), Z( LDZ, * ) 317* .. 318* 319* ===================================================================== 320* 321* .. Parameters .. 322* $ SAFETY = 1.0E+0 ) 323 DOUBLE PRECISION HALF, ZERO, ONE, SAFETY 324 PARAMETER ( HALF = 0.5D+0, ZERO = 0.0D+0, ONE = 1.0D+0, 325 $ SAFETY = 1.0D+2 ) 326* .. 327* .. Local Scalars .. 328 LOGICAL ILAZR2, ILAZRO, ILPIVT, ILQ, ILSCHR, ILZ, 329 $ LQUERY 330 INTEGER ICOMPQ, ICOMPZ, IFIRST, IFRSTM, IITER, ILAST, 331 $ ILASTM, IN, ISCHUR, ISTART, J, JC, JCH, JITER, 332 $ JR, MAXIT 333 DOUBLE PRECISION A11, A12, A1I, A1R, A21, A22, A2I, A2R, AD11, 334 $ AD11L, AD12, AD12L, AD21, AD21L, AD22, AD22L, 335 $ AD32L, AN, ANORM, ASCALE, ATOL, B11, B1A, B1I, 336 $ B1R, B22, B2A, B2I, B2R, BN, BNORM, BSCALE, 337 $ BTOL, C, C11I, C11R, C12, C21, C22I, C22R, CL, 338 $ CQ, CR, CZ, ESHIFT, S, S1, S1INV, S2, SAFMAX, 339 $ SAFMIN, SCALE, SL, SQI, SQR, SR, SZI, SZR, T1, 340 $ TAU, TEMP, TEMP2, TEMPI, TEMPR, U1, U12, U12L, 341 $ U2, ULP, VS, W11, W12, W21, W22, WABS, WI, WR, 342 $ WR2 343* .. 344* .. Local Arrays .. 345 DOUBLE PRECISION V( 3 ) 346* .. 347* .. External Functions .. 348 LOGICAL LSAME 349 DOUBLE PRECISION DLAMCH, DLANHS, DLAPY2, DLAPY3 350 EXTERNAL LSAME, DLAMCH, DLANHS, DLAPY2, DLAPY3 351* .. 352* .. External Subroutines .. 353 EXTERNAL DLAG2, DLARFG, DLARTG, DLASET, DLASV2, DROT, 354 $ XERBLA 355* .. 356* .. Intrinsic Functions .. 357 INTRINSIC ABS, DBLE, MAX, MIN, SQRT 358* .. 359* .. Executable Statements .. 360* 361* Decode JOB, COMPQ, COMPZ 362* 363 IF( LSAME( JOB, 'E' ) ) THEN 364 ILSCHR = .FALSE. 365 ISCHUR = 1 366 ELSE IF( LSAME( JOB, 'S' ) ) THEN 367 ILSCHR = .TRUE. 368 ISCHUR = 2 369 ELSE 370 ISCHUR = 0 371 END IF 372* 373 IF( LSAME( COMPQ, 'N' ) ) THEN 374 ILQ = .FALSE. 375 ICOMPQ = 1 376 ELSE IF( LSAME( COMPQ, 'V' ) ) THEN 377 ILQ = .TRUE. 378 ICOMPQ = 2 379 ELSE IF( LSAME( COMPQ, 'I' ) ) THEN 380 ILQ = .TRUE. 381 ICOMPQ = 3 382 ELSE 383 ICOMPQ = 0 384 END IF 385* 386 IF( LSAME( COMPZ, 'N' ) ) THEN 387 ILZ = .FALSE. 388 ICOMPZ = 1 389 ELSE IF( LSAME( COMPZ, 'V' ) ) THEN 390 ILZ = .TRUE. 391 ICOMPZ = 2 392 ELSE IF( LSAME( COMPZ, 'I' ) ) THEN 393 ILZ = .TRUE. 394 ICOMPZ = 3 395 ELSE 396 ICOMPZ = 0 397 END IF 398* 399* Check Argument Values 400* 401 INFO = 0 402 WORK( 1 ) = MAX( 1, N ) 403 LQUERY = ( LWORK.EQ.-1 ) 404 IF( ISCHUR.EQ.0 ) THEN 405 INFO = -1 406 ELSE IF( ICOMPQ.EQ.0 ) THEN 407 INFO = -2 408 ELSE IF( ICOMPZ.EQ.0 ) THEN 409 INFO = -3 410 ELSE IF( N.LT.0 ) THEN 411 INFO = -4 412 ELSE IF( ILO.LT.1 ) THEN 413 INFO = -5 414 ELSE IF( IHI.GT.N .OR. IHI.LT.ILO-1 ) THEN 415 INFO = -6 416 ELSE IF( LDH.LT.N ) THEN 417 INFO = -8 418 ELSE IF( LDT.LT.N ) THEN 419 INFO = -10 420 ELSE IF( LDQ.LT.1 .OR. ( ILQ .AND. LDQ.LT.N ) ) THEN 421 INFO = -15 422 ELSE IF( LDZ.LT.1 .OR. ( ILZ .AND. LDZ.LT.N ) ) THEN 423 INFO = -17 424 ELSE IF( LWORK.LT.MAX( 1, N ) .AND. .NOT.LQUERY ) THEN 425 INFO = -19 426 END IF 427 IF( INFO.NE.0 ) THEN 428 CALL XERBLA( 'DHGEQZ', -INFO ) 429 RETURN 430 ELSE IF( LQUERY ) THEN 431 RETURN 432 END IF 433* 434* Quick return if possible 435* 436 IF( N.LE.0 ) THEN 437 WORK( 1 ) = DBLE( 1 ) 438 RETURN 439 END IF 440* 441* Initialize Q and Z 442* 443 IF( ICOMPQ.EQ.3 ) 444 $ CALL DLASET( 'Full', N, N, ZERO, ONE, Q, LDQ ) 445 IF( ICOMPZ.EQ.3 ) 446 $ CALL DLASET( 'Full', N, N, ZERO, ONE, Z, LDZ ) 447* 448* Machine Constants 449* 450 IN = IHI + 1 - ILO 451 SAFMIN = DLAMCH( 'S' ) 452 SAFMAX = ONE / SAFMIN 453 ULP = DLAMCH( 'E' )*DLAMCH( 'B' ) 454 ANORM = DLANHS( 'F', IN, H( ILO, ILO ), LDH, WORK ) 455 BNORM = DLANHS( 'F', IN, T( ILO, ILO ), LDT, WORK ) 456 ATOL = MAX( SAFMIN, ULP*ANORM ) 457 BTOL = MAX( SAFMIN, ULP*BNORM ) 458 ASCALE = ONE / MAX( SAFMIN, ANORM ) 459 BSCALE = ONE / MAX( SAFMIN, BNORM ) 460* 461* Set Eigenvalues IHI+1:N 462* 463 DO 30 J = IHI + 1, N 464 IF( T( J, J ).LT.ZERO ) THEN 465 IF( ILSCHR ) THEN 466 DO 10 JR = 1, J 467 H( JR, J ) = -H( JR, J ) 468 T( JR, J ) = -T( JR, J ) 469 10 CONTINUE 470 ELSE 471 H( J, J ) = -H( J, J ) 472 T( J, J ) = -T( J, J ) 473 END IF 474 IF( ILZ ) THEN 475 DO 20 JR = 1, N 476 Z( JR, J ) = -Z( JR, J ) 477 20 CONTINUE 478 END IF 479 END IF 480 ALPHAR( J ) = H( J, J ) 481 ALPHAI( J ) = ZERO 482 BETA( J ) = T( J, J ) 483 30 CONTINUE 484* 485* If IHI < ILO, skip QZ steps 486* 487 IF( IHI.LT.ILO ) 488 $ GO TO 380 489* 490* MAIN QZ ITERATION LOOP 491* 492* Initialize dynamic indices 493* 494* Eigenvalues ILAST+1:N have been found. 495* Column operations modify rows IFRSTM:whatever. 496* Row operations modify columns whatever:ILASTM. 497* 498* If only eigenvalues are being computed, then 499* IFRSTM is the row of the last splitting row above row ILAST; 500* this is always at least ILO. 501* IITER counts iterations since the last eigenvalue was found, 502* to tell when to use an extraordinary shift. 503* MAXIT is the maximum number of QZ sweeps allowed. 504* 505 ILAST = IHI 506 IF( ILSCHR ) THEN 507 IFRSTM = 1 508 ILASTM = N 509 ELSE 510 IFRSTM = ILO 511 ILASTM = IHI 512 END IF 513 IITER = 0 514 ESHIFT = ZERO 515 MAXIT = 30*( IHI-ILO+1 ) 516* 517 DO 360 JITER = 1, MAXIT 518* 519* Split the matrix if possible. 520* 521* Two tests: 522* 1: H(j,j-1)=0 or j=ILO 523* 2: T(j,j)=0 524* 525 IF( ILAST.EQ.ILO ) THEN 526* 527* Special case: j=ILAST 528* 529 GO TO 80 530 ELSE 531 IF( ABS( H( ILAST, ILAST-1 ) ).LE.MAX( SAFMIN, ULP*( 532 $ ABS( H( ILAST, ILAST ) ) + ABS( H( ILAST-1, ILAST-1 ) ) 533 $ ) ) ) THEN 534 H( ILAST, ILAST-1 ) = ZERO 535 GO TO 80 536 END IF 537 END IF 538* 539 IF( ABS( T( ILAST, ILAST ) ).LE.MAX( SAFMIN, ULP*( 540 $ ABS( T( ILAST - 1, ILAST ) ) + ABS( T( ILAST-1, ILAST-1 ) 541 $ ) ) ) ) THEN 542 T( ILAST, ILAST ) = ZERO 543 GO TO 70 544 END IF 545* 546* General case: j<ILAST 547* 548 DO 60 J = ILAST - 1, ILO, -1 549* 550* Test 1: for H(j,j-1)=0 or j=ILO 551* 552 IF( J.EQ.ILO ) THEN 553 ILAZRO = .TRUE. 554 ELSE 555 IF( ABS( H( J, J-1 ) ).LE.MAX( SAFMIN, ULP*( 556 $ ABS( H( J, J ) ) + ABS( H( J-1, J-1 ) ) 557 $ ) ) ) THEN 558 H( J, J-1 ) = ZERO 559 ILAZRO = .TRUE. 560 ELSE 561 ILAZRO = .FALSE. 562 END IF 563 END IF 564* 565* Test 2: for T(j,j)=0 566* 567 TEMP = ABS ( T( J, J + 1 ) ) 568 IF ( J .GT. ILO ) 569 $ TEMP = TEMP + ABS ( T( J - 1, J ) ) 570 IF( ABS( T( J, J ) ).LT.MAX( SAFMIN,ULP*TEMP ) ) THEN 571 T( J, J ) = ZERO 572* 573* Test 1a: Check for 2 consecutive small subdiagonals in A 574* 575 ILAZR2 = .FALSE. 576 IF( .NOT.ILAZRO ) THEN 577 TEMP = ABS( H( J, J-1 ) ) 578 TEMP2 = ABS( H( J, J ) ) 579 TEMPR = MAX( TEMP, TEMP2 ) 580 IF( TEMPR.LT.ONE .AND. TEMPR.NE.ZERO ) THEN 581 TEMP = TEMP / TEMPR 582 TEMP2 = TEMP2 / TEMPR 583 END IF 584 IF( TEMP*( ASCALE*ABS( H( J+1, J ) ) ).LE.TEMP2* 585 $ ( ASCALE*ATOL ) )ILAZR2 = .TRUE. 586 END IF 587* 588* If both tests pass (1 & 2), i.e., the leading diagonal 589* element of B in the block is zero, split a 1x1 block off 590* at the top. (I.e., at the J-th row/column) The leading 591* diagonal element of the remainder can also be zero, so 592* this may have to be done repeatedly. 593* 594 IF( ILAZRO .OR. ILAZR2 ) THEN 595 DO 40 JCH = J, ILAST - 1 596 TEMP = H( JCH, JCH ) 597 CALL DLARTG( TEMP, H( JCH+1, JCH ), C, S, 598 $ H( JCH, JCH ) ) 599 H( JCH+1, JCH ) = ZERO 600 CALL DROT( ILASTM-JCH, H( JCH, JCH+1 ), LDH, 601 $ H( JCH+1, JCH+1 ), LDH, C, S ) 602 CALL DROT( ILASTM-JCH, T( JCH, JCH+1 ), LDT, 603 $ T( JCH+1, JCH+1 ), LDT, C, S ) 604 IF( ILQ ) 605 $ CALL DROT( N, Q( 1, JCH ), 1, Q( 1, JCH+1 ), 1, 606 $ C, S ) 607 IF( ILAZR2 ) 608 $ H( JCH, JCH-1 ) = H( JCH, JCH-1 )*C 609 ILAZR2 = .FALSE. 610 IF( ABS( T( JCH+1, JCH+1 ) ).GE.BTOL ) THEN 611 IF( JCH+1.GE.ILAST ) THEN 612 GO TO 80 613 ELSE 614 IFIRST = JCH + 1 615 GO TO 110 616 END IF 617 END IF 618 T( JCH+1, JCH+1 ) = ZERO 619 40 CONTINUE 620 GO TO 70 621 ELSE 622* 623* Only test 2 passed -- chase the zero to T(ILAST,ILAST) 624* Then process as in the case T(ILAST,ILAST)=0 625* 626 DO 50 JCH = J, ILAST - 1 627 TEMP = T( JCH, JCH+1 ) 628 CALL DLARTG( TEMP, T( JCH+1, JCH+1 ), C, S, 629 $ T( JCH, JCH+1 ) ) 630 T( JCH+1, JCH+1 ) = ZERO 631 IF( JCH.LT.ILASTM-1 ) 632 $ CALL DROT( ILASTM-JCH-1, T( JCH, JCH+2 ), LDT, 633 $ T( JCH+1, JCH+2 ), LDT, C, S ) 634 CALL DROT( ILASTM-JCH+2, H( JCH, JCH-1 ), LDH, 635 $ H( JCH+1, JCH-1 ), LDH, C, S ) 636 IF( ILQ ) 637 $ CALL DROT( N, Q( 1, JCH ), 1, Q( 1, JCH+1 ), 1, 638 $ C, S ) 639 TEMP = H( JCH+1, JCH ) 640 CALL DLARTG( TEMP, H( JCH+1, JCH-1 ), C, S, 641 $ H( JCH+1, JCH ) ) 642 H( JCH+1, JCH-1 ) = ZERO 643 CALL DROT( JCH+1-IFRSTM, H( IFRSTM, JCH ), 1, 644 $ H( IFRSTM, JCH-1 ), 1, C, S ) 645 CALL DROT( JCH-IFRSTM, T( IFRSTM, JCH ), 1, 646 $ T( IFRSTM, JCH-1 ), 1, C, S ) 647 IF( ILZ ) 648 $ CALL DROT( N, Z( 1, JCH ), 1, Z( 1, JCH-1 ), 1, 649 $ C, S ) 650 50 CONTINUE 651 GO TO 70 652 END IF 653 ELSE IF( ILAZRO ) THEN 654* 655* Only test 1 passed -- work on J:ILAST 656* 657 IFIRST = J 658 GO TO 110 659 END IF 660* 661* Neither test passed -- try next J 662* 663 60 CONTINUE 664* 665* (Drop-through is "impossible") 666* 667 INFO = N + 1 668 GO TO 420 669* 670* T(ILAST,ILAST)=0 -- clear H(ILAST,ILAST-1) to split off a 671* 1x1 block. 672* 673 70 CONTINUE 674 TEMP = H( ILAST, ILAST ) 675 CALL DLARTG( TEMP, H( ILAST, ILAST-1 ), C, S, 676 $ H( ILAST, ILAST ) ) 677 H( ILAST, ILAST-1 ) = ZERO 678 CALL DROT( ILAST-IFRSTM, H( IFRSTM, ILAST ), 1, 679 $ H( IFRSTM, ILAST-1 ), 1, C, S ) 680 CALL DROT( ILAST-IFRSTM, T( IFRSTM, ILAST ), 1, 681 $ T( IFRSTM, ILAST-1 ), 1, C, S ) 682 IF( ILZ ) 683 $ CALL DROT( N, Z( 1, ILAST ), 1, Z( 1, ILAST-1 ), 1, C, S ) 684* 685* H(ILAST,ILAST-1)=0 -- Standardize B, set ALPHAR, ALPHAI, 686* and BETA 687* 688 80 CONTINUE 689 IF( T( ILAST, ILAST ).LT.ZERO ) THEN 690 IF( ILSCHR ) THEN 691 DO 90 J = IFRSTM, ILAST 692 H( J, ILAST ) = -H( J, ILAST ) 693 T( J, ILAST ) = -T( J, ILAST ) 694 90 CONTINUE 695 ELSE 696 H( ILAST, ILAST ) = -H( ILAST, ILAST ) 697 T( ILAST, ILAST ) = -T( ILAST, ILAST ) 698 END IF 699 IF( ILZ ) THEN 700 DO 100 J = 1, N 701 Z( J, ILAST ) = -Z( J, ILAST ) 702 100 CONTINUE 703 END IF 704 END IF 705 ALPHAR( ILAST ) = H( ILAST, ILAST ) 706 ALPHAI( ILAST ) = ZERO 707 BETA( ILAST ) = T( ILAST, ILAST ) 708* 709* Go to next block -- exit if finished. 710* 711 ILAST = ILAST - 1 712 IF( ILAST.LT.ILO ) 713 $ GO TO 380 714* 715* Reset counters 716* 717 IITER = 0 718 ESHIFT = ZERO 719 IF( .NOT.ILSCHR ) THEN 720 ILASTM = ILAST 721 IF( IFRSTM.GT.ILAST ) 722 $ IFRSTM = ILO 723 END IF 724 GO TO 350 725* 726* QZ step 727* 728* This iteration only involves rows/columns IFIRST:ILAST. We 729* assume IFIRST < ILAST, and that the diagonal of B is non-zero. 730* 731 110 CONTINUE 732 IITER = IITER + 1 733 IF( .NOT.ILSCHR ) THEN 734 IFRSTM = IFIRST 735 END IF 736* 737* Compute single shifts. 738* 739* At this point, IFIRST < ILAST, and the diagonal elements of 740* T(IFIRST:ILAST,IFIRST,ILAST) are larger than BTOL (in 741* magnitude) 742* 743 IF( ( IITER / 10 )*10.EQ.IITER ) THEN 744* 745* Exceptional shift. Chosen for no particularly good reason. 746* (Single shift only.) 747* 748 IF( ( DBLE( MAXIT )*SAFMIN )*ABS( H( ILAST, ILAST-1 ) ).LT. 749 $ ABS( T( ILAST-1, ILAST-1 ) ) ) THEN 750 ESHIFT = H( ILAST, ILAST-1 ) / 751 $ T( ILAST-1, ILAST-1 ) 752 ELSE 753 ESHIFT = ESHIFT + ONE / ( SAFMIN*DBLE( MAXIT ) ) 754 END IF 755 S1 = ONE 756 WR = ESHIFT 757* 758 ELSE 759* 760* Shifts based on the generalized eigenvalues of the 761* bottom-right 2x2 block of A and B. The first eigenvalue 762* returned by DLAG2 is the Wilkinson shift (AEP p.512), 763* 764 CALL DLAG2( H( ILAST-1, ILAST-1 ), LDH, 765 $ T( ILAST-1, ILAST-1 ), LDT, SAFMIN*SAFETY, S1, 766 $ S2, WR, WR2, WI ) 767* 768 IF ( ABS( (WR/S1)*T( ILAST, ILAST ) - H( ILAST, ILAST ) ) 769 $ .GT. ABS( (WR2/S2)*T( ILAST, ILAST ) 770 $ - H( ILAST, ILAST ) ) ) THEN 771 TEMP = WR 772 WR = WR2 773 WR2 = TEMP 774 TEMP = S1 775 S1 = S2 776 S2 = TEMP 777 END IF 778 TEMP = MAX( S1, SAFMIN*MAX( ONE, ABS( WR ), ABS( WI ) ) ) 779 IF( WI.NE.ZERO ) 780 $ GO TO 200 781 END IF 782* 783* Fiddle with shift to avoid overflow 784* 785 TEMP = MIN( ASCALE, ONE )*( HALF*SAFMAX ) 786 IF( S1.GT.TEMP ) THEN 787 SCALE = TEMP / S1 788 ELSE 789 SCALE = ONE 790 END IF 791* 792 TEMP = MIN( BSCALE, ONE )*( HALF*SAFMAX ) 793 IF( ABS( WR ).GT.TEMP ) 794 $ SCALE = MIN( SCALE, TEMP / ABS( WR ) ) 795 S1 = SCALE*S1 796 WR = SCALE*WR 797* 798* Now check for two consecutive small subdiagonals. 799* 800 DO 120 J = ILAST - 1, IFIRST + 1, -1 801 ISTART = J 802 TEMP = ABS( S1*H( J, J-1 ) ) 803 TEMP2 = ABS( S1*H( J, J )-WR*T( J, J ) ) 804 TEMPR = MAX( TEMP, TEMP2 ) 805 IF( TEMPR.LT.ONE .AND. TEMPR.NE.ZERO ) THEN 806 TEMP = TEMP / TEMPR 807 TEMP2 = TEMP2 / TEMPR 808 END IF 809 IF( ABS( ( ASCALE*H( J+1, J ) )*TEMP ).LE.( ASCALE*ATOL )* 810 $ TEMP2 )GO TO 130 811 120 CONTINUE 812* 813 ISTART = IFIRST 814 130 CONTINUE 815* 816* Do an implicit single-shift QZ sweep. 817* 818* Initial Q 819* 820 TEMP = S1*H( ISTART, ISTART ) - WR*T( ISTART, ISTART ) 821 TEMP2 = S1*H( ISTART+1, ISTART ) 822 CALL DLARTG( TEMP, TEMP2, C, S, TEMPR ) 823* 824* Sweep 825* 826 DO 190 J = ISTART, ILAST - 1 827 IF( J.GT.ISTART ) THEN 828 TEMP = H( J, J-1 ) 829 CALL DLARTG( TEMP, H( J+1, J-1 ), C, S, H( J, J-1 ) ) 830 H( J+1, J-1 ) = ZERO 831 END IF 832* 833 DO 140 JC = J, ILASTM 834 TEMP = C*H( J, JC ) + S*H( J+1, JC ) 835 H( J+1, JC ) = -S*H( J, JC ) + C*H( J+1, JC ) 836 H( J, JC ) = TEMP 837 TEMP2 = C*T( J, JC ) + S*T( J+1, JC ) 838 T( J+1, JC ) = -S*T( J, JC ) + C*T( J+1, JC ) 839 T( J, JC ) = TEMP2 840 140 CONTINUE 841 IF( ILQ ) THEN 842 DO 150 JR = 1, N 843 TEMP = C*Q( JR, J ) + S*Q( JR, J+1 ) 844 Q( JR, J+1 ) = -S*Q( JR, J ) + C*Q( JR, J+1 ) 845 Q( JR, J ) = TEMP 846 150 CONTINUE 847 END IF 848* 849 TEMP = T( J+1, J+1 ) 850 CALL DLARTG( TEMP, T( J+1, J ), C, S, T( J+1, J+1 ) ) 851 T( J+1, J ) = ZERO 852* 853 DO 160 JR = IFRSTM, MIN( J+2, ILAST ) 854 TEMP = C*H( JR, J+1 ) + S*H( JR, J ) 855 H( JR, J ) = -S*H( JR, J+1 ) + C*H( JR, J ) 856 H( JR, J+1 ) = TEMP 857 160 CONTINUE 858 DO 170 JR = IFRSTM, J 859 TEMP = C*T( JR, J+1 ) + S*T( JR, J ) 860 T( JR, J ) = -S*T( JR, J+1 ) + C*T( JR, J ) 861 T( JR, J+1 ) = TEMP 862 170 CONTINUE 863 IF( ILZ ) THEN 864 DO 180 JR = 1, N 865 TEMP = C*Z( JR, J+1 ) + S*Z( JR, J ) 866 Z( JR, J ) = -S*Z( JR, J+1 ) + C*Z( JR, J ) 867 Z( JR, J+1 ) = TEMP 868 180 CONTINUE 869 END IF 870 190 CONTINUE 871* 872 GO TO 350 873* 874* Use Francis double-shift 875* 876* Note: the Francis double-shift should work with real shifts, 877* but only if the block is at least 3x3. 878* This code may break if this point is reached with 879* a 2x2 block with real eigenvalues. 880* 881 200 CONTINUE 882 IF( IFIRST+1.EQ.ILAST ) THEN 883* 884* Special case -- 2x2 block with complex eigenvectors 885* 886* Step 1: Standardize, that is, rotate so that 887* 888* ( B11 0 ) 889* B = ( ) with B11 non-negative. 890* ( 0 B22 ) 891* 892 CALL DLASV2( T( ILAST-1, ILAST-1 ), T( ILAST-1, ILAST ), 893 $ T( ILAST, ILAST ), B22, B11, SR, CR, SL, CL ) 894* 895 IF( B11.LT.ZERO ) THEN 896 CR = -CR 897 SR = -SR 898 B11 = -B11 899 B22 = -B22 900 END IF 901* 902 CALL DROT( ILASTM+1-IFIRST, H( ILAST-1, ILAST-1 ), LDH, 903 $ H( ILAST, ILAST-1 ), LDH, CL, SL ) 904 CALL DROT( ILAST+1-IFRSTM, H( IFRSTM, ILAST-1 ), 1, 905 $ H( IFRSTM, ILAST ), 1, CR, SR ) 906* 907 IF( ILAST.LT.ILASTM ) 908 $ CALL DROT( ILASTM-ILAST, T( ILAST-1, ILAST+1 ), LDT, 909 $ T( ILAST, ILAST+1 ), LDT, CL, SL ) 910 IF( IFRSTM.LT.ILAST-1 ) 911 $ CALL DROT( IFIRST-IFRSTM, T( IFRSTM, ILAST-1 ), 1, 912 $ T( IFRSTM, ILAST ), 1, CR, SR ) 913* 914 IF( ILQ ) 915 $ CALL DROT( N, Q( 1, ILAST-1 ), 1, Q( 1, ILAST ), 1, CL, 916 $ SL ) 917 IF( ILZ ) 918 $ CALL DROT( N, Z( 1, ILAST-1 ), 1, Z( 1, ILAST ), 1, CR, 919 $ SR ) 920* 921 T( ILAST-1, ILAST-1 ) = B11 922 T( ILAST-1, ILAST ) = ZERO 923 T( ILAST, ILAST-1 ) = ZERO 924 T( ILAST, ILAST ) = B22 925* 926* If B22 is negative, negate column ILAST 927* 928 IF( B22.LT.ZERO ) THEN 929 DO 210 J = IFRSTM, ILAST 930 H( J, ILAST ) = -H( J, ILAST ) 931 T( J, ILAST ) = -T( J, ILAST ) 932 210 CONTINUE 933* 934 IF( ILZ ) THEN 935 DO 220 J = 1, N 936 Z( J, ILAST ) = -Z( J, ILAST ) 937 220 CONTINUE 938 END IF 939 B22 = -B22 940 END IF 941* 942* Step 2: Compute ALPHAR, ALPHAI, and BETA (see refs.) 943* 944* Recompute shift 945* 946 CALL DLAG2( H( ILAST-1, ILAST-1 ), LDH, 947 $ T( ILAST-1, ILAST-1 ), LDT, SAFMIN*SAFETY, S1, 948 $ TEMP, WR, TEMP2, WI ) 949* 950* If standardization has perturbed the shift onto real line, 951* do another (real single-shift) QR step. 952* 953 IF( WI.EQ.ZERO ) 954 $ GO TO 350 955 S1INV = ONE / S1 956* 957* Do EISPACK (QZVAL) computation of alpha and beta 958* 959 A11 = H( ILAST-1, ILAST-1 ) 960 A21 = H( ILAST, ILAST-1 ) 961 A12 = H( ILAST-1, ILAST ) 962 A22 = H( ILAST, ILAST ) 963* 964* Compute complex Givens rotation on right 965* (Assume some element of C = (sA - wB) > unfl ) 966* __ 967* (sA - wB) ( CZ -SZ ) 968* ( SZ CZ ) 969* 970 C11R = S1*A11 - WR*B11 971 C11I = -WI*B11 972 C12 = S1*A12 973 C21 = S1*A21 974 C22R = S1*A22 - WR*B22 975 C22I = -WI*B22 976* 977 IF( ABS( C11R )+ABS( C11I )+ABS( C12 ).GT.ABS( C21 )+ 978 $ ABS( C22R )+ABS( C22I ) ) THEN 979 T1 = DLAPY3( C12, C11R, C11I ) 980 CZ = C12 / T1 981 SZR = -C11R / T1 982 SZI = -C11I / T1 983 ELSE 984 CZ = DLAPY2( C22R, C22I ) 985 IF( CZ.LE.SAFMIN ) THEN 986 CZ = ZERO 987 SZR = ONE 988 SZI = ZERO 989 ELSE 990 TEMPR = C22R / CZ 991 TEMPI = C22I / CZ 992 T1 = DLAPY2( CZ, C21 ) 993 CZ = CZ / T1 994 SZR = -C21*TEMPR / T1 995 SZI = C21*TEMPI / T1 996 END IF 997 END IF 998* 999* Compute Givens rotation on left 1000* 1001* ( CQ SQ ) 1002* ( __ ) A or B 1003* ( -SQ CQ ) 1004* 1005 AN = ABS( A11 ) + ABS( A12 ) + ABS( A21 ) + ABS( A22 ) 1006 BN = ABS( B11 ) + ABS( B22 ) 1007 WABS = ABS( WR ) + ABS( WI ) 1008 IF( S1*AN.GT.WABS*BN ) THEN 1009 CQ = CZ*B11 1010 SQR = SZR*B22 1011 SQI = -SZI*B22 1012 ELSE 1013 A1R = CZ*A11 + SZR*A12 1014 A1I = SZI*A12 1015 A2R = CZ*A21 + SZR*A22 1016 A2I = SZI*A22 1017 CQ = DLAPY2( A1R, A1I ) 1018 IF( CQ.LE.SAFMIN ) THEN 1019 CQ = ZERO 1020 SQR = ONE 1021 SQI = ZERO 1022 ELSE 1023 TEMPR = A1R / CQ 1024 TEMPI = A1I / CQ 1025 SQR = TEMPR*A2R + TEMPI*A2I 1026 SQI = TEMPI*A2R - TEMPR*A2I 1027 END IF 1028 END IF 1029 T1 = DLAPY3( CQ, SQR, SQI ) 1030 CQ = CQ / T1 1031 SQR = SQR / T1 1032 SQI = SQI / T1 1033* 1034* Compute diagonal elements of QBZ 1035* 1036 TEMPR = SQR*SZR - SQI*SZI 1037 TEMPI = SQR*SZI + SQI*SZR 1038 B1R = CQ*CZ*B11 + TEMPR*B22 1039 B1I = TEMPI*B22 1040 B1A = DLAPY2( B1R, B1I ) 1041 B2R = CQ*CZ*B22 + TEMPR*B11 1042 B2I = -TEMPI*B11 1043 B2A = DLAPY2( B2R, B2I ) 1044* 1045* Normalize so beta > 0, and Im( alpha1 ) > 0 1046* 1047 BETA( ILAST-1 ) = B1A 1048 BETA( ILAST ) = B2A 1049 ALPHAR( ILAST-1 ) = ( WR*B1A )*S1INV 1050 ALPHAI( ILAST-1 ) = ( WI*B1A )*S1INV 1051 ALPHAR( ILAST ) = ( WR*B2A )*S1INV 1052 ALPHAI( ILAST ) = -( WI*B2A )*S1INV 1053* 1054* Step 3: Go to next block -- exit if finished. 1055* 1056 ILAST = IFIRST - 1 1057 IF( ILAST.LT.ILO ) 1058 $ GO TO 380 1059* 1060* Reset counters 1061* 1062 IITER = 0 1063 ESHIFT = ZERO 1064 IF( .NOT.ILSCHR ) THEN 1065 ILASTM = ILAST 1066 IF( IFRSTM.GT.ILAST ) 1067 $ IFRSTM = ILO 1068 END IF 1069 GO TO 350 1070 ELSE 1071* 1072* Usual case: 3x3 or larger block, using Francis implicit 1073* double-shift 1074* 1075* 2 1076* Eigenvalue equation is w - c w + d = 0, 1077* 1078* -1 2 -1 1079* so compute 1st column of (A B ) - c A B + d 1080* using the formula in QZIT (from EISPACK) 1081* 1082* We assume that the block is at least 3x3 1083* 1084 AD11 = ( ASCALE*H( ILAST-1, ILAST-1 ) ) / 1085 $ ( BSCALE*T( ILAST-1, ILAST-1 ) ) 1086 AD21 = ( ASCALE*H( ILAST, ILAST-1 ) ) / 1087 $ ( BSCALE*T( ILAST-1, ILAST-1 ) ) 1088 AD12 = ( ASCALE*H( ILAST-1, ILAST ) ) / 1089 $ ( BSCALE*T( ILAST, ILAST ) ) 1090 AD22 = ( ASCALE*H( ILAST, ILAST ) ) / 1091 $ ( BSCALE*T( ILAST, ILAST ) ) 1092 U12 = T( ILAST-1, ILAST ) / T( ILAST, ILAST ) 1093 AD11L = ( ASCALE*H( IFIRST, IFIRST ) ) / 1094 $ ( BSCALE*T( IFIRST, IFIRST ) ) 1095 AD21L = ( ASCALE*H( IFIRST+1, IFIRST ) ) / 1096 $ ( BSCALE*T( IFIRST, IFIRST ) ) 1097 AD12L = ( ASCALE*H( IFIRST, IFIRST+1 ) ) / 1098 $ ( BSCALE*T( IFIRST+1, IFIRST+1 ) ) 1099 AD22L = ( ASCALE*H( IFIRST+1, IFIRST+1 ) ) / 1100 $ ( BSCALE*T( IFIRST+1, IFIRST+1 ) ) 1101 AD32L = ( ASCALE*H( IFIRST+2, IFIRST+1 ) ) / 1102 $ ( BSCALE*T( IFIRST+1, IFIRST+1 ) ) 1103 U12L = T( IFIRST, IFIRST+1 ) / T( IFIRST+1, IFIRST+1 ) 1104* 1105 V( 1 ) = ( AD11-AD11L )*( AD22-AD11L ) - AD12*AD21 + 1106 $ AD21*U12*AD11L + ( AD12L-AD11L*U12L )*AD21L 1107 V( 2 ) = ( ( AD22L-AD11L )-AD21L*U12L-( AD11-AD11L )- 1108 $ ( AD22-AD11L )+AD21*U12 )*AD21L 1109 V( 3 ) = AD32L*AD21L 1110* 1111 ISTART = IFIRST 1112* 1113 CALL DLARFG( 3, V( 1 ), V( 2 ), 1, TAU ) 1114 V( 1 ) = ONE 1115* 1116* Sweep 1117* 1118 DO 290 J = ISTART, ILAST - 2 1119* 1120* All but last elements: use 3x3 Householder transforms. 1121* 1122* Zero (j-1)st column of A 1123* 1124 IF( J.GT.ISTART ) THEN 1125 V( 1 ) = H( J, J-1 ) 1126 V( 2 ) = H( J+1, J-1 ) 1127 V( 3 ) = H( J+2, J-1 ) 1128* 1129 CALL DLARFG( 3, H( J, J-1 ), V( 2 ), 1, TAU ) 1130 V( 1 ) = ONE 1131 H( J+1, J-1 ) = ZERO 1132 H( J+2, J-1 ) = ZERO 1133 END IF 1134* 1135 DO 230 JC = J, ILASTM 1136 TEMP = TAU*( H( J, JC )+V( 2 )*H( J+1, JC )+V( 3 )* 1137 $ H( J+2, JC ) ) 1138 H( J, JC ) = H( J, JC ) - TEMP 1139 H( J+1, JC ) = H( J+1, JC ) - TEMP*V( 2 ) 1140 H( J+2, JC ) = H( J+2, JC ) - TEMP*V( 3 ) 1141 TEMP2 = TAU*( T( J, JC )+V( 2 )*T( J+1, JC )+V( 3 )* 1142 $ T( J+2, JC ) ) 1143 T( J, JC ) = T( J, JC ) - TEMP2 1144 T( J+1, JC ) = T( J+1, JC ) - TEMP2*V( 2 ) 1145 T( J+2, JC ) = T( J+2, JC ) - TEMP2*V( 3 ) 1146 230 CONTINUE 1147 IF( ILQ ) THEN 1148 DO 240 JR = 1, N 1149 TEMP = TAU*( Q( JR, J )+V( 2 )*Q( JR, J+1 )+V( 3 )* 1150 $ Q( JR, J+2 ) ) 1151 Q( JR, J ) = Q( JR, J ) - TEMP 1152 Q( JR, J+1 ) = Q( JR, J+1 ) - TEMP*V( 2 ) 1153 Q( JR, J+2 ) = Q( JR, J+2 ) - TEMP*V( 3 ) 1154 240 CONTINUE 1155 END IF 1156* 1157* Zero j-th column of B (see DLAGBC for details) 1158* 1159* Swap rows to pivot 1160* 1161 ILPIVT = .FALSE. 1162 TEMP = MAX( ABS( T( J+1, J+1 ) ), ABS( T( J+1, J+2 ) ) ) 1163 TEMP2 = MAX( ABS( T( J+2, J+1 ) ), ABS( T( J+2, J+2 ) ) ) 1164 IF( MAX( TEMP, TEMP2 ).LT.SAFMIN ) THEN 1165 SCALE = ZERO 1166 U1 = ONE 1167 U2 = ZERO 1168 GO TO 250 1169 ELSE IF( TEMP.GE.TEMP2 ) THEN 1170 W11 = T( J+1, J+1 ) 1171 W21 = T( J+2, J+1 ) 1172 W12 = T( J+1, J+2 ) 1173 W22 = T( J+2, J+2 ) 1174 U1 = T( J+1, J ) 1175 U2 = T( J+2, J ) 1176 ELSE 1177 W21 = T( J+1, J+1 ) 1178 W11 = T( J+2, J+1 ) 1179 W22 = T( J+1, J+2 ) 1180 W12 = T( J+2, J+2 ) 1181 U2 = T( J+1, J ) 1182 U1 = T( J+2, J ) 1183 END IF 1184* 1185* Swap columns if nec. 1186* 1187 IF( ABS( W12 ).GT.ABS( W11 ) ) THEN 1188 ILPIVT = .TRUE. 1189 TEMP = W12 1190 TEMP2 = W22 1191 W12 = W11 1192 W22 = W21 1193 W11 = TEMP 1194 W21 = TEMP2 1195 END IF 1196* 1197* LU-factor 1198* 1199 TEMP = W21 / W11 1200 U2 = U2 - TEMP*U1 1201 W22 = W22 - TEMP*W12 1202 W21 = ZERO 1203* 1204* Compute SCALE 1205* 1206 SCALE = ONE 1207 IF( ABS( W22 ).LT.SAFMIN ) THEN 1208 SCALE = ZERO 1209 U2 = ONE 1210 U1 = -W12 / W11 1211 GO TO 250 1212 END IF 1213 IF( ABS( W22 ).LT.ABS( U2 ) ) 1214 $ SCALE = ABS( W22 / U2 ) 1215 IF( ABS( W11 ).LT.ABS( U1 ) ) 1216 $ SCALE = MIN( SCALE, ABS( W11 / U1 ) ) 1217* 1218* Solve 1219* 1220 U2 = ( SCALE*U2 ) / W22 1221 U1 = ( SCALE*U1-W12*U2 ) / W11 1222* 1223 250 CONTINUE 1224 IF( ILPIVT ) THEN 1225 TEMP = U2 1226 U2 = U1 1227 U1 = TEMP 1228 END IF 1229* 1230* Compute Householder Vector 1231* 1232 T1 = SQRT( SCALE**2+U1**2+U2**2 ) 1233 TAU = ONE + SCALE / T1 1234 VS = -ONE / ( SCALE+T1 ) 1235 V( 1 ) = ONE 1236 V( 2 ) = VS*U1 1237 V( 3 ) = VS*U2 1238* 1239* Apply transformations from the right. 1240* 1241 DO 260 JR = IFRSTM, MIN( J+3, ILAST ) 1242 TEMP = TAU*( H( JR, J )+V( 2 )*H( JR, J+1 )+V( 3 )* 1243 $ H( JR, J+2 ) ) 1244 H( JR, J ) = H( JR, J ) - TEMP 1245 H( JR, J+1 ) = H( JR, J+1 ) - TEMP*V( 2 ) 1246 H( JR, J+2 ) = H( JR, J+2 ) - TEMP*V( 3 ) 1247 260 CONTINUE 1248 DO 270 JR = IFRSTM, J + 2 1249 TEMP = TAU*( T( JR, J )+V( 2 )*T( JR, J+1 )+V( 3 )* 1250 $ T( JR, J+2 ) ) 1251 T( JR, J ) = T( JR, J ) - TEMP 1252 T( JR, J+1 ) = T( JR, J+1 ) - TEMP*V( 2 ) 1253 T( JR, J+2 ) = T( JR, J+2 ) - TEMP*V( 3 ) 1254 270 CONTINUE 1255 IF( ILZ ) THEN 1256 DO 280 JR = 1, N 1257 TEMP = TAU*( Z( JR, J )+V( 2 )*Z( JR, J+1 )+V( 3 )* 1258 $ Z( JR, J+2 ) ) 1259 Z( JR, J ) = Z( JR, J ) - TEMP 1260 Z( JR, J+1 ) = Z( JR, J+1 ) - TEMP*V( 2 ) 1261 Z( JR, J+2 ) = Z( JR, J+2 ) - TEMP*V( 3 ) 1262 280 CONTINUE 1263 END IF 1264 T( J+1, J ) = ZERO 1265 T( J+2, J ) = ZERO 1266 290 CONTINUE 1267* 1268* Last elements: Use Givens rotations 1269* 1270* Rotations from the left 1271* 1272 J = ILAST - 1 1273 TEMP = H( J, J-1 ) 1274 CALL DLARTG( TEMP, H( J+1, J-1 ), C, S, H( J, J-1 ) ) 1275 H( J+1, J-1 ) = ZERO 1276* 1277 DO 300 JC = J, ILASTM 1278 TEMP = C*H( J, JC ) + S*H( J+1, JC ) 1279 H( J+1, JC ) = -S*H( J, JC ) + C*H( J+1, JC ) 1280 H( J, JC ) = TEMP 1281 TEMP2 = C*T( J, JC ) + S*T( J+1, JC ) 1282 T( J+1, JC ) = -S*T( J, JC ) + C*T( J+1, JC ) 1283 T( J, JC ) = TEMP2 1284 300 CONTINUE 1285 IF( ILQ ) THEN 1286 DO 310 JR = 1, N 1287 TEMP = C*Q( JR, J ) + S*Q( JR, J+1 ) 1288 Q( JR, J+1 ) = -S*Q( JR, J ) + C*Q( JR, J+1 ) 1289 Q( JR, J ) = TEMP 1290 310 CONTINUE 1291 END IF 1292* 1293* Rotations from the right. 1294* 1295 TEMP = T( J+1, J+1 ) 1296 CALL DLARTG( TEMP, T( J+1, J ), C, S, T( J+1, J+1 ) ) 1297 T( J+1, J ) = ZERO 1298* 1299 DO 320 JR = IFRSTM, ILAST 1300 TEMP = C*H( JR, J+1 ) + S*H( JR, J ) 1301 H( JR, J ) = -S*H( JR, J+1 ) + C*H( JR, J ) 1302 H( JR, J+1 ) = TEMP 1303 320 CONTINUE 1304 DO 330 JR = IFRSTM, ILAST - 1 1305 TEMP = C*T( JR, J+1 ) + S*T( JR, J ) 1306 T( JR, J ) = -S*T( JR, J+1 ) + C*T( JR, J ) 1307 T( JR, J+1 ) = TEMP 1308 330 CONTINUE 1309 IF( ILZ ) THEN 1310 DO 340 JR = 1, N 1311 TEMP = C*Z( JR, J+1 ) + S*Z( JR, J ) 1312 Z( JR, J ) = -S*Z( JR, J+1 ) + C*Z( JR, J ) 1313 Z( JR, J+1 ) = TEMP 1314 340 CONTINUE 1315 END IF 1316* 1317* End of Double-Shift code 1318* 1319 END IF 1320* 1321 GO TO 350 1322* 1323* End of iteration loop 1324* 1325 350 CONTINUE 1326 360 CONTINUE 1327* 1328* Drop-through = non-convergence 1329* 1330 INFO = ILAST 1331 GO TO 420 1332* 1333* Successful completion of all QZ steps 1334* 1335 380 CONTINUE 1336* 1337* Set Eigenvalues 1:ILO-1 1338* 1339 DO 410 J = 1, ILO - 1 1340 IF( T( J, J ).LT.ZERO ) THEN 1341 IF( ILSCHR ) THEN 1342 DO 390 JR = 1, J 1343 H( JR, J ) = -H( JR, J ) 1344 T( JR, J ) = -T( JR, J ) 1345 390 CONTINUE 1346 ELSE 1347 H( J, J ) = -H( J, J ) 1348 T( J, J ) = -T( J, J ) 1349 END IF 1350 IF( ILZ ) THEN 1351 DO 400 JR = 1, N 1352 Z( JR, J ) = -Z( JR, J ) 1353 400 CONTINUE 1354 END IF 1355 END IF 1356 ALPHAR( J ) = H( J, J ) 1357 ALPHAI( J ) = ZERO 1358 BETA( J ) = T( J, J ) 1359 410 CONTINUE 1360* 1361* Normal Termination 1362* 1363 INFO = 0 1364* 1365* Exit (other than argument error) -- return optimal workspace size 1366* 1367 420 CONTINUE 1368 WORK( 1 ) = DBLE( N ) 1369 RETURN 1370* 1371* End of DHGEQZ 1372* 1373 END 1374