3 * =========== DOCUMENTATION ===========
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
9 *> Download STGEVC + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/stgevc.f">
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/stgevc.f">
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/stgevc.f">
21 * SUBROUTINE STGEVC( SIDE, HOWMNY, SELECT, N, S, LDS, P, LDP, VL,
22 * LDVL, VR, LDVR, MM, M, WORK, INFO )
24 * .. Scalar Arguments ..
25 * CHARACTER HOWMNY, SIDE
26 * INTEGER INFO, LDP, LDS, LDVL, LDVR, M, MM, N
28 * .. Array Arguments ..
30 * REAL P( LDP, * ), S( LDS, * ), VL( LDVL, * ),
31 * $ VR( LDVR, * ), WORK( * )
41 *> STGEVC computes some or all of the right and/or left eigenvectors of
42 *> a pair of real matrices (S,P), where S is a quasi-triangular matrix
43 *> and P is upper triangular. Matrix pairs of this type are produced by
44 *> the generalized Schur factorization of a matrix pair (A,B):
46 *> A = Q*S*Z**T, B = Q*P*Z**T
48 *> as computed by SGGHRD + SHGEQZ.
50 *> The right eigenvector x and the left eigenvector y of (S,P)
51 *> corresponding to an eigenvalue w are defined by:
53 *> S*x = w*P*x, (y**H)*S = w*(y**H)*P,
55 *> where y**H denotes the conjugate tranpose of y.
56 *> The eigenvalues are not input to this routine, but are computed
57 *> directly from the diagonal blocks of S and P.
59 *> This routine returns the matrices X and/or Y of right and left
60 *> eigenvectors of (S,P), or the products Z*X and/or Q*Y,
61 *> where Z and Q are input matrices.
62 *> If Q and Z are the orthogonal factors from the generalized Schur
63 *> factorization of a matrix pair (A,B), then Z*X and Q*Y
64 *> are the matrices of right and left eigenvectors of (A,B).
73 *> SIDE is CHARACTER*1
74 *> = 'R': compute right eigenvectors only;
75 *> = 'L': compute left eigenvectors only;
76 *> = 'B': compute both right and left eigenvectors.
81 *> HOWMNY is CHARACTER*1
82 *> = 'A': compute all right and/or left eigenvectors;
83 *> = 'B': compute all right and/or left eigenvectors,
84 *> backtransformed by the matrices in VR and/or VL;
85 *> = 'S': compute selected right and/or left eigenvectors,
86 *> specified by the logical array SELECT.
91 *> SELECT is LOGICAL array, dimension (N)
92 *> If HOWMNY='S', SELECT specifies the eigenvectors to be
93 *> computed. If w(j) is a real eigenvalue, the corresponding
94 *> real eigenvector is computed if SELECT(j) is .TRUE..
95 *> If w(j) and w(j+1) are the real and imaginary parts of a
96 *> complex eigenvalue, the corresponding complex eigenvector
97 *> is computed if either SELECT(j) or SELECT(j+1) is .TRUE.,
98 *> and on exit SELECT(j) is set to .TRUE. and SELECT(j+1) is
100 *> Not referenced if HOWMNY = 'A' or 'B'.
106 *> The order of the matrices S and P. N >= 0.
111 *> S is REAL array, dimension (LDS,N)
112 *> The upper quasi-triangular matrix S from a generalized Schur
113 *> factorization, as computed by SHGEQZ.
119 *> The leading dimension of array S. LDS >= max(1,N).
124 *> P is REAL array, dimension (LDP,N)
125 *> The upper triangular matrix P from a generalized Schur
126 *> factorization, as computed by SHGEQZ.
127 *> 2-by-2 diagonal blocks of P corresponding to 2-by-2 blocks
128 *> of S must be in positive diagonal form.
134 *> The leading dimension of array P. LDP >= max(1,N).
139 *> VL is REAL array, dimension (LDVL,MM)
140 *> On entry, if SIDE = 'L' or 'B' and HOWMNY = 'B', VL must
141 *> contain an N-by-N matrix Q (usually the orthogonal matrix Q
142 *> of left Schur vectors returned by SHGEQZ).
143 *> On exit, if SIDE = 'L' or 'B', VL contains:
144 *> if HOWMNY = 'A', the matrix Y of left eigenvectors of (S,P);
145 *> if HOWMNY = 'B', the matrix Q*Y;
146 *> if HOWMNY = 'S', the left eigenvectors of (S,P) specified by
147 *> SELECT, stored consecutively in the columns of
148 *> VL, in the same order as their eigenvalues.
150 *> A complex eigenvector corresponding to a complex eigenvalue
151 *> is stored in two consecutive columns, the first holding the
152 *> real part, and the second the imaginary part.
154 *> Not referenced if SIDE = 'R'.
160 *> The leading dimension of array VL. LDVL >= 1, and if
161 *> SIDE = 'L' or 'B', LDVL >= N.
166 *> VR is REAL array, dimension (LDVR,MM)
167 *> On entry, if SIDE = 'R' or 'B' and HOWMNY = 'B', VR must
168 *> contain an N-by-N matrix Z (usually the orthogonal matrix Z
169 *> of right Schur vectors returned by SHGEQZ).
171 *> On exit, if SIDE = 'R' or 'B', VR contains:
172 *> if HOWMNY = 'A', the matrix X of right eigenvectors of (S,P);
173 *> if HOWMNY = 'B' or 'b', the matrix Z*X;
174 *> if HOWMNY = 'S' or 's', the right eigenvectors of (S,P)
175 *> specified by SELECT, stored consecutively in the
176 *> columns of VR, in the same order as their
179 *> A complex eigenvector corresponding to a complex eigenvalue
180 *> is stored in two consecutive columns, the first holding the
181 *> real part and the second the imaginary part.
183 *> Not referenced if SIDE = 'L'.
189 *> The leading dimension of the array VR. LDVR >= 1, and if
190 *> SIDE = 'R' or 'B', LDVR >= N.
196 *> The number of columns in the arrays VL and/or VR. MM >= M.
202 *> The number of columns in the arrays VL and/or VR actually
203 *> used to store the eigenvectors. If HOWMNY = 'A' or 'B', M
204 *> is set to N. Each selected real eigenvector occupies one
205 *> column and each selected complex eigenvector occupies two
211 *> WORK is REAL array, dimension (6*N)
217 *> = 0: successful exit.
218 *> < 0: if INFO = -i, the i-th argument had an illegal value.
219 *> > 0: the 2-by-2 block (INFO:INFO+1) does not have a complex
226 *> \author Univ. of Tennessee
227 *> \author Univ. of California Berkeley
228 *> \author Univ. of Colorado Denver
231 *> \date November 2011
233 *> \ingroup realGEcomputational
235 *> \par Further Details:
236 * =====================
240 *> Allocation of workspace:
241 *> ---------- -- ---------
243 *> WORK( j ) = 1-norm of j-th column of A, above the diagonal
244 *> WORK( N+j ) = 1-norm of j-th column of B, above the diagonal
245 *> WORK( 2*N+1:3*N ) = real part of eigenvector
246 *> WORK( 3*N+1:4*N ) = imaginary part of eigenvector
247 *> WORK( 4*N+1:5*N ) = real part of back-transformed eigenvector
248 *> WORK( 5*N+1:6*N ) = imaginary part of back-transformed eigenvector
250 *> Rowwise vs. columnwise solution methods:
251 *> ------- -- ---------- -------- -------
253 *> Finding a generalized eigenvector consists basically of solving the
254 *> singular triangular system
256 *> (A - w B) x = 0 (for right) or: (A - w B)**H y = 0 (for left)
258 *> Consider finding the i-th right eigenvector (assume all eigenvalues
259 *> are real). The equation to be solved is:
261 *> 0 = sum C(j,k) v(k) = sum C(j,k) v(k) for j = i,. . .,1
264 *> where C = (A - w B) (The components v(i+1:n) are 0.)
266 *> The "rowwise" method is:
269 *> for j = i-1,. . .,1:
271 *> (2) compute s = - sum C(j,k) v(k) and
274 *> (3) v(j) := s / C(j,j)
276 *> Step 2 is sometimes called the "dot product" step, since it is an
277 *> inner product between the j-th row and the portion of the eigenvector
278 *> that has been computed so far.
280 *> The "columnwise" method consists basically in doing the sums
281 *> for all the rows in parallel. As each v(j) is computed, the
282 *> contribution of v(j) times the j-th column of C is added to the
283 *> partial sums. Since FORTRAN arrays are stored columnwise, this has
284 *> the advantage that at each step, the elements of C that are accessed
285 *> are adjacent to one another, whereas with the rowwise method, the
286 *> elements accessed at a step are spaced LDS (and LDP) words apart.
288 *> When finding left eigenvectors, the matrix in question is the
289 *> transpose of the one in storage, so the rowwise method then
290 *> actually accesses columns of A and B at each step, and so is the
294 * =====================================================================
295 SUBROUTINE STGEVC( SIDE, HOWMNY, SELECT, N, S, LDS, P, LDP, VL,
296 $ LDVL, VR, LDVR, MM, M, WORK, INFO )
298 * -- LAPACK computational routine (version 3.4.0) --
299 * -- LAPACK is a software package provided by Univ. of Tennessee, --
300 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
303 * .. Scalar Arguments ..
304 CHARACTER HOWMNY, SIDE
305 INTEGER INFO, LDP, LDS, LDVL, LDVR, M, MM, N
307 * .. Array Arguments ..
309 REAL P( LDP, * ), S( LDS, * ), VL( LDVL, * ),
310 $ VR( LDVR, * ), WORK( * )
314 * =====================================================================
317 REAL ZERO, ONE, SAFETY
318 PARAMETER ( ZERO = 0.0E+0, ONE = 1.0E+0,
321 * .. Local Scalars ..
322 LOGICAL COMPL, COMPR, IL2BY2, ILABAD, ILALL, ILBACK,
323 $ ILBBAD, ILCOMP, ILCPLX, LSA, LSB
324 INTEGER I, IBEG, IEIG, IEND, IHWMNY, IINFO, IM, ISIDE,
325 $ J, JA, JC, JE, JR, JW, NA, NW
326 REAL ACOEF, ACOEFA, ANORM, ASCALE, BCOEFA, BCOEFI,
327 $ BCOEFR, BIG, BIGNUM, BNORM, BSCALE, CIM2A,
328 $ CIM2B, CIMAGA, CIMAGB, CRE2A, CRE2B, CREALA,
329 $ CREALB, DMIN, SAFMIN, SALFAR, SBETA, SCALE,
330 $ SMALL, TEMP, TEMP2, TEMP2I, TEMP2R, ULP, XMAX,
334 REAL BDIAG( 2 ), SUM( 2, 2 ), SUMS( 2, 2 ),
337 * .. External Functions ..
340 EXTERNAL LSAME, SLAMCH
342 * .. External Subroutines ..
343 EXTERNAL SGEMV, SLABAD, SLACPY, SLAG2, SLALN2, XERBLA
345 * .. Intrinsic Functions ..
346 INTRINSIC ABS, MAX, MIN
348 * .. Executable Statements ..
350 * Decode and Test the input parameters
352 IF( LSAME( HOWMNY, 'A' ) ) THEN
356 ELSE IF( LSAME( HOWMNY, 'S' ) ) THEN
360 ELSE IF( LSAME( HOWMNY, 'B' ) ) THEN
369 IF( LSAME( SIDE, 'R' ) ) THEN
373 ELSE IF( LSAME( SIDE, 'L' ) ) THEN
377 ELSE IF( LSAME( SIDE, 'B' ) ) THEN
386 IF( ISIDE.LT.0 ) THEN
388 ELSE IF( IHWMNY.LT.0 ) THEN
390 ELSE IF( N.LT.0 ) THEN
392 ELSE IF( LDS.LT.MAX( 1, N ) ) THEN
394 ELSE IF( LDP.LT.MAX( 1, N ) ) THEN
398 CALL XERBLA( 'STGEVC', -INFO )
402 * Count the number of eigenvectors to be computed
404 IF( .NOT.ILALL ) THEN
413 IF( S( J+1, J ).NE.ZERO )
417 IF( SELECT( J ) .OR. SELECT( J+1 ) )
428 * Check 2-by-2 diagonal blocks of A, B
433 IF( S( J+1, J ).NE.ZERO ) THEN
434 IF( P( J, J ).EQ.ZERO .OR. P( J+1, J+1 ).EQ.ZERO .OR.
435 $ P( J, J+1 ).NE.ZERO )ILBBAD = .TRUE.
437 IF( S( J+2, J+1 ).NE.ZERO )
445 ELSE IF( ILBBAD ) THEN
447 ELSE IF( COMPL .AND. LDVL.LT.N .OR. LDVL.LT.1 ) THEN
449 ELSE IF( COMPR .AND. LDVR.LT.N .OR. LDVR.LT.1 ) THEN
451 ELSE IF( MM.LT.IM ) THEN
455 CALL XERBLA( 'STGEVC', -INFO )
459 * Quick return if possible
467 SAFMIN = SLAMCH( 'Safe minimum' )
469 CALL SLABAD( SAFMIN, BIG )
470 ULP = SLAMCH( 'Epsilon' )*SLAMCH( 'Base' )
471 SMALL = SAFMIN*N / ULP
473 BIGNUM = ONE / ( SAFMIN*N )
475 * Compute the 1-norm of each column of the strictly upper triangular
476 * part (i.e., excluding all elements belonging to the diagonal
477 * blocks) of A and B to check for possible overflow in the
480 ANORM = ABS( S( 1, 1 ) )
482 $ ANORM = ANORM + ABS( S( 2, 1 ) )
483 BNORM = ABS( P( 1, 1 ) )
490 IF( S( J, J-1 ).EQ.ZERO ) THEN
496 TEMP = TEMP + ABS( S( I, J ) )
497 TEMP2 = TEMP2 + ABS( P( I, J ) )
501 DO 40 I = IEND + 1, MIN( J+1, N )
502 TEMP = TEMP + ABS( S( I, J ) )
503 TEMP2 = TEMP2 + ABS( P( I, J ) )
505 ANORM = MAX( ANORM, TEMP )
506 BNORM = MAX( BNORM, TEMP2 )
509 ASCALE = ONE / MAX( ANORM, SAFMIN )
510 BSCALE = ONE / MAX( BNORM, SAFMIN )
517 * Main loop over eigenvalues
522 * Skip this iteration if (a) HOWMNY='S' and SELECT=.FALSE., or
523 * (b) this would be the second of a complex pair.
524 * Check for complex eigenvalue, so as to be sure of which
525 * entry(-ies) of SELECT to look at.
533 IF( S( JE+1, JE ).NE.ZERO ) THEN
540 ELSE IF( ILCPLX ) THEN
541 ILCOMP = SELECT( JE ) .OR. SELECT( JE+1 )
543 ILCOMP = SELECT( JE )
548 * Decide if (a) singular pencil, (b) real eigenvalue, or
549 * (c) complex eigenvalue.
551 IF( .NOT.ILCPLX ) THEN
552 IF( ABS( S( JE, JE ) ).LE.SAFMIN .AND.
553 $ ABS( P( JE, JE ) ).LE.SAFMIN ) THEN
555 * Singular matrix pencil -- return unit eigenvector
559 VL( JR, IEIG ) = ZERO
561 VL( IEIG, IEIG ) = ONE
569 WORK( 2*N+JR ) = ZERO
572 * Compute coefficients in ( a A - b B ) y = 0
574 * b is BCOEFR + i*BCOEFI
576 IF( .NOT.ILCPLX ) THEN
580 TEMP = ONE / MAX( ABS( S( JE, JE ) )*ASCALE,
581 $ ABS( P( JE, JE ) )*BSCALE, SAFMIN )
582 SALFAR = ( TEMP*S( JE, JE ) )*ASCALE
583 SBETA = ( TEMP*P( JE, JE ) )*BSCALE
585 BCOEFR = SALFAR*BSCALE
588 * Scale to avoid underflow
591 LSA = ABS( SBETA ).GE.SAFMIN .AND. ABS( ACOEF ).LT.SMALL
592 LSB = ABS( SALFAR ).GE.SAFMIN .AND. ABS( BCOEFR ).LT.
595 $ SCALE = ( SMALL / ABS( SBETA ) )*MIN( ANORM, BIG )
597 $ SCALE = MAX( SCALE, ( SMALL / ABS( SALFAR ) )*
598 $ MIN( BNORM, BIG ) )
599 IF( LSA .OR. LSB ) THEN
600 SCALE = MIN( SCALE, ONE /
601 $ ( SAFMIN*MAX( ONE, ABS( ACOEF ),
602 $ ABS( BCOEFR ) ) ) )
604 ACOEF = ASCALE*( SCALE*SBETA )
609 BCOEFR = BSCALE*( SCALE*SALFAR )
611 BCOEFR = SCALE*BCOEFR
614 ACOEFA = ABS( ACOEF )
615 BCOEFA = ABS( BCOEFR )
617 * First component is 1
625 CALL SLAG2( S( JE, JE ), LDS, P( JE, JE ), LDP,
626 $ SAFMIN*SAFETY, ACOEF, TEMP, BCOEFR, TEMP2,
629 IF( BCOEFI.EQ.ZERO ) THEN
634 * Scale to avoid over/underflow
636 ACOEFA = ABS( ACOEF )
637 BCOEFA = ABS( BCOEFR ) + ABS( BCOEFI )
639 IF( ACOEFA*ULP.LT.SAFMIN .AND. ACOEFA.GE.SAFMIN )
640 $ SCALE = ( SAFMIN / ULP ) / ACOEFA
641 IF( BCOEFA*ULP.LT.SAFMIN .AND. BCOEFA.GE.SAFMIN )
642 $ SCALE = MAX( SCALE, ( SAFMIN / ULP ) / BCOEFA )
643 IF( SAFMIN*ACOEFA.GT.ASCALE )
644 $ SCALE = ASCALE / ( SAFMIN*ACOEFA )
645 IF( SAFMIN*BCOEFA.GT.BSCALE )
646 $ SCALE = MIN( SCALE, BSCALE / ( SAFMIN*BCOEFA ) )
647 IF( SCALE.NE.ONE ) THEN
649 ACOEFA = ABS( ACOEF )
650 BCOEFR = SCALE*BCOEFR
651 BCOEFI = SCALE*BCOEFI
652 BCOEFA = ABS( BCOEFR ) + ABS( BCOEFI )
655 * Compute first two components of eigenvector
657 TEMP = ACOEF*S( JE+1, JE )
658 TEMP2R = ACOEF*S( JE, JE ) - BCOEFR*P( JE, JE )
659 TEMP2I = -BCOEFI*P( JE, JE )
660 IF( ABS( TEMP ).GT.ABS( TEMP2R )+ABS( TEMP2I ) ) THEN
662 WORK( 3*N+JE ) = ZERO
663 WORK( 2*N+JE+1 ) = -TEMP2R / TEMP
664 WORK( 3*N+JE+1 ) = -TEMP2I / TEMP
666 WORK( 2*N+JE+1 ) = ONE
667 WORK( 3*N+JE+1 ) = ZERO
668 TEMP = ACOEF*S( JE, JE+1 )
669 WORK( 2*N+JE ) = ( BCOEFR*P( JE+1, JE+1 )-ACOEF*
670 $ S( JE+1, JE+1 ) ) / TEMP
671 WORK( 3*N+JE ) = BCOEFI*P( JE+1, JE+1 ) / TEMP
673 XMAX = MAX( ABS( WORK( 2*N+JE ) )+ABS( WORK( 3*N+JE ) ),
674 $ ABS( WORK( 2*N+JE+1 ) )+ABS( WORK( 3*N+JE+1 ) ) )
677 DMIN = MAX( ULP*ACOEFA*ANORM, ULP*BCOEFA*BNORM, SAFMIN )
680 * Triangular solve of (a A - b B) y = 0
683 * (rowwise in (a A - b B) , or columnwise in (a A - b B) )
687 DO 160 J = JE + NW, N
694 BDIAG( 1 ) = P( J, J )
696 IF( S( J+1, J ).NE.ZERO ) THEN
698 BDIAG( 2 ) = P( J+1, J+1 )
703 * Check whether scaling is necessary for dot products
705 XSCALE = ONE / MAX( ONE, XMAX )
706 TEMP = MAX( WORK( J ), WORK( N+J ),
707 $ ACOEFA*WORK( J )+BCOEFA*WORK( N+J ) )
709 $ TEMP = MAX( TEMP, WORK( J+1 ), WORK( N+J+1 ),
710 $ ACOEFA*WORK( J+1 )+BCOEFA*WORK( N+J+1 ) )
711 IF( TEMP.GT.BIGNUM*XSCALE ) THEN
714 WORK( ( JW+2 )*N+JR ) = XSCALE*
715 $ WORK( ( JW+2 )*N+JR )
721 * Compute dot products
724 * SUM = sum conjg( a*S(k,j) - b*P(k,j) )*x(k)
727 * To reduce the op count, this is done as
730 * a*conjg( sum S(k,j)*x(k) ) - b*conjg( sum P(k,j)*x(k) )
733 * which may cause underflow problems if A or B are close
734 * to underflow. (E.g., less than SMALL.)
739 SUMS( JA, JW ) = ZERO
740 SUMP( JA, JW ) = ZERO
742 DO 100 JR = JE, J - 1
743 SUMS( JA, JW ) = SUMS( JA, JW ) +
745 $ WORK( ( JW+1 )*N+JR )
746 SUMP( JA, JW ) = SUMP( JA, JW ) +
748 $ WORK( ( JW+1 )*N+JR )
755 SUM( JA, 1 ) = -ACOEF*SUMS( JA, 1 ) +
756 $ BCOEFR*SUMP( JA, 1 ) -
757 $ BCOEFI*SUMP( JA, 2 )
758 SUM( JA, 2 ) = -ACOEF*SUMS( JA, 2 ) +
759 $ BCOEFR*SUMP( JA, 2 ) +
760 $ BCOEFI*SUMP( JA, 1 )
762 SUM( JA, 1 ) = -ACOEF*SUMS( JA, 1 ) +
763 $ BCOEFR*SUMP( JA, 1 )
768 * Solve ( a A - b B ) y = SUM(,)
769 * with scaling and perturbation of the denominator
771 CALL SLALN2( .TRUE., NA, NW, DMIN, ACOEF, S( J, J ), LDS,
772 $ BDIAG( 1 ), BDIAG( 2 ), SUM, 2, BCOEFR,
773 $ BCOEFI, WORK( 2*N+J ), N, SCALE, TEMP,
775 IF( SCALE.LT.ONE ) THEN
776 DO 150 JW = 0, NW - 1
777 DO 140 JR = JE, J - 1
778 WORK( ( JW+2 )*N+JR ) = SCALE*
779 $ WORK( ( JW+2 )*N+JR )
784 XMAX = MAX( XMAX, TEMP )
787 * Copy eigenvector to VL, back transforming if
792 DO 170 JW = 0, NW - 1
793 CALL SGEMV( 'N', N, N+1-JE, ONE, VL( 1, JE ), LDVL,
794 $ WORK( ( JW+2 )*N+JE ), 1, ZERO,
795 $ WORK( ( JW+4 )*N+1 ), 1 )
797 CALL SLACPY( ' ', N, NW, WORK( 4*N+1 ), N, VL( 1, JE ),
801 CALL SLACPY( ' ', N, NW, WORK( 2*N+1 ), N, VL( 1, IEIG ),
811 XMAX = MAX( XMAX, ABS( VL( J, IEIG ) )+
812 $ ABS( VL( J, IEIG+1 ) ) )
816 XMAX = MAX( XMAX, ABS( VL( J, IEIG ) ) )
820 IF( XMAX.GT.SAFMIN ) THEN
823 DO 210 JW = 0, NW - 1
825 VL( JR, IEIG+JW ) = XSCALE*VL( JR, IEIG+JW )
839 * Main loop over eigenvalues
844 * Skip this iteration if (a) HOWMNY='S' and SELECT=.FALSE., or
845 * (b) this would be the second of a complex pair.
846 * Check for complex eigenvalue, so as to be sure of which
847 * entry(-ies) of SELECT to look at -- if complex, SELECT(JE)
849 * If this is a complex pair, the 2-by-2 diagonal block
850 * corresponding to the eigenvalue is in rows/columns JE-1:JE
858 IF( S( JE, JE-1 ).NE.ZERO ) THEN
865 ELSE IF( ILCPLX ) THEN
866 ILCOMP = SELECT( JE ) .OR. SELECT( JE-1 )
868 ILCOMP = SELECT( JE )
873 * Decide if (a) singular pencil, (b) real eigenvalue, or
874 * (c) complex eigenvalue.
876 IF( .NOT.ILCPLX ) THEN
877 IF( ABS( S( JE, JE ) ).LE.SAFMIN .AND.
878 $ ABS( P( JE, JE ) ).LE.SAFMIN ) THEN
880 * Singular matrix pencil -- unit eigenvector
884 VR( JR, IEIG ) = ZERO
886 VR( IEIG, IEIG ) = ONE
893 DO 250 JW = 0, NW - 1
895 WORK( ( JW+2 )*N+JR ) = ZERO
899 * Compute coefficients in ( a A - b B ) x = 0
901 * b is BCOEFR + i*BCOEFI
903 IF( .NOT.ILCPLX ) THEN
907 TEMP = ONE / MAX( ABS( S( JE, JE ) )*ASCALE,
908 $ ABS( P( JE, JE ) )*BSCALE, SAFMIN )
909 SALFAR = ( TEMP*S( JE, JE ) )*ASCALE
910 SBETA = ( TEMP*P( JE, JE ) )*BSCALE
912 BCOEFR = SALFAR*BSCALE
915 * Scale to avoid underflow
918 LSA = ABS( SBETA ).GE.SAFMIN .AND. ABS( ACOEF ).LT.SMALL
919 LSB = ABS( SALFAR ).GE.SAFMIN .AND. ABS( BCOEFR ).LT.
922 $ SCALE = ( SMALL / ABS( SBETA ) )*MIN( ANORM, BIG )
924 $ SCALE = MAX( SCALE, ( SMALL / ABS( SALFAR ) )*
925 $ MIN( BNORM, BIG ) )
926 IF( LSA .OR. LSB ) THEN
927 SCALE = MIN( SCALE, ONE /
928 $ ( SAFMIN*MAX( ONE, ABS( ACOEF ),
929 $ ABS( BCOEFR ) ) ) )
931 ACOEF = ASCALE*( SCALE*SBETA )
936 BCOEFR = BSCALE*( SCALE*SALFAR )
938 BCOEFR = SCALE*BCOEFR
941 ACOEFA = ABS( ACOEF )
942 BCOEFA = ABS( BCOEFR )
944 * First component is 1
949 * Compute contribution from column JE of A and B to sum
950 * (See "Further Details", above.)
952 DO 260 JR = 1, JE - 1
953 WORK( 2*N+JR ) = BCOEFR*P( JR, JE ) -
960 CALL SLAG2( S( JE-1, JE-1 ), LDS, P( JE-1, JE-1 ), LDP,
961 $ SAFMIN*SAFETY, ACOEF, TEMP, BCOEFR, TEMP2,
963 IF( BCOEFI.EQ.ZERO ) THEN
968 * Scale to avoid over/underflow
970 ACOEFA = ABS( ACOEF )
971 BCOEFA = ABS( BCOEFR ) + ABS( BCOEFI )
973 IF( ACOEFA*ULP.LT.SAFMIN .AND. ACOEFA.GE.SAFMIN )
974 $ SCALE = ( SAFMIN / ULP ) / ACOEFA
975 IF( BCOEFA*ULP.LT.SAFMIN .AND. BCOEFA.GE.SAFMIN )
976 $ SCALE = MAX( SCALE, ( SAFMIN / ULP ) / BCOEFA )
977 IF( SAFMIN*ACOEFA.GT.ASCALE )
978 $ SCALE = ASCALE / ( SAFMIN*ACOEFA )
979 IF( SAFMIN*BCOEFA.GT.BSCALE )
980 $ SCALE = MIN( SCALE, BSCALE / ( SAFMIN*BCOEFA ) )
981 IF( SCALE.NE.ONE ) THEN
983 ACOEFA = ABS( ACOEF )
984 BCOEFR = SCALE*BCOEFR
985 BCOEFI = SCALE*BCOEFI
986 BCOEFA = ABS( BCOEFR ) + ABS( BCOEFI )
989 * Compute first two components of eigenvector
990 * and contribution to sums
992 TEMP = ACOEF*S( JE, JE-1 )
993 TEMP2R = ACOEF*S( JE, JE ) - BCOEFR*P( JE, JE )
994 TEMP2I = -BCOEFI*P( JE, JE )
995 IF( ABS( TEMP ).GE.ABS( TEMP2R )+ABS( TEMP2I ) ) THEN
997 WORK( 3*N+JE ) = ZERO
998 WORK( 2*N+JE-1 ) = -TEMP2R / TEMP
999 WORK( 3*N+JE-1 ) = -TEMP2I / TEMP
1001 WORK( 2*N+JE-1 ) = ONE
1002 WORK( 3*N+JE-1 ) = ZERO
1003 TEMP = ACOEF*S( JE-1, JE )
1004 WORK( 2*N+JE ) = ( BCOEFR*P( JE-1, JE-1 )-ACOEF*
1005 $ S( JE-1, JE-1 ) ) / TEMP
1006 WORK( 3*N+JE ) = BCOEFI*P( JE-1, JE-1 ) / TEMP
1009 XMAX = MAX( ABS( WORK( 2*N+JE ) )+ABS( WORK( 3*N+JE ) ),
1010 $ ABS( WORK( 2*N+JE-1 ) )+ABS( WORK( 3*N+JE-1 ) ) )
1012 * Compute contribution from columns JE and JE-1
1013 * of A and B to the sums.
1015 CREALA = ACOEF*WORK( 2*N+JE-1 )
1016 CIMAGA = ACOEF*WORK( 3*N+JE-1 )
1017 CREALB = BCOEFR*WORK( 2*N+JE-1 ) -
1018 $ BCOEFI*WORK( 3*N+JE-1 )
1019 CIMAGB = BCOEFI*WORK( 2*N+JE-1 ) +
1020 $ BCOEFR*WORK( 3*N+JE-1 )
1021 CRE2A = ACOEF*WORK( 2*N+JE )
1022 CIM2A = ACOEF*WORK( 3*N+JE )
1023 CRE2B = BCOEFR*WORK( 2*N+JE ) - BCOEFI*WORK( 3*N+JE )
1024 CIM2B = BCOEFI*WORK( 2*N+JE ) + BCOEFR*WORK( 3*N+JE )
1025 DO 270 JR = 1, JE - 2
1026 WORK( 2*N+JR ) = -CREALA*S( JR, JE-1 ) +
1027 $ CREALB*P( JR, JE-1 ) -
1028 $ CRE2A*S( JR, JE ) + CRE2B*P( JR, JE )
1029 WORK( 3*N+JR ) = -CIMAGA*S( JR, JE-1 ) +
1030 $ CIMAGB*P( JR, JE-1 ) -
1031 $ CIM2A*S( JR, JE ) + CIM2B*P( JR, JE )
1035 DMIN = MAX( ULP*ACOEFA*ANORM, ULP*BCOEFA*BNORM, SAFMIN )
1037 * Columnwise triangular solve of (a A - b B) x = 0
1040 DO 370 J = JE - NW, 1, -1
1042 * If a 2-by-2 block, is in position j-1:j, wait until
1043 * next iteration to process it (when it will be j:j+1)
1045 IF( .NOT.IL2BY2 .AND. J.GT.1 ) THEN
1046 IF( S( J, J-1 ).NE.ZERO ) THEN
1051 BDIAG( 1 ) = P( J, J )
1054 BDIAG( 2 ) = P( J+1, J+1 )
1059 * Compute x(j) (and x(j+1), if 2-by-2 block)
1061 CALL SLALN2( .FALSE., NA, NW, DMIN, ACOEF, S( J, J ),
1062 $ LDS, BDIAG( 1 ), BDIAG( 2 ), WORK( 2*N+J ),
1063 $ N, BCOEFR, BCOEFI, SUM, 2, SCALE, TEMP,
1065 IF( SCALE.LT.ONE ) THEN
1067 DO 290 JW = 0, NW - 1
1069 WORK( ( JW+2 )*N+JR ) = SCALE*
1070 $ WORK( ( JW+2 )*N+JR )
1074 XMAX = MAX( SCALE*XMAX, TEMP )
1078 WORK( ( JW+1 )*N+J+JA-1 ) = SUM( JA, JW )
1082 * w = w + x(j)*(a S(*,j) - b P(*,j) ) with scaling
1086 * Check whether scaling is necessary for sum.
1088 XSCALE = ONE / MAX( ONE, XMAX )
1089 TEMP = ACOEFA*WORK( J ) + BCOEFA*WORK( N+J )
1091 $ TEMP = MAX( TEMP, ACOEFA*WORK( J+1 )+BCOEFA*
1093 TEMP = MAX( TEMP, ACOEFA, BCOEFA )
1094 IF( TEMP.GT.BIGNUM*XSCALE ) THEN
1096 DO 330 JW = 0, NW - 1
1098 WORK( ( JW+2 )*N+JR ) = XSCALE*
1099 $ WORK( ( JW+2 )*N+JR )
1105 * Compute the contributions of the off-diagonals of
1106 * column j (and j+1, if 2-by-2 block) of A and B to the
1112 CREALA = ACOEF*WORK( 2*N+J+JA-1 )
1113 CIMAGA = ACOEF*WORK( 3*N+J+JA-1 )
1114 CREALB = BCOEFR*WORK( 2*N+J+JA-1 ) -
1115 $ BCOEFI*WORK( 3*N+J+JA-1 )
1116 CIMAGB = BCOEFI*WORK( 2*N+J+JA-1 ) +
1117 $ BCOEFR*WORK( 3*N+J+JA-1 )
1118 DO 340 JR = 1, J - 1
1119 WORK( 2*N+JR ) = WORK( 2*N+JR ) -
1120 $ CREALA*S( JR, J+JA-1 ) +
1121 $ CREALB*P( JR, J+JA-1 )
1122 WORK( 3*N+JR ) = WORK( 3*N+JR ) -
1123 $ CIMAGA*S( JR, J+JA-1 ) +
1124 $ CIMAGB*P( JR, J+JA-1 )
1127 CREALA = ACOEF*WORK( 2*N+J+JA-1 )
1128 CREALB = BCOEFR*WORK( 2*N+J+JA-1 )
1129 DO 350 JR = 1, J - 1
1130 WORK( 2*N+JR ) = WORK( 2*N+JR ) -
1131 $ CREALA*S( JR, J+JA-1 ) +
1132 $ CREALB*P( JR, J+JA-1 )
1141 * Copy eigenvector to VR, back transforming if
1147 DO 410 JW = 0, NW - 1
1149 WORK( ( JW+4 )*N+JR ) = WORK( ( JW+2 )*N+1 )*
1153 * A series of compiler directives to defeat
1154 * vectorization for the next loop
1159 WORK( ( JW+4 )*N+JR ) = WORK( ( JW+4 )*N+JR ) +
1160 $ WORK( ( JW+2 )*N+JC )*VR( JR, JC )
1165 DO 430 JW = 0, NW - 1
1167 VR( JR, IEIG+JW ) = WORK( ( JW+4 )*N+JR )
1173 DO 450 JW = 0, NW - 1
1175 VR( JR, IEIG+JW ) = WORK( ( JW+2 )*N+JR )
1187 XMAX = MAX( XMAX, ABS( VR( J, IEIG ) )+
1188 $ ABS( VR( J, IEIG+1 ) ) )
1192 XMAX = MAX( XMAX, ABS( VR( J, IEIG ) ) )
1196 IF( XMAX.GT.SAFMIN ) THEN
1198 DO 490 JW = 0, NW - 1
1200 VR( JR, IEIG+JW ) = XSCALE*VR( JR, IEIG+JW )