SUBROUTINE CHGEQZ( JOB, COMPQ, COMPZ, N, ILO, IHI, A, LDA, B, LDB, $ ALPHA, BETA, Q, LDQ, Z, LDZ, WORK, LWORK, $ RWORK, INFO ) * * -- LAPACK routine (version 3.0) -- * Univ. of Tennessee, Univ. of California Berkeley, NAG Ltd., * Courant Institute, Argonne National Lab, and Rice University * June 30, 1999 * * .. Scalar Arguments .. CHARACTER COMPQ, COMPZ, JOB INTEGER IHI, ILO, INFO, LDA, LDB, LDQ, LDZ, LWORK, N * .. * .. Array Arguments .. REAL RWORK( * ) COMPLEX A( LDA, * ), ALPHA( * ), B( LDB, * ), $ BETA( * ), Q( LDQ, * ), WORK( * ), Z( LDZ, * ) * .. * * Purpose * ======= * * CHGEQZ implements a single-shift version of the QZ * method for finding the generalized eigenvalues w(i)=ALPHA(i)/BETA(i) * of the equation * * det( A - w(i) B ) = 0 * * If JOB='S', then the pair (A,B) is simultaneously * reduced to Schur form (i.e., A and B are both upper triangular) by * applying one unitary tranformation (usually called Q) on the left and * another (usually called Z) on the right. The diagonal elements of * A are then ALPHA(1),...,ALPHA(N), and of B are BETA(1),...,BETA(N). * * If JOB='S' and COMPQ and COMPZ are 'V' or 'I', then the unitary * transformations used to reduce (A,B) are accumulated into the arrays * Q and Z s.t.: * * Q(in) A(in) Z(in)* = Q(out) A(out) Z(out)* * Q(in) B(in) Z(in)* = Q(out) B(out) Z(out)* * * Ref: C.B. Moler & G.W. Stewart, "An Algorithm for Generalized Matrix * Eigenvalue Problems", SIAM J. Numer. Anal., 10(1973), * pp. 241--256. * * Arguments * ========= * * JOB (input) CHARACTER*1 * = 'E': compute only ALPHA and BETA. A and B will not * necessarily be put into generalized Schur form. * = 'S': put A and B into generalized Schur form, as well * as computing ALPHA and BETA. * * COMPQ (input) CHARACTER*1 * = 'N': do not modify Q. * = 'V': multiply the array Q on the right by the conjugate * transpose of the unitary tranformation that is * applied to the left side of A and B to reduce them * to Schur form. * = 'I': like COMPQ='V', except that Q will be initialized to * the identity first. * * COMPZ (input) CHARACTER*1 * = 'N': do not modify Z. * = 'V': multiply the array Z on the right by the unitary * tranformation that is applied to the right side of * A and B to reduce them to Schur form. * = 'I': like COMPZ='V', except that Z will be initialized to * the identity first. * * N (input) INTEGER * The order of the matrices A, B, Q, and Z. N >= 0. * * ILO (input) INTEGER * IHI (input) INTEGER * It is assumed that A is already upper triangular in rows and * columns 1:ILO-1 and IHI+1:N. * 1 <= ILO <= IHI <= N, if N > 0; ILO=1 and IHI=0, if N=0. * * A (input/output) COMPLEX array, dimension (LDA, N) * On entry, the N-by-N upper Hessenberg matrix A. Elements * below the subdiagonal must be zero. * If JOB='S', then on exit A and B will have been * simultaneously reduced to upper triangular form. * If JOB='E', then on exit A will have been destroyed. * * LDA (input) INTEGER * The leading dimension of the array A. LDA >= max( 1, N ). * * B (input/output) COMPLEX array, dimension (LDB, N) * On entry, the N-by-N upper triangular matrix B. Elements * below the diagonal must be zero. * If JOB='S', then on exit A and B will have been * simultaneously reduced to upper triangular form. * If JOB='E', then on exit B will have been destroyed. * * LDB (input) INTEGER * The leading dimension of the array B. LDB >= max( 1, N ). * * ALPHA (output) COMPLEX array, dimension (N) * The diagonal elements of A when the pair (A,B) has been * reduced to Schur form. ALPHA(i)/BETA(i) i=1,...,N * are the generalized eigenvalues. * * BETA (output) COMPLEX array, dimension (N) * The diagonal elements of B when the pair (A,B) has been * reduced to Schur form. ALPHA(i)/BETA(i) i=1,...,N * are the generalized eigenvalues. A and B are normalized * so that BETA(1),...,BETA(N) are non-negative real numbers. * * Q (input/output) COMPLEX array, dimension (LDQ, N) * If COMPQ='N', then Q will not be referenced. * If COMPQ='V' or 'I', then the conjugate transpose of the * unitary transformations which are applied to A and B on * the left will be applied to the array Q on the right. * * LDQ (input) INTEGER * The leading dimension of the array Q. LDQ >= 1. * If COMPQ='V' or 'I', then LDQ >= N. * * Z (input/output) COMPLEX array, dimension (LDZ, N) * If COMPZ='N', then Z will not be referenced. * If COMPZ='V' or 'I', then the unitary transformations which * are applied to A and B on the right will be applied to the * array Z on the right. * * LDZ (input) INTEGER * The leading dimension of the array Z. LDZ >= 1. * If COMPZ='V' or 'I', then LDZ >= N. * * WORK (workspace/output) COMPLEX array, dimension (LWORK) * On exit, if INFO >= 0, WORK(1) returns the optimal LWORK. * * LWORK (input) INTEGER * The dimension of the array WORK. LWORK >= max(1,N). * * If LWORK = -1, then a workspace query is assumed; the routine * only calculates the optimal size of the WORK array, returns * this value as the first entry of the WORK array, and no error * message related to LWORK is issued by XERBLA. * * RWORK (workspace) REAL array, dimension (N) * * INFO (output) INTEGER * = 0: successful exit * < 0: if INFO = -i, the i-th argument had an illegal value * = 1,...,N: the QZ iteration did not converge. (A,B) is not * in Schur form, but ALPHA(i) and BETA(i), * i=INFO+1,...,N should be correct. * = N+1,...,2*N: the shift calculation failed. (A,B) is not * in Schur form, but ALPHA(i) and BETA(i), * i=INFO-N+1,...,N should be correct. * > 2*N: various "impossible" errors. * * Further Details * =============== * * We assume that complex ABS works as long as its value is less than * overflow. * * ===================================================================== * * .. Parameters .. COMPLEX CZERO, CONE PARAMETER ( CZERO = ( 0.0E+0, 0.0E+0 ), $ CONE = ( 1.0E+0, 0.0E+0 ) ) REAL ZERO, ONE PARAMETER ( ZERO = 0.0E+0, ONE = 1.0E+0 ) REAL HALF PARAMETER ( HALF = 0.5E+0 ) * .. * .. Local Scalars .. LOGICAL ILAZR2, ILAZRO, ILQ, ILSCHR, ILZ, LQUERY INTEGER ICOMPQ, ICOMPZ, IFIRST, IFRSTM, IITER, ILAST, $ ILASTM, IN, ISCHUR, ISTART, J, JC, JCH, JITER, $ JR, MAXIT REAL ABSB, ANORM, ASCALE, ATOL, BNORM, BSCALE, BTOL, $ C, SAFMIN, TEMP, TEMP2, TEMPR, ULP COMPLEX ABI22, AD11, AD12, AD21, AD22, CTEMP, CTEMP2, $ CTEMP3, ESHIFT, RTDISC, S, SHIFT, SIGNBC, T, $ U12, X * .. * .. External Functions .. LOGICAL LSAME REAL CLANHS, SLAMCH EXTERNAL LSAME, CLANHS, SLAMCH * .. * .. External Subroutines .. EXTERNAL CLARTG, CLASET, CROT, CSCAL, XERBLA * .. * .. Intrinsic Functions .. INTRINSIC ABS, AIMAG, CMPLX, CONJG, MAX, MIN, REAL, SQRT * .. * .. Statement Functions .. REAL ABS1 * .. * .. Statement Function definitions .. ABS1( X ) = ABS( REAL( X ) ) + ABS( AIMAG( X ) ) * .. * .. Executable Statements .. * * Decode JOB, COMPQ, COMPZ * IF( LSAME( JOB, 'E' ) ) THEN ILSCHR = .FALSE. ISCHUR = 1 ELSE IF( LSAME( JOB, 'S' ) ) THEN ILSCHR = .TRUE. ISCHUR = 2 ELSE ISCHUR = 0 END IF * IF( LSAME( COMPQ, 'N' ) ) THEN ILQ = .FALSE. ICOMPQ = 1 ELSE IF( LSAME( COMPQ, 'V' ) ) THEN ILQ = .TRUE. ICOMPQ = 2 ELSE IF( LSAME( COMPQ, 'I' ) ) THEN ILQ = .TRUE. ICOMPQ = 3 ELSE ICOMPQ = 0 END IF * IF( LSAME( COMPZ, 'N' ) ) THEN ILZ = .FALSE. ICOMPZ = 1 ELSE IF( LSAME( COMPZ, 'V' ) ) THEN ILZ = .TRUE. ICOMPZ = 2 ELSE IF( LSAME( COMPZ, 'I' ) ) THEN ILZ = .TRUE. ICOMPZ = 3 ELSE ICOMPZ = 0 END IF * * Check Argument Values * INFO = 0 WORK( 1 ) = MAX( 1, N ) LQUERY = ( LWORK.EQ.-1 ) IF( ISCHUR.EQ.0 ) THEN INFO = -1 ELSE IF( ICOMPQ.EQ.0 ) THEN INFO = -2 ELSE IF( ICOMPZ.EQ.0 ) THEN INFO = -3 ELSE IF( N.LT.0 ) THEN INFO = -4 ELSE IF( ILO.LT.1 ) THEN INFO = -5 ELSE IF( IHI.GT.N .OR. IHI.LT.ILO-1 ) THEN INFO = -6 ELSE IF( LDA.LT.N ) THEN INFO = -8 ELSE IF( LDB.LT.N ) THEN INFO = -10 ELSE IF( LDQ.LT.1 .OR. ( ILQ .AND. LDQ.LT.N ) ) THEN INFO = -14 ELSE IF( LDZ.LT.1 .OR. ( ILZ .AND. LDZ.LT.N ) ) THEN INFO = -16 ELSE IF( LWORK.LT.MAX( 1, N ) .AND. .NOT.LQUERY ) THEN INFO = -18 END IF IF( INFO.NE.0 ) THEN CALL XERBLA( 'CHGEQZ', -INFO ) RETURN ELSE IF( LQUERY ) THEN RETURN END IF * * Quick return if possible * c WORK( 1 ) = CMPLX( 1 ) IF( N.LE.0 ) THEN WORK( 1 ) = CMPLX( 1 ) RETURN END IF * * Initialize Q and Z * IF( ICOMPQ.EQ.3 ) $ CALL CLASET( 'Full', N, N, CZERO, CONE, Q, LDQ ) IF( ICOMPZ.EQ.3 ) $ CALL CLASET( 'Full', N, N, CZERO, CONE, Z, LDZ ) * * Machine Constants * IN = IHI + 1 - ILO SAFMIN = SLAMCH( 'S' ) ULP = SLAMCH( 'E' )*SLAMCH( 'B' ) ANORM = CLANHS( 'F', IN, A( ILO, ILO ), LDA, RWORK ) BNORM = CLANHS( 'F', IN, B( ILO, ILO ), LDB, RWORK ) ATOL = MAX( SAFMIN, ULP*ANORM ) BTOL = MAX( SAFMIN, ULP*BNORM ) ASCALE = ONE / MAX( SAFMIN, ANORM ) BSCALE = ONE / MAX( SAFMIN, BNORM ) * * * Set Eigenvalues IHI+1:N * DO 10 J = IHI + 1, N ABSB = ABS( B( J, J ) ) IF( ABSB.GT.SAFMIN ) THEN SIGNBC = CONJG( B( J, J ) / ABSB ) B( J, J ) = ABSB IF( ILSCHR ) THEN CALL CSCAL( J-1, SIGNBC, B( 1, J ), 1 ) CALL CSCAL( J, SIGNBC, A( 1, J ), 1 ) ELSE A( J, J ) = A( J, J )*SIGNBC END IF IF( ILZ ) $ CALL CSCAL( N, SIGNBC, Z( 1, J ), 1 ) ELSE B( J, J ) = CZERO END IF ALPHA( J ) = A( J, J ) BETA( J ) = B( J, J ) 10 CONTINUE * * If IHI < ILO, skip QZ steps * IF( IHI.LT.ILO ) $ GO TO 190 * * MAIN QZ ITERATION LOOP * * Initialize dynamic indices * * Eigenvalues ILAST+1:N have been found. * Column operations modify rows IFRSTM:whatever * Row operations modify columns whatever:ILASTM * * If only eigenvalues are being computed, then * IFRSTM is the row of the last splitting row above row ILAST; * this is always at least ILO. * IITER counts iterations since the last eigenvalue was found, * to tell when to use an extraordinary shift. * MAXIT is the maximum number of QZ sweeps allowed. * ILAST = IHI IF( ILSCHR ) THEN IFRSTM = 1 ILASTM = N ELSE IFRSTM = ILO ILASTM = IHI END IF IITER = 0 ESHIFT = CZERO MAXIT = 30*( IHI-ILO+1 ) * DO 170 JITER = 1, MAXIT * * Check for too many iterations. * IF( JITER.GT.MAXIT ) $ GO TO 180 * * Split the matrix if possible. * * Two tests: * 1: A(j,j-1)=0 or j=ILO * 2: B(j,j)=0 * * Special case: j=ILAST * IF( ILAST.EQ.ILO ) THEN GO TO 60 ELSE IF( ABS1( A( ILAST, ILAST-1 ) ).LE.ATOL ) THEN A( ILAST, ILAST-1 ) = CZERO GO TO 60 END IF END IF * IF( ABS( B( ILAST, ILAST ) ).LE.BTOL ) THEN B( ILAST, ILAST ) = CZERO GO TO 50 END IF * * General case: j