1       SUBROUTINE ZHGEQZ( JOB, COMPQ, COMPZ, N, ILO, IHI, H, LDH, T, LDT,
  2      $                   ALPHA, BETA, Q, LDQ, Z, LDZ, WORK, LWORK,
  3      $                   RWORK, INFO )
  4 *
  5 *  -- LAPACK routine (version 3.2) --
  6 *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
  7 *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  8 *     November 2006
  9 *
 10 *     .. Scalar Arguments ..
 11       CHARACTER          COMPQ, COMPZ, JOB
 12       INTEGER            IHI, ILO, INFO, LDH, LDQ, LDT, LDZ, LWORK, N
 13 *     ..
 14 *     .. Array Arguments ..
 15       DOUBLE PRECISION   RWORK( * )
 16       COMPLEX*16         ALPHA( * ), BETA( * ), H( LDH, * ),
 17      $                   Q( LDQ, * ), T( LDT, * ), WORK( * ),
 18      $                   Z( LDZ, * )
 19 *     ..
 20 *
 21 *  Purpose
 22 *  =======
 23 *
 24 *  ZHGEQZ computes the eigenvalues of a complex matrix pair (H,T),
 25 *  where H is an upper Hessenberg matrix and T is upper triangular,
 26 *  using the single-shift QZ method.
 27 *  Matrix pairs of this type are produced by the reduction to
 28 *  generalized upper Hessenberg form of a complex matrix pair (A,B):
 29 *  
 30 *     A = Q1*H*Z1**H,  B = Q1*T*Z1**H,
 31 *  
 32 *  as computed by ZGGHRD.
 33 *  
 34 *  If JOB='S', then the Hessenberg-triangular pair (H,T) is
 35 *  also reduced to generalized Schur form,
 36 *  
 37 *     H = Q*S*Z**H,  T = Q*P*Z**H,
 38 *  
 39 *  where Q and Z are unitary matrices and S and P are upper triangular.
 40 *  
 41 *  Optionally, the unitary matrix Q from the generalized Schur
 42 *  factorization may be postmultiplied into an input matrix Q1, and the
 43 *  unitary matrix Z may be postmultiplied into an input matrix Z1.
 44 *  If Q1 and Z1 are the unitary matrices from ZGGHRD that reduced
 45 *  the matrix pair (A,B) to generalized Hessenberg form, then the output
 46 *  matrices Q1*Q and Z1*Z are the unitary factors from the generalized
 47 *  Schur factorization of (A,B):
 48 *  
 49 *     A = (Q1*Q)*S*(Z1*Z)**H,  B = (Q1*Q)*P*(Z1*Z)**H.
 50 *  
 51 *  To avoid overflow, eigenvalues of the matrix pair (H,T)
 52 *  (equivalently, of (A,B)) are computed as a pair of complex values
 53 *  (alpha,beta).  If beta is nonzero, lambda = alpha / beta is an
 54 *  eigenvalue of the generalized nonsymmetric eigenvalue problem (GNEP)
 55 *     A*x = lambda*B*x
 56 *  and if alpha is nonzero, mu = beta / alpha is an eigenvalue of the
 57 *  alternate form of the GNEP
 58 *     mu*A*y = B*y.
 59 *  The values of alpha and beta for the i-th eigenvalue can be read
 60 *  directly from the generalized Schur form:  alpha = S(i,i),
 61 *  beta = P(i,i).
 62 *
 63 *  Ref: C.B. Moler & G.W. Stewart, "An Algorithm for Generalized Matrix
 64 *       Eigenvalue Problems", SIAM J. Numer. Anal., 10(1973),
 65 *       pp. 241--256.
 66 *
 67 *  Arguments
 68 *  =========
 69 *
 70 *  JOB     (input) CHARACTER*1
 71 *          = 'E': Compute eigenvalues only;
 72 *          = 'S': Computer eigenvalues and the Schur form.
 73 *
 74 *  COMPQ   (input) CHARACTER*1
 75 *          = 'N': Left Schur vectors (Q) are not computed;
 76 *          = 'I': Q is initialized to the unit matrix and the matrix Q
 77 *                 of left Schur vectors of (H,T) is returned;
 78 *          = 'V': Q must contain a unitary matrix Q1 on entry and
 79 *                 the product Q1*Q is returned.
 80 *
 81 *  COMPZ   (input) CHARACTER*1
 82 *          = 'N': Right Schur vectors (Z) are not computed;
 83 *          = 'I': Q is initialized to the unit matrix and the matrix Z
 84 *                 of right Schur vectors of (H,T) is returned;
 85 *          = 'V': Z must contain a unitary matrix Z1 on entry and
 86 *                 the product Z1*Z is returned.
 87 *
 88 *  N       (input) INTEGER
 89 *          The order of the matrices H, T, Q, and Z.  N >= 0.
 90 *
 91 *  ILO     (input) INTEGER
 92 *  IHI     (input) INTEGER
 93 *          ILO and IHI mark the rows and columns of H which are in
 94 *          Hessenberg form.  It is assumed that A is already upper
 95 *          triangular in rows and columns 1:ILO-1 and IHI+1:N.
 96 *          If N > 0, 1 <= ILO <= IHI <= N; if N = 0, ILO=1 and IHI=0.
 97 *
 98 *  H       (input/output) COMPLEX*16 array, dimension (LDH, N)
 99 *          On entry, the N-by-N upper Hessenberg matrix H.
100 *          On exit, if JOB = 'S', H contains the upper triangular
101 *          matrix S from the generalized Schur factorization.
102 *          If JOB = 'E', the diagonal of H matches that of S, but
103 *          the rest of H is unspecified.
104 *
105 *  LDH     (input) INTEGER
106 *          The leading dimension of the array H.  LDH >= max( 1, N ).
107 *
108 *  T       (input/output) COMPLEX*16 array, dimension (LDT, N)
109 *          On entry, the N-by-N upper triangular matrix T.
110 *          On exit, if JOB = 'S', T contains the upper triangular
111 *          matrix P from the generalized Schur factorization.
112 *          If JOB = 'E', the diagonal of T matches that of P, but
113 *          the rest of T is unspecified.
114 *
115 *  LDT     (input) INTEGER
116 *          The leading dimension of the array T.  LDT >= max( 1, N ).
117 *
118 *  ALPHA   (output) COMPLEX*16 array, dimension (N)
119 *          The complex scalars alpha that define the eigenvalues of
120 *          GNEP.  ALPHA(i) = S(i,i) in the generalized Schur
121 *          factorization.
122 *
123 *  BETA    (output) COMPLEX*16 array, dimension (N)
124 *          The real non-negative scalars beta that define the
125 *          eigenvalues of GNEP.  BETA(i) = P(i,i) in the generalized
126 *          Schur factorization.
127 *
128 *          Together, the quantities alpha = ALPHA(j) and beta = BETA(j)
129 *          represent the j-th eigenvalue of the matrix pair (A,B), in
130 *          one of the forms lambda = alpha/beta or mu = beta/alpha.
131 *          Since either lambda or mu may overflow, they should not,
132 *          in general, be computed.
133 *
134 *  Q       (input/output) COMPLEX*16 array, dimension (LDQ, N)
135 *          On entry, if COMPZ = 'V', the unitary matrix Q1 used in the
136 *          reduction of (A,B) to generalized Hessenberg form.
137 *          On exit, if COMPZ = 'I', the unitary matrix of left Schur
138 *          vectors of (H,T), and if COMPZ = 'V', the unitary matrix of
139 *          left Schur vectors of (A,B).
140 *          Not referenced if COMPZ = 'N'.
141 *
142 *  LDQ     (input) INTEGER
143 *          The leading dimension of the array Q.  LDQ >= 1.
144 *          If COMPQ='V' or 'I', then LDQ >= N.
145 *
146 *  Z       (input/output) COMPLEX*16 array, dimension (LDZ, N)
147 *          On entry, if COMPZ = 'V', the unitary matrix Z1 used in the
148 *          reduction of (A,B) to generalized Hessenberg form.
149 *          On exit, if COMPZ = 'I', the unitary matrix of right Schur
150 *          vectors of (H,T), and if COMPZ = 'V', the unitary matrix of
151 *          right Schur vectors of (A,B).
152 *          Not referenced if COMPZ = 'N'.
153 *
154 *  LDZ     (input) INTEGER
155 *          The leading dimension of the array Z.  LDZ >= 1.
156 *          If COMPZ='V' or 'I', then LDZ >= N.
157 *
158 *  WORK    (workspace/output) COMPLEX*16 array, dimension (MAX(1,LWORK))
159 *          On exit, if INFO >= 0, WORK(1) returns the optimal LWORK.
160 *
161 *  LWORK   (input) INTEGER
162 *          The dimension of the array WORK.  LWORK >= max(1,N).
163 *
164 *          If LWORK = -1, then a workspace query is assumed; the routine
165 *          only calculates the optimal size of the WORK array, returns
166 *          this value as the first entry of the WORK array, and no error
167 *          message related to LWORK is issued by XERBLA.
168 *
169 *  RWORK   (workspace) DOUBLE PRECISION array, dimension (N)
170 *
171 *  INFO    (output) INTEGER
172 *          = 0: successful exit
173 *          < 0: if INFO = -i, the i-th argument had an illegal value
174 *          = 1,...,N: the QZ iteration did not converge.  (H,T) is not
175 *                     in Schur form, but ALPHA(i) and BETA(i),
176 *                     i=INFO+1,...,N should be correct.
177 *          = N+1,...,2*N: the shift calculation failed.  (H,T) is not
178 *                     in Schur form, but ALPHA(i) and BETA(i),
179 *                     i=INFO-N+1,...,N should be correct.
180 *
181 *  Further Details
182 *  ===============
183 *
184 *  We assume that complex ABS works as long as its value is less than
185 *  overflow.
186 *
187 *  =====================================================================
188 *
189 *     .. Parameters ..
190       COMPLEX*16         CZERO, CONE
191       PARAMETER          ( CZERO = ( 0.0D+00.0D+0 ),
192      $                   CONE = ( 1.0D+00.0D+0 ) )
193       DOUBLE PRECISION   ZERO, ONE
194       PARAMETER          ( ZERO = 0.0D+0, ONE = 1.0D+0 )
195       DOUBLE PRECISION   HALF
196       PARAMETER          ( HALF = 0.5D+0 )
197 *     ..
198 *     .. Local Scalars ..
199       LOGICAL            ILAZR2, ILAZRO, ILQ, ILSCHR, ILZ, LQUERY
200       INTEGER            ICOMPQ, ICOMPZ, IFIRST, IFRSTM, IITER, ILAST,
201      $                   ILASTM, IN, ISCHUR, ISTART, J, JC, JCH, JITER,
202      $                   JR, MAXIT
203       DOUBLE PRECISION   ABSB, ANORM, ASCALE, ATOL, BNORM, BSCALE, BTOL,
204      $                   C, SAFMIN, TEMP, TEMP2, TEMPR, ULP
205       COMPLEX*16         ABI22, AD11, AD12, AD21, AD22, CTEMP, CTEMP2,
206      $                   CTEMP3, ESHIFT, RTDISC, S, SHIFT, SIGNBC, T1,
207      $                   U12, X
208 *     ..
209 *     .. External Functions ..
210       LOGICAL            LSAME
211       DOUBLE PRECISION   DLAMCH, ZLANHS
212       EXTERNAL           LSAME, DLAMCH, ZLANHS
213 *     ..
214 *     .. External Subroutines ..
215       EXTERNAL           XERBLA, ZLARTG, ZLASET, ZROT, ZSCAL
216 *     ..
217 *     .. Intrinsic Functions ..
218       INTRINSIC          ABSDBLEDCMPLXDCONJGDIMAGMAXMIN,
219      $                   SQRT
220 *     ..
221 *     .. Statement Functions ..
222       DOUBLE PRECISION   ABS1
223 *     ..
224 *     .. Statement Function definitions ..
225       ABS1( X ) = ABSDBLE( X ) ) + ABSDIMAG( X ) )
226 *     ..
227 *     .. Executable Statements ..
228 *
229 *     Decode JOB, COMPQ, COMPZ
230 *
231       IF( LSAME( JOB, 'E' ) ) THEN
232          ILSCHR = .FALSE.
233          ISCHUR = 1
234       ELSE IF( LSAME( JOB, 'S' ) ) THEN
235          ILSCHR = .TRUE.
236          ISCHUR = 2
237       ELSE
238          ISCHUR = 0
239       END IF
240 *
241       IF( LSAME( COMPQ, 'N' ) ) THEN
242          ILQ = .FALSE.
243          ICOMPQ = 1
244       ELSE IF( LSAME( COMPQ, 'V' ) ) THEN
245          ILQ = .TRUE.
246          ICOMPQ = 2
247       ELSE IF( LSAME( COMPQ, 'I' ) ) THEN
248          ILQ = .TRUE.
249          ICOMPQ = 3
250       ELSE
251          ICOMPQ = 0
252       END IF
253 *
254       IF( LSAME( COMPZ, 'N' ) ) THEN
255          ILZ = .FALSE.
256          ICOMPZ = 1
257       ELSE IF( LSAME( COMPZ, 'V' ) ) THEN
258          ILZ = .TRUE.
259          ICOMPZ = 2
260       ELSE IF( LSAME( COMPZ, 'I' ) ) THEN
261          ILZ = .TRUE.
262          ICOMPZ = 3
263       ELSE
264          ICOMPZ = 0
265       END IF
266 *
267 *     Check Argument Values
268 *
269       INFO = 0
270       WORK( 1 ) = MAX1, N )
271       LQUERY = ( LWORK.EQ.-1 )
272       IF( ISCHUR.EQ.0 ) THEN
273          INFO = -1
274       ELSE IF( ICOMPQ.EQ.0 ) THEN
275          INFO = -2
276       ELSE IF( ICOMPZ.EQ.0 ) THEN
277          INFO = -3
278       ELSE IF( N.LT.0 ) THEN
279          INFO = -4
280       ELSE IF( ILO.LT.1 ) THEN
281          INFO = -5
282       ELSE IF( IHI.GT..OR. IHI.LT.ILO-1 ) THEN
283          INFO = -6
284       ELSE IF( LDH.LT.N ) THEN
285          INFO = -8
286       ELSE IF( LDT.LT.N ) THEN
287          INFO = -10
288       ELSE IF( LDQ.LT.1 .OR. ( ILQ .AND. LDQ.LT.N ) ) THEN
289          INFO = -14
290       ELSE IF( LDZ.LT.1 .OR. ( ILZ .AND. LDZ.LT.N ) ) THEN
291          INFO = -16
292       ELSE IF( LWORK.LT.MAX1, N ) .AND. .NOT.LQUERY ) THEN
293          INFO = -18
294       END IF
295       IF( INFO.NE.0 ) THEN
296          CALL XERBLA( 'ZHGEQZ'-INFO )
297          RETURN
298       ELSE IF( LQUERY ) THEN
299          RETURN
300       END IF
301 *
302 *     Quick return if possible
303 *
304 *     WORK( 1 ) = CMPLX( 1 )
305       IF( N.LE.0 ) THEN
306          WORK( 1 ) = DCMPLX1 )
307          RETURN
308       END IF
309 *
310 *     Initialize Q and Z
311 *
312       IF( ICOMPQ.EQ.3 )
313      $   CALL ZLASET( 'Full', N, N, CZERO, CONE, Q, LDQ )
314       IF( ICOMPZ.EQ.3 )
315      $   CALL ZLASET( 'Full', N, N, CZERO, CONE, Z, LDZ )
316 *
317 *     Machine Constants
318 *
319       IN = IHI + 1 - ILO
320       SAFMIN = DLAMCH( 'S' )
321       ULP = DLAMCH( 'E' )*DLAMCH( 'B' )
322       ANORM = ZLANHS( 'F'IN, H( ILO, ILO ), LDH, RWORK )
323       BNORM = ZLANHS( 'F'IN, T( ILO, ILO ), LDT, RWORK )
324       ATOL = MAX( SAFMIN, ULP*ANORM )
325       BTOL = MAX( SAFMIN, ULP*BNORM )
326       ASCALE = ONE / MAX( SAFMIN, ANORM )
327       BSCALE = ONE / MAX( SAFMIN, BNORM )
328 *
329 *
330 *     Set Eigenvalues IHI+1:N
331 *
332       DO 10 J = IHI + 1, N
333          ABSB = ABS( T( J, J ) )
334          IF( ABSB.GT.SAFMIN ) THEN
335             SIGNBC = DCONJG( T( J, J ) / ABSB )
336             T( J, J ) = ABSB
337             IF( ILSCHR ) THEN
338                CALL ZSCAL( J-1, SIGNBC, T( 1, J ), 1 )
339                CALL ZSCAL( J, SIGNBC, H( 1, J ), 1 )
340             ELSE
341                H( J, J ) = H( J, J )*SIGNBC
342             END IF
343             IF( ILZ )
344      $         CALL ZSCAL( N, SIGNBC, Z( 1, J ), 1 )
345          ELSE
346             T( J, J ) = CZERO
347          END IF
348          ALPHA( J ) = H( J, J )
349          BETA( J ) = T( J, J )
350    10 CONTINUE
351 *
352 *     If IHI < ILO, skip QZ steps
353 *
354       IF( IHI.LT.ILO )
355      $   GO TO 190
356 *
357 *     MAIN QZ ITERATION LOOP
358 *
359 *     Initialize dynamic indices
360 *
361 *     Eigenvalues ILAST+1:N have been found.
362 *        Column operations modify rows IFRSTM:whatever
363 *        Row operations modify columns whatever:ILASTM
364 *
365 *     If only eigenvalues are being computed, then
366 *        IFRSTM is the row of the last splitting row above row ILAST;
367 *        this is always at least ILO.
368 *     IITER counts iterations since the last eigenvalue was found,
369 *        to tell when to use an extraordinary shift.
370 *     MAXIT is the maximum number of QZ sweeps allowed.
371 *
372       ILAST = IHI
373       IF( ILSCHR ) THEN
374          IFRSTM = 1
375          ILASTM = N
376       ELSE
377          IFRSTM = ILO
378          ILASTM = IHI
379       END IF
380       IITER = 0
381       ESHIFT = CZERO
382       MAXIT = 30*( IHI-ILO+1 )
383 *
384       DO 170 JITER = 1, MAXIT
385 *
386 *        Check for too many iterations.
387 *
388          IF( JITER.GT.MAXIT )
389      $      GO TO 180
390 *
391 *        Split the matrix if possible.
392 *
393 *        Two tests:
394 *           1: H(j,j-1)=0  or  j=ILO
395 *           2: T(j,j)=0
396 *
397 *        Special case: j=ILAST
398 *
399          IF( ILAST.EQ.ILO ) THEN
400             GO TO 60
401          ELSE
402             IF( ABS1( H( ILAST, ILAST-1 ) ).LE.ATOL ) THEN
403                H( ILAST, ILAST-1 ) = CZERO
404                GO TO 60
405             END IF
406          END IF
407 *
408          IFABS( T( ILAST, ILAST ) ).LE.BTOL ) THEN
409             T( ILAST, ILAST ) = CZERO
410             GO TO 50
411          END IF
412 *
413 *        General case: j<ILAST
414 *
415          DO 40 J = ILAST - 1, ILO, -1
416 *
417 *           Test 1: for H(j,j-1)=0 or j=ILO
418 *
419             IF( J.EQ.ILO ) THEN
420                ILAZRO = .TRUE.
421             ELSE
422                IF( ABS1( H( J, J-1 ) ).LE.ATOL ) THEN
423                   H( J, J-1 ) = CZERO
424                   ILAZRO = .TRUE.
425                ELSE
426                   ILAZRO = .FALSE.
427                END IF
428             END IF
429 *
430 *           Test 2: for T(j,j)=0
431 *
432             IFABS( T( J, J ) ).LT.BTOL ) THEN
433                T( J, J ) = CZERO
434 *
435 *              Test 1a: Check for 2 consecutive small subdiagonals in A
436 *
437                ILAZR2 = .FALSE.
438                IF.NOT.ILAZRO ) THEN
439                   IF( ABS1( H( J, J-1 ) )*( ASCALE*ABS1( H( J+1,
440      $                J ) ) ).LE.ABS1( H( J, J ) )*( ASCALE*ATOL ) )
441      $                ILAZR2 = .TRUE.
442                END IF
443 *
444 *              If both tests pass (1 & 2), i.e., the leading diagonal
445 *              element of B in the block is zero, split a 1x1 block off
446 *              at the top. (I.e., at the J-th row/column) The leading
447 *              diagonal element of the remainder can also be zero, so
448 *              this may have to be done repeatedly.
449 *
450                IF( ILAZRO .OR. ILAZR2 ) THEN
451                   DO 20 JCH = J, ILAST - 1
452                      CTEMP = H( JCH, JCH )
453                      CALL ZLARTG( CTEMP, H( JCH+1, JCH ), C, S,
454      $                            H( JCH, JCH ) )
455                      H( JCH+1, JCH ) = CZERO
456                      CALL ZROT( ILASTM-JCH, H( JCH, JCH+1 ), LDH,
457      $                          H( JCH+1, JCH+1 ), LDH, C, S )
458                      CALL ZROT( ILASTM-JCH, T( JCH, JCH+1 ), LDT,
459      $                          T( JCH+1, JCH+1 ), LDT, C, S )
460                      IF( ILQ )
461      $                  CALL ZROT( N, Q( 1, JCH ), 1, Q( 1, JCH+1 ), 1,
462      $                             C, DCONJG( S ) )
463                      IF( ILAZR2 )
464      $                  H( JCH, JCH-1 ) = H( JCH, JCH-1 )*C
465                      ILAZR2 = .FALSE.
466                      IF( ABS1( T( JCH+1, JCH+1 ) ).GE.BTOL ) THEN
467                         IF( JCH+1.GE.ILAST ) THEN
468                            GO TO 60
469                         ELSE
470                            IFIRST = JCH + 1
471                            GO TO 70
472                         END IF
473                      END IF
474                      T( JCH+1, JCH+1 ) = CZERO
475    20             CONTINUE
476                   GO TO 50
477                ELSE
478 *
479 *                 Only test 2 passed -- chase the zero to T(ILAST,ILAST)
480 *                 Then process as in the case T(ILAST,ILAST)=0
481 *
482                   DO 30 JCH = J, ILAST - 1
483                      CTEMP = T( JCH, JCH+1 )
484                      CALL ZLARTG( CTEMP, T( JCH+1, JCH+1 ), C, S,
485      $                            T( JCH, JCH+1 ) )
486                      T( JCH+1, JCH+1 ) = CZERO
487                      IF( JCH.LT.ILASTM-1 )
488      $                  CALL ZROT( ILASTM-JCH-1, T( JCH, JCH+2 ), LDT,
489      $                             T( JCH+1, JCH+2 ), LDT, C, S )
490                      CALL ZROT( ILASTM-JCH+2, H( JCH, JCH-1 ), LDH,
491      $                          H( JCH+1, JCH-1 ), LDH, C, S )
492                      IF( ILQ )
493      $                  CALL ZROT( N, Q( 1, JCH ), 1, Q( 1, JCH+1 ), 1,
494      $                             C, DCONJG( S ) )
495                      CTEMP = H( JCH+1, JCH )
496                      CALL ZLARTG( CTEMP, H( JCH+1, JCH-1 ), C, S,
497      $                            H( JCH+1, JCH ) )
498                      H( JCH+1, JCH-1 ) = CZERO
499                      CALL ZROT( JCH+1-IFRSTM, H( IFRSTM, JCH ), 1,
500      $                          H( IFRSTM, JCH-1 ), 1, C, S )
501                      CALL ZROT( JCH-IFRSTM, T( IFRSTM, JCH ), 1,
502      $                          T( IFRSTM, JCH-1 ), 1, C, S )
503                      IF( ILZ )
504      $                  CALL ZROT( N, Z( 1, JCH ), 1, Z( 1, JCH-1 ), 1,
505      $                             C, S )
506    30             CONTINUE
507                   GO TO 50
508                END IF
509             ELSE IF( ILAZRO ) THEN
510 *
511 *              Only test 1 passed -- work on J:ILAST
512 *
513                IFIRST = J
514                GO TO 70
515             END IF
516 *
517 *           Neither test passed -- try next J
518 *
519    40    CONTINUE
520 *
521 *        (Drop-through is "impossible")
522 *
523          INFO = 2*+ 1
524          GO TO 210
525 *
526 *        T(ILAST,ILAST)=0 -- clear H(ILAST,ILAST-1) to split off a
527 *        1x1 block.
528 *
529    50    CONTINUE
530          CTEMP = H( ILAST, ILAST )
531          CALL ZLARTG( CTEMP, H( ILAST, ILAST-1 ), C, S,
532      $                H( ILAST, ILAST ) )
533          H( ILAST, ILAST-1 ) = CZERO
534          CALL ZROT( ILAST-IFRSTM, H( IFRSTM, ILAST ), 1,
535      $              H( IFRSTM, ILAST-1 ), 1, C, S )
536          CALL ZROT( ILAST-IFRSTM, T( IFRSTM, ILAST ), 1,
537      $              T( IFRSTM, ILAST-1 ), 1, C, S )
538          IF( ILZ )
539      $      CALL ZROT( N, Z( 1, ILAST ), 1, Z( 1, ILAST-1 ), 1, C, S )
540 *
541 *        H(ILAST,ILAST-1)=0 -- Standardize B, set ALPHA and BETA
542 *
543    60    CONTINUE
544          ABSB = ABS( T( ILAST, ILAST ) )
545          IF( ABSB.GT.SAFMIN ) THEN
546             SIGNBC = DCONJG( T( ILAST, ILAST ) / ABSB )
547             T( ILAST, ILAST ) = ABSB
548             IF( ILSCHR ) THEN
549                CALL ZSCAL( ILAST-IFRSTM, SIGNBC, T( IFRSTM, ILAST ), 1 )
550                CALL ZSCAL( ILAST+1-IFRSTM, SIGNBC, H( IFRSTM, ILAST ),
551      $                     1 )
552             ELSE
553                H( ILAST, ILAST ) = H( ILAST, ILAST )*SIGNBC
554             END IF
555             IF( ILZ )
556      $         CALL ZSCAL( N, SIGNBC, Z( 1, ILAST ), 1 )
557          ELSE
558             T( ILAST, ILAST ) = CZERO
559          END IF
560          ALPHA( ILAST ) = H( ILAST, ILAST )
561          BETA( ILAST ) = T( ILAST, ILAST )
562 *
563 *        Go to next block -- exit if finished.
564 *
565          ILAST = ILAST - 1
566          IF( ILAST.LT.ILO )
567      $      GO TO 190
568 *
569 *        Reset counters
570 *
571          IITER = 0
572          ESHIFT = CZERO
573          IF.NOT.ILSCHR ) THEN
574             ILASTM = ILAST
575             IF( IFRSTM.GT.ILAST )
576      $         IFRSTM = ILO
577          END IF
578          GO TO 160
579 *
580 *        QZ step
581 *
582 *        This iteration only involves rows/columns IFIRST:ILAST.  We
583 *        assume IFIRST < ILAST, and that the diagonal of B is non-zero.
584 *
585    70    CONTINUE
586          IITER = IITER + 1
587          IF.NOT.ILSCHR ) THEN
588             IFRSTM = IFIRST
589          END IF
590 *
591 *        Compute the Shift.
592 *
593 *        At this point, IFIRST < ILAST, and the diagonal elements of
594 *        T(IFIRST:ILAST,IFIRST,ILAST) are larger than BTOL (in
595 *        magnitude)
596 *
597          IF( ( IITER / 10 )*10.NE.IITER ) THEN
598 *
599 *           The Wilkinson shift (AEP p.512), i.e., the eigenvalue of
600 *           the bottom-right 2x2 block of A inv(B) which is nearest to
601 *           the bottom-right element.
602 *
603 *           We factor B as U*D, where U has unit diagonals, and
604 *           compute (A*inv(D))*inv(U).
605 *
606             U12 = ( BSCALE*T( ILAST-1, ILAST ) ) /
607      $            ( BSCALE*T( ILAST, ILAST ) )
608             AD11 = ( ASCALE*H( ILAST-1, ILAST-1 ) ) /
609      $             ( BSCALE*T( ILAST-1, ILAST-1 ) )
610             AD21 = ( ASCALE*H( ILAST, ILAST-1 ) ) /
611      $             ( BSCALE*T( ILAST-1, ILAST-1 ) )
612             AD12 = ( ASCALE*H( ILAST-1, ILAST ) ) /
613      $             ( BSCALE*T( ILAST, ILAST ) )
614             AD22 = ( ASCALE*H( ILAST, ILAST ) ) /
615      $             ( BSCALE*T( ILAST, ILAST ) )
616             ABI22 = AD22 - U12*AD21
617 *
618             T1 = HALF*( AD11+ABI22 )
619             RTDISC = SQRT( T1**2+AD12*AD21-AD11*AD22 )
620             TEMP = DBLE( T1-ABI22 )*DBLE( RTDISC ) +
621      $             DIMAG( T1-ABI22 )*DIMAG( RTDISC )
622             IF( TEMP.LE.ZERO ) THEN
623                SHIFT = T1 + RTDISC
624             ELSE
625                SHIFT = T1 - RTDISC
626             END IF
627          ELSE
628 *
629 *           Exceptional shift.  Chosen for no particularly good reason.
630 *
631             ESHIFT = ESHIFT + DCONJG( ( ASCALE*H( ILAST-1, ILAST ) ) /
632      $               ( BSCALE*T( ILAST-1, ILAST-1 ) ) )
633             SHIFT = ESHIFT
634          END IF
635 *
636 *        Now check for two consecutive small subdiagonals.
637 *
638          DO 80 J = ILAST - 1, IFIRST + 1-1
639             ISTART = J
640             CTEMP = ASCALE*H( J, J ) - SHIFT*( BSCALE*T( J, J ) )
641             TEMP = ABS1( CTEMP )
642             TEMP2 = ASCALE*ABS1( H( J+1, J ) )
643             TEMPR = MAX( TEMP, TEMP2 )
644             IF( TEMPR.LT.ONE .AND. TEMPR.NE.ZERO ) THEN
645                TEMP = TEMP / TEMPR
646                TEMP2 = TEMP2 / TEMPR
647             END IF
648             IF( ABS1( H( J, J-1 ) )*TEMP2.LE.TEMP*ATOL )
649      $         GO TO 90
650    80    CONTINUE
651 *
652          ISTART = IFIRST
653          CTEMP = ASCALE*H( IFIRST, IFIRST ) -
654      $           SHIFT*( BSCALE*T( IFIRST, IFIRST ) )
655    90    CONTINUE
656 *
657 *        Do an implicit-shift QZ sweep.
658 *
659 *        Initial Q
660 *
661          CTEMP2 = ASCALE*H( ISTART+1, ISTART )
662          CALL ZLARTG( CTEMP, CTEMP2, C, S, CTEMP3 )
663 *
664 *        Sweep
665 *
666          DO 150 J = ISTART, ILAST - 1
667             IF( J.GT.ISTART ) THEN
668                CTEMP = H( J, J-1 )
669                CALL ZLARTG( CTEMP, H( J+1, J-1 ), C, S, H( J, J-1 ) )
670                H( J+1, J-1 ) = CZERO
671             END IF
672 *
673             DO 100 JC = J, ILASTM
674                CTEMP = C*H( J, JC ) + S*H( J+1, JC )
675                H( J+1, JC ) = -DCONJG( S )*H( J, JC ) + C*H( J+1, JC )
676                H( J, JC ) = CTEMP
677                CTEMP2 = C*T( J, JC ) + S*T( J+1, JC )
678                T( J+1, JC ) = -DCONJG( S )*T( J, JC ) + C*T( J+1, JC )
679                T( J, JC ) = CTEMP2
680   100       CONTINUE
681             IF( ILQ ) THEN
682                DO 110 JR = 1, N
683                   CTEMP = C*Q( JR, J ) + DCONJG( S )*Q( JR, J+1 )
684                   Q( JR, J+1 ) = -S*Q( JR, J ) + C*Q( JR, J+1 )
685                   Q( JR, J ) = CTEMP
686   110          CONTINUE
687             END IF
688 *
689             CTEMP = T( J+1, J+1 )
690             CALL ZLARTG( CTEMP, T( J+1, J ), C, S, T( J+1, J+1 ) )
691             T( J+1, J ) = CZERO
692 *
693             DO 120 JR = IFRSTM, MIN( J+2, ILAST )
694                CTEMP = C*H( JR, J+1 ) + S*H( JR, J )
695                H( JR, J ) = -DCONJG( S )*H( JR, J+1 ) + C*H( JR, J )
696                H( JR, J+1 ) = CTEMP
697   120       CONTINUE
698             DO 130 JR = IFRSTM, J
699                CTEMP = C*T( JR, J+1 ) + S*T( JR, J )
700                T( JR, J ) = -DCONJG( S )*T( JR, J+1 ) + C*T( JR, J )
701                T( JR, J+1 ) = CTEMP
702   130       CONTINUE
703             IF( ILZ ) THEN
704                DO 140 JR = 1, N
705                   CTEMP = C*Z( JR, J+1 ) + S*Z( JR, J )
706                   Z( JR, J ) = -DCONJG( S )*Z( JR, J+1 ) + C*Z( JR, J )
707                   Z( JR, J+1 ) = CTEMP
708   140          CONTINUE
709             END IF
710   150    CONTINUE
711 *
712   160    CONTINUE
713 *
714   170 CONTINUE
715 *
716 *     Drop-through = non-convergence
717 *
718   180 CONTINUE
719       INFO = ILAST
720       GO TO 210
721 *
722 *     Successful completion of all QZ steps
723 *
724   190 CONTINUE
725 *
726 *     Set Eigenvalues 1:ILO-1
727 *
728       DO 200 J = 1, ILO - 1
729          ABSB = ABS( T( J, J ) )
730          IF( ABSB.GT.SAFMIN ) THEN
731             SIGNBC = DCONJG( T( J, J ) / ABSB )
732             T( J, J ) = ABSB
733             IF( ILSCHR ) THEN
734                CALL ZSCAL( J-1, SIGNBC, T( 1, J ), 1 )
735                CALL ZSCAL( J, SIGNBC, H( 1, J ), 1 )
736             ELSE
737                H( J, J ) = H( J, J )*SIGNBC
738             END IF
739             IF( ILZ )
740      $         CALL ZSCAL( N, SIGNBC, Z( 1, J ), 1 )
741          ELSE
742             T( J, J ) = CZERO
743          END IF
744          ALPHA( J ) = H( J, J )
745          BETA( J ) = T( J, J )
746   200 CONTINUE
747 *
748 *     Normal Termination
749 *
750       INFO = 0
751 *
752 *     Exit (other than argument error) -- return optimal workspace size
753 *
754   210 CONTINUE
755       WORK( 1 ) = DCMPLX( N )
756       RETURN
757 *
758 *     End of ZHGEQZ
759 *
760       END