SUBROUTINE SGEGS( JOBVSL, JOBVSR, N, A, LDA, B, LDB, ALPHAR, $ ALPHAI, BETA, VSL, LDVSL, VSR, LDVSR, WORK, $ LWORK, INFO ) * * -- LAPACK driver routine (version 2.0) -- * Univ. of Tennessee, Univ. of California Berkeley, NAG Ltd., * Courant Institute, Argonne National Lab, and Rice University * September 30, 1994 * * .. Scalar Arguments .. CHARACTER JOBVSL, JOBVSR INTEGER INFO, LDA, LDB, LDVSL, LDVSR, LWORK, N * .. * .. Array Arguments .. REAL A( LDA, * ), ALPHAI( * ), ALPHAR( * ), $ B( LDB, * ), BETA( * ), VSL( LDVSL, * ), $ VSR( LDVSR, * ), WORK( * ) * .. * * Purpose * ======= * * SGEGS computes for a pair of N-by-N real nonsymmetric matrices A, B: * the generalized eigenvalues (alphar +/- alphai*i, beta), the real * Schur form (A, B), and optionally left and/or right Schur vectors * (VSL and VSR). * * (If only the generalized eigenvalues are needed, use the driver SGEGV * instead.) * * A generalized eigenvalue for a pair of matrices (A,B) is, roughly * speaking, a scalar w or a ratio alpha/beta = w, such that A - w*B * is singular. It is usually represented as the pair (alpha,beta), * as there is a reasonable interpretation for beta=0, and even for * both being zero. A good beginning reference is the book, "Matrix * Computations", by G. Golub & C. van Loan (Johns Hopkins U. Press) * * The (generalized) Schur form of a pair of matrices is the result of * multiplying both matrices on the left by one orthogonal matrix and * both on the right by another orthogonal matrix, these two orthogonal * matrices being chosen so as to bring the pair of matrices into * (real) Schur form. * * A pair of matrices A, B is in generalized real Schur form if B is * upper triangular with non-negative diagonal and A is block upper * triangular with 1-by-1 and 2-by-2 blocks. 1-by-1 blocks correspond * to real generalized eigenvalues, while 2-by-2 blocks of A will be * "standardized" by making the corresponding elements of B have the * form: * [ a 0 ] * [ 0 b ] * * and the pair of corresponding 2-by-2 blocks in A and B will * have a complex conjugate pair of generalized eigenvalues. * * The left and right Schur vectors are the columns of VSL and VSR, * respectively, where VSL and VSR are the orthogonal matrices * which reduce A and B to Schur form: * * Schur form of (A,B) = ( (VSL)**T A (VSR), (VSL)**T B (VSR) ) * * Arguments * ========= * * JOBVSL (input) CHARACTER*1 * = 'N': do not compute the left Schur vectors; * = 'V': compute the left Schur vectors. * * JOBVSR (input) CHARACTER*1 * = 'N': do not compute the right Schur vectors; * = 'V': compute the right Schur vectors. * * N (input) INTEGER * The order of the matrices A, B, VSL, and VSR. N >= 0. * * A (input/output) REAL array, dimension (LDA, N) * On entry, the first of the pair of matrices whose generalized * eigenvalues and (optionally) Schur vectors are to be * computed. * On exit, the generalized Schur form of A. * Note: to avoid overflow, the Frobenius norm of the matrix * A should be less than the overflow threshold. * * LDA (input) INTEGER * The leading dimension of A. LDA >= max(1,N). * * B (input/output) REAL array, dimension (LDB, N) * On entry, the second of the pair of matrices whose * generalized eigenvalues and (optionally) Schur vectors are * to be computed. * On exit, the generalized Schur form of B. * Note: to avoid overflow, the Frobenius norm of the matrix * B should be less than the overflow threshold. * * LDB (input) INTEGER * The leading dimension of B. LDB >= max(1,N). * * ALPHAR (output) REAL array, dimension (N) * ALPHAI (output) REAL array, dimension (N) * BETA (output) REAL array, dimension (N) * On exit, (ALPHAR(j) + ALPHAI(j)*i)/BETA(j), j=1,...,N, will * be the generalized eigenvalues. ALPHAR(j) + ALPHAI(j)*i, * j=1,...,N and BETA(j),j=1,...,N are the diagonals of the * complex Schur form (A,B) that would result if the 2-by-2 * diagonal blocks of the real Schur form of (A,B) were further * reduced to triangular form using 2-by-2 complex unitary * transformations. If ALPHAI(j) is zero, then the j-th * eigenvalue is real; if positive, then the j-th and (j+1)-st * eigenvalues are a complex conjugate pair, with ALPHAI(j+1) * negative. * * Note: the quotients ALPHAR(j)/BETA(j) and ALPHAI(j)/BETA(j) * may easily over- or underflow, and BETA(j) may even be zero. * Thus, the user should avoid naively computing the ratio * alpha/beta. However, ALPHAR and ALPHAI will be always less * than and usually comparable with norm(A) in magnitude, and * BETA always less than and usually comparable with norm(B). * * VSL (output) REAL array, dimension (LDVSL,N) * If JOBVSL = 'V', VSL will contain the left Schur vectors. * (See "Purpose", above.) * Not referenced if JOBVSL = 'N'. * * LDVSL (input) INTEGER * The leading dimension of the matrix VSL. LDVSL >=1, and * if JOBVSL = 'V', LDVSL >= N. * * VSR (output) REAL array, dimension (LDVSR,N) * If JOBVSR = 'V', VSR will contain the right Schur vectors. * (See "Purpose", above.) * Not referenced if JOBVSR = 'N'. * * LDVSR (input) INTEGER * The leading dimension of the matrix VSR. LDVSR >= 1, and * if JOBVSR = 'V', LDVSR >= N. * * WORK (workspace/output) REAL 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,4*N). * For good performance, LWORK must generally be larger. * To compute the optimal value of LWORK, call ILAENV to get * blocksizes (for SGEQRF, SORMQR, and SORGQR.) Then compute: * NB -- MAX of the blocksizes for SGEQRF, SORMQR, and SORGQR * The optimal LWORK is 2*N + N*(NB+1). * * INFO (output) INTEGER * = 0: successful exit * < 0: if INFO = -i, the i-th argument had an illegal value. * = 1,...,N: * The QZ iteration failed. (A,B) are not in Schur * form, but ALPHAR(j), ALPHAI(j), and BETA(j) should * be correct for j=INFO+1,...,N. * > N: errors that usually indicate LAPACK problems: * =N+1: error return from SGGBAL * =N+2: error return from SGEQRF * =N+3: error return from SORMQR * =N+4: error return from SORGQR * =N+5: error return from SGGHRD * =N+6: error return from SHGEQZ (other than failed * iteration) * =N+7: error return from SGGBAK (computing VSL) * =N+8: error return from SGGBAK (computing VSR) * =N+9: error return from SLASCL (various places) * * ===================================================================== * * .. Parameters .. REAL ZERO, ONE PARAMETER ( ZERO = 0.0E0, ONE = 1.0E0 ) * .. * .. Local Scalars .. LOGICAL ILASCL, ILBSCL, ILVSL, ILVSR INTEGER ICOLS, IHI, IINFO, IJOBVL, IJOBVR, ILEFT, $ ILO, IRIGHT, IROWS, ITAU, IWORK, LWKMIN, LWKOPT REAL ANRM, ANRMTO, BIGNUM, BNRM, BNRMTO, EPS, $ SAFMIN, SMLNUM * .. * .. External Subroutines .. EXTERNAL SGEQRF, SGGBAK, SGGBAL, SGGHRD, SHGEQZ, SLACPY, $ SLASCL, SLASET, SORGQR, SORMQR, XERBLA * .. * .. External Functions .. LOGICAL LSAME REAL SLAMCH, SLANGE EXTERNAL LSAME, SLAMCH, SLANGE * .. * .. Intrinsic Functions .. INTRINSIC INT, MAX * .. * .. Executable Statements .. * * Decode the input arguments * IF( LSAME( JOBVSL, 'N' ) ) THEN IJOBVL = 1 ILVSL = .FALSE. ELSE IF( LSAME( JOBVSL, 'V' ) ) THEN IJOBVL = 2 ILVSL = .TRUE. ELSE IJOBVL = -1 ILVSL = .FALSE. END IF * IF( LSAME( JOBVSR, 'N' ) ) THEN IJOBVR = 1 ILVSR = .FALSE. ELSE IF( LSAME( JOBVSR, 'V' ) ) THEN IJOBVR = 2 ILVSR = .TRUE. ELSE IJOBVR = -1 ILVSR = .FALSE. END IF * * Test the input arguments * LWKMIN = MAX( 4*N, 1 ) LWKOPT = LWKMIN INFO = 0 IF( IJOBVL.LE.0 ) THEN INFO = -1 ELSE IF( IJOBVR.LE.0 ) THEN INFO = -2 ELSE IF( N.LT.0 ) THEN INFO = -3 ELSE IF( LDA.LT.MAX( 1, N ) ) THEN INFO = -5 ELSE IF( LDB.LT.MAX( 1, N ) ) THEN INFO = -7 ELSE IF( LDVSL.LT.1 .OR. ( ILVSL .AND. LDVSL.LT.N ) ) THEN INFO = -12 ELSE IF( LDVSR.LT.1 .OR. ( ILVSR .AND. LDVSR.LT.N ) ) THEN INFO = -14 ELSE IF( LWORK.LT.LWKMIN ) THEN INFO = -16 END IF * IF( INFO.NE.0 ) THEN CALL XERBLA( 'SGEGS ', -INFO ) RETURN END IF * * Quick return if possible * WORK( 1 ) = LWKOPT IF( N.EQ.0 ) $ RETURN * * Get machine constants * EPS = SLAMCH( 'E' )*SLAMCH( 'B' ) SAFMIN = SLAMCH( 'S' ) SMLNUM = N*SAFMIN / EPS BIGNUM = ONE / SMLNUM * * Scale A if max element outside range [SMLNUM,BIGNUM] * ANRM = SLANGE( 'M', N, N, A, LDA, WORK ) ILASCL = .FALSE. IF( ANRM.GT.ZERO .AND. ANRM.LT.SMLNUM ) THEN ANRMTO = SMLNUM ILASCL = .TRUE. ELSE IF( ANRM.GT.BIGNUM ) THEN ANRMTO = BIGNUM ILASCL = .TRUE. END IF * IF( ILASCL ) THEN CALL SLASCL( 'G', -1, -1, ANRM, ANRMTO, N, N, A, LDA, IINFO ) IF( IINFO.NE.0 ) THEN INFO = N + 9 RETURN END IF END IF * * Scale B if max element outside range [SMLNUM,BIGNUM] * BNRM = SLANGE( 'M', N, N, B, LDB, WORK ) ILBSCL = .FALSE. IF( BNRM.GT.ZERO .AND. BNRM.LT.SMLNUM ) THEN BNRMTO = SMLNUM ILBSCL = .TRUE. ELSE IF( BNRM.GT.BIGNUM ) THEN BNRMTO = BIGNUM ILBSCL = .TRUE. END IF * IF( ILBSCL ) THEN CALL SLASCL( 'G', -1, -1, BNRM, BNRMTO, N, N, B, LDB, IINFO ) IF( IINFO.NE.0 ) THEN INFO = N + 9 RETURN END IF END IF * * Permute the matrix to make it more nearly triangular * Workspace layout: (2*N words -- "work..." not actually used) * left_permutation, right_permutation, work... * ILEFT = 1 IRIGHT = N + 1 IWORK = IRIGHT + N CALL SGGBAL( 'P', N, A, LDA, B, LDB, ILO, IHI, WORK( ILEFT ), $ WORK( IRIGHT ), WORK( IWORK ), IINFO ) IF( IINFO.NE.0 ) THEN INFO = N + 1 GO TO 10 END IF * * Reduce B to triangular form, and initialize VSL and/or VSR * Workspace layout: ("work..." must have at least N words) * left_permutation, right_permutation, tau, work... * IROWS = IHI + 1 - ILO ICOLS = N + 1 - ILO ITAU = IWORK IWORK = ITAU + IROWS CALL SGEQRF( IROWS, ICOLS, B( ILO, ILO ), LDB, WORK( ITAU ), $ WORK( IWORK ), LWORK+1-IWORK, IINFO ) IF( IINFO.GE.0 ) $ LWKOPT = MAX( LWKOPT, INT( WORK( IWORK ) )+IWORK-1 ) IF( IINFO.NE.0 ) THEN INFO = N + 2 GO TO 10 END IF * CALL SORMQR( 'L', 'T', IROWS, ICOLS, IROWS, B( ILO, ILO ), LDB, $ WORK( ITAU ), A( ILO, ILO ), LDA, WORK( IWORK ), $ LWORK+1-IWORK, IINFO ) IF( IINFO.GE.0 ) $ LWKOPT = MAX( LWKOPT, INT( WORK( IWORK ) )+IWORK-1 ) IF( IINFO.NE.0 ) THEN INFO = N + 3 GO TO 10 END IF * IF( ILVSL ) THEN CALL SLASET( 'Full', N, N, ZERO, ONE, VSL, LDVSL ) CALL SLACPY( 'L', IROWS-1, IROWS-1, B( ILO+1, ILO ), LDB, $ VSL( ILO+1, ILO ), LDVSL ) CALL SORGQR( IROWS, IROWS, IROWS, VSL( ILO, ILO ), LDVSL, $ WORK( ITAU ), WORK( IWORK ), LWORK+1-IWORK, $ IINFO ) IF( IINFO.GE.0 ) $ LWKOPT = MAX( LWKOPT, INT( WORK( IWORK ) )+IWORK-1 ) IF( IINFO.NE.0 ) THEN INFO = N + 4 GO TO 10 END IF END IF * IF( ILVSR ) $ CALL SLASET( 'Full', N, N, ZERO, ONE, VSR, LDVSR ) * * Reduce to generalized Hessenberg form * CALL SGGHRD( JOBVSL, JOBVSR, N, ILO, IHI, A, LDA, B, LDB, VSL, $ LDVSL, VSR, LDVSR, IINFO ) IF( IINFO.NE.0 ) THEN INFO = N + 5 GO TO 10 END IF * * Perform QZ algorithm, computing Schur vectors if desired * Workspace layout: ("work..." must have at least 1 word) * left_permutation, right_permutation, work... * IWORK = ITAU CALL SHGEQZ( 'S', JOBVSL, JOBVSR, N, ILO, IHI, A, LDA, B, LDB, $ ALPHAR, ALPHAI, BETA, VSL, LDVSL, VSR, LDVSR, $ WORK( IWORK ), LWORK+1-IWORK, IINFO ) IF( IINFO.GE.0 ) $ LWKOPT = MAX( LWKOPT, INT( WORK( IWORK ) )+IWORK-1 ) IF( IINFO.NE.0 ) THEN IF( IINFO.GT.0 .AND. IINFO.LE.N ) THEN INFO = IINFO ELSE IF( IINFO.GT.N .AND. IINFO.LE.2*N ) THEN INFO = IINFO - N ELSE INFO = N + 6 END IF GO TO 10 END IF * * Apply permutation to VSL and VSR * IF( ILVSL ) THEN CALL SGGBAK( 'P', 'L', N, ILO, IHI, WORK( ILEFT ), $ WORK( IRIGHT ), N, VSL, LDVSL, IINFO ) IF( IINFO.NE.0 ) THEN INFO = N + 7 GO TO 10 END IF END IF IF( ILVSR ) THEN CALL SGGBAK( 'P', 'R', N, ILO, IHI, WORK( ILEFT ), $ WORK( IRIGHT ), N, VSR, LDVSR, IINFO ) IF( IINFO.NE.0 ) THEN INFO = N + 8 GO TO 10 END IF END IF * * Undo scaling * IF( ILASCL ) THEN CALL SLASCL( 'H', -1, -1, ANRMTO, ANRM, N, N, A, LDA, IINFO ) IF( IINFO.NE.0 ) THEN INFO = N + 9 RETURN END IF CALL SLASCL( 'G', -1, -1, ANRMTO, ANRM, N, 1, ALPHAR, N, $ IINFO ) IF( IINFO.NE.0 ) THEN INFO = N + 9 RETURN END IF CALL SLASCL( 'G', -1, -1, ANRMTO, ANRM, N, 1, ALPHAI, N, $ IINFO ) IF( IINFO.NE.0 ) THEN INFO = N + 9 RETURN END IF END IF * IF( ILBSCL ) THEN CALL SLASCL( 'U', -1, -1, BNRMTO, BNRM, N, N, B, LDB, IINFO ) IF( IINFO.NE.0 ) THEN INFO = N + 9 RETURN END IF CALL SLASCL( 'G', -1, -1, BNRMTO, BNRM, N, 1, BETA, N, IINFO ) IF( IINFO.NE.0 ) THEN INFO = N + 9 RETURN END IF END IF * 10 CONTINUE WORK( 1 ) = LWKOPT * RETURN * * End of SGEGS * END