3 * =========== DOCUMENTATION ===========
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
9 *> Download STREVC3 + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/strevc3.f">
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/strevc3.f">
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/strevc3.f">
21 * SUBROUTINE STREVC3( SIDE, HOWMNY, SELECT, N, T, LDT, VL, LDVL,
22 * VR, LDVR, MM, M, WORK, LWORK, INFO )
24 * .. Scalar Arguments ..
25 * CHARACTER HOWMNY, SIDE
26 * INTEGER INFO, LDT, LDVL, LDVR, LWORK, M, MM, N
28 * .. Array Arguments ..
30 * REAL T( LDT, * ), VL( LDVL, * ), VR( LDVR, * ),
40 *> STREVC3 computes some or all of the right and/or left eigenvectors of
41 *> a real upper quasi-triangular matrix T.
42 *> Matrices of this type are produced by the Schur factorization of
43 *> a real general matrix: A = Q*T*Q**T, as computed by SHSEQR.
45 *> The right eigenvector x and the left eigenvector y of T corresponding
46 *> to an eigenvalue w are defined by:
48 *> T*x = w*x, (y**T)*T = w*(y**T)
50 *> where y**T denotes the transpose of the vector y.
51 *> The eigenvalues are not input to this routine, but are read directly
52 *> from the diagonal blocks of T.
54 *> This routine returns the matrices X and/or Y of right and left
55 *> eigenvectors of T, or the products Q*X and/or Q*Y, where Q is an
56 *> input matrix. If Q is the orthogonal factor that reduces a matrix
57 *> A to Schur form T, then Q*X and Q*Y are the matrices of right and
58 *> left eigenvectors of A.
60 *> This uses a Level 3 BLAS version of the back transformation.
68 *> SIDE is CHARACTER*1
69 *> = 'R': compute right eigenvectors only;
70 *> = 'L': compute left eigenvectors only;
71 *> = 'B': compute both right and left eigenvectors.
76 *> HOWMNY is CHARACTER*1
77 *> = 'A': compute all right and/or left eigenvectors;
78 *> = 'B': compute all right and/or left eigenvectors,
79 *> backtransformed by the matrices in VR and/or VL;
80 *> = 'S': compute selected right and/or left eigenvectors,
81 *> as indicated by the logical array SELECT.
84 *> \param[in,out] SELECT
86 *> SELECT is LOGICAL array, dimension (N)
87 *> If HOWMNY = 'S', SELECT specifies the eigenvectors to be
89 *> If w(j) is a real eigenvalue, the corresponding real
90 *> eigenvector is computed if SELECT(j) is .TRUE..
91 *> If w(j) and w(j+1) are the real and imaginary parts of a
92 *> complex eigenvalue, the corresponding complex eigenvector is
93 *> computed if either SELECT(j) or SELECT(j+1) is .TRUE., and
94 *> on exit SELECT(j) is set to .TRUE. and SELECT(j+1) is set to
96 *> Not referenced if HOWMNY = 'A' or 'B'.
102 *> The order of the matrix T. N >= 0.
107 *> T is REAL array, dimension (LDT,N)
108 *> The upper quasi-triangular matrix T in Schur canonical form.
114 *> The leading dimension of the array T. LDT >= max(1,N).
119 *> VL is REAL array, dimension (LDVL,MM)
120 *> On entry, if SIDE = 'L' or 'B' and HOWMNY = 'B', VL must
121 *> contain an N-by-N matrix Q (usually the orthogonal matrix Q
122 *> of Schur vectors returned by SHSEQR).
123 *> On exit, if SIDE = 'L' or 'B', VL contains:
124 *> if HOWMNY = 'A', the matrix Y of left eigenvectors of T;
125 *> if HOWMNY = 'B', the matrix Q*Y;
126 *> if HOWMNY = 'S', the left eigenvectors of T specified by
127 *> SELECT, stored consecutively in the columns
128 *> of VL, in the same order as their
130 *> A complex eigenvector corresponding to a complex eigenvalue
131 *> is stored in two consecutive columns, the first holding the
132 *> real part, and the second the imaginary part.
133 *> Not referenced if SIDE = 'R'.
139 *> The leading dimension of the array VL.
140 *> LDVL >= 1, and if SIDE = 'L' or 'B', LDVL >= N.
145 *> VR is REAL array, dimension (LDVR,MM)
146 *> On entry, if SIDE = 'R' or 'B' and HOWMNY = 'B', VR must
147 *> contain an N-by-N matrix Q (usually the orthogonal matrix Q
148 *> of Schur vectors returned by SHSEQR).
149 *> On exit, if SIDE = 'R' or 'B', VR contains:
150 *> if HOWMNY = 'A', the matrix X of right eigenvectors of T;
151 *> if HOWMNY = 'B', the matrix Q*X;
152 *> if HOWMNY = 'S', the right eigenvectors of T specified by
153 *> SELECT, stored consecutively in the columns
154 *> of VR, in the same order as their
156 *> A complex eigenvector corresponding to a complex eigenvalue
157 *> is stored in two consecutive columns, the first holding the
158 *> real part and the second the imaginary part.
159 *> Not referenced if SIDE = 'L'.
165 *> The leading dimension of the array VR.
166 *> LDVR >= 1, and if SIDE = 'R' or 'B', LDVR >= N.
172 *> The number of columns in the arrays VL and/or VR. MM >= M.
178 *> The number of columns in the arrays VL and/or VR actually
179 *> used to store the eigenvectors.
180 *> If HOWMNY = 'A' or 'B', M is set to N.
181 *> Each selected real eigenvector occupies one column and each
182 *> selected complex eigenvector occupies two columns.
187 *> WORK is REAL array, dimension (MAX(1,LWORK))
193 *> The dimension of array WORK. LWORK >= max(1,3*N).
194 *> For optimum performance, LWORK >= N + 2*N*NB, where NB is
195 *> the optimal blocksize.
197 *> If LWORK = -1, then a workspace query is assumed; the routine
198 *> only calculates the optimal size of the WORK array, returns
199 *> this value as the first entry of the WORK array, and no error
200 *> message related to LWORK is issued by XERBLA.
206 *> = 0: successful exit
207 *> < 0: if INFO = -i, the i-th argument had an illegal value
213 *> \author Univ. of Tennessee
214 *> \author Univ. of California Berkeley
215 *> \author Univ. of Colorado Denver
218 *> \date November 2011
220 * @generated from dtrevc3.f, fortran d -> s, Tue Apr 19 01:47:44 2016
222 *> \ingroup realOTHERcomputational
224 *> \par Further Details:
225 * =====================
229 *> The algorithm used in this program is basically backward (forward)
230 *> substitution, with scaling to make the the code robust against
231 *> possible overflow.
233 *> Each eigenvector is normalized so that the element of largest
234 *> magnitude has magnitude 1; here the magnitude of a complex number
235 *> (x,y) is taken to be |x| + |y|.
238 * =====================================================================
239 SUBROUTINE STREVC3( SIDE, HOWMNY, SELECT, N, T, LDT, VL, LDVL,
240 $ VR, LDVR, MM, M, WORK, LWORK, INFO )
243 * -- LAPACK computational routine (version 3.4.0) --
244 * -- LAPACK is a software package provided by Univ. of Tennessee, --
245 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
248 * .. Scalar Arguments ..
249 CHARACTER HOWMNY, SIDE
250 INTEGER INFO, LDT, LDVL, LDVR, LWORK, M, MM, N
252 * .. Array Arguments ..
254 REAL T( LDT, * ), VL( LDVL, * ), VR( LDVR, * ),
258 * =====================================================================
262 PARAMETER ( ZERO = 0.0E+0, ONE = 1.0E+0 )
264 PARAMETER ( NBMIN = 8, NBMAX = 128 )
266 * .. Local Scalars ..
267 LOGICAL ALLV, BOTHV, LEFTV, LQUERY, OVER, PAIR,
269 INTEGER I, IERR, II, IP, IS, J, J1, J2, JNXT, K, KI,
270 $ IV, MAXWRK, NB, KI2
271 REAL BETA, BIGNUM, EMAX, OVFL, REC, REMAX, SCALE,
272 $ SMIN, SMLNUM, ULP, UNFL, VCRIT, VMAX, WI, WR,
275 * .. External Functions ..
277 INTEGER ISAMAX, ILAENV
279 EXTERNAL LSAME, ISAMAX, ILAENV, SDOT, SLAMCH
281 * .. External Subroutines ..
282 EXTERNAL SAXPY, SCOPY, SGEMV, SLALN2, SSCAL, XERBLA
284 * .. Intrinsic Functions ..
285 INTRINSIC ABS, MAX, SQRT
289 INTEGER ISCOMPLEX( NBMAX )
291 * .. Executable Statements ..
293 * Decode and test the input parameters
295 BOTHV = LSAME( SIDE, 'B' )
296 RIGHTV = LSAME( SIDE, 'R' ) .OR. BOTHV
297 LEFTV = LSAME( SIDE, 'L' ) .OR. BOTHV
299 ALLV = LSAME( HOWMNY, 'A' )
300 OVER = LSAME( HOWMNY, 'B' )
301 SOMEV = LSAME( HOWMNY, 'S' )
304 NB = ILAENV( 1, 'STREVC', SIDE // HOWMNY, N, -1, -1, -1 )
307 LQUERY = ( LWORK.EQ.-1 )
308 IF( .NOT.RIGHTV .AND. .NOT.LEFTV ) THEN
310 ELSE IF( .NOT.ALLV .AND. .NOT.OVER .AND. .NOT.SOMEV ) THEN
312 ELSE IF( N.LT.0 ) THEN
314 ELSE IF( LDT.LT.MAX( 1, N ) ) THEN
316 ELSE IF( LDVL.LT.1 .OR. ( LEFTV .AND. LDVL.LT.N ) ) THEN
318 ELSE IF( LDVR.LT.1 .OR. ( RIGHTV .AND. LDVR.LT.N ) ) THEN
320 ELSE IF( LWORK.LT.MAX( 1, 3*N ) .AND. .NOT.LQUERY ) THEN
324 * Set M to the number of columns required to store the selected
325 * eigenvectors, standardize the array SELECT if necessary, and
334 SELECT( J ) = .FALSE.
337 IF( T( J+1, J ).EQ.ZERO ) THEN
342 IF( SELECT( J ) .OR. SELECT( J+1 ) ) THEN
362 CALL XERBLA( 'STREVC3', -INFO )
364 ELSE IF( LQUERY ) THEN
368 * Quick return if possible.
373 * Use blocked version of back-transformation if sufficient workspace.
374 * Zero-out the workspace to avoid potential NaN propagation.
376 IF( OVER .AND. LWORK .GE. N + 2*N*NBMIN ) THEN
377 NB = (LWORK - N) / (2*N)
378 NB = MIN( NB, NBMAX )
379 CALL SLASET( 'F', N, 1+2*NB, ZERO, ZERO, WORK, N )
384 * Set the constants to control overflow.
386 UNFL = SLAMCH( 'Safe minimum' )
388 CALL SLABAD( UNFL, OVFL )
389 ULP = SLAMCH( 'Precision' )
390 SMLNUM = UNFL*( N / ULP )
391 BIGNUM = ( ONE-ULP ) / SMLNUM
393 * Compute 1-norm of each column of strictly upper triangular
394 * part of T to control overflow in triangular solver.
400 WORK( J ) = WORK( J ) + ABS( T( I, J ) )
404 * Index IP is used to specify the real or complex eigenvalue:
405 * IP = 0, real eigenvalue,
406 * 1, first of conjugate complex pair: (wr,wi)
407 * -1, second of conjugate complex pair: (wr,wi)
408 * ISCOMPLEX array stores IP for each column in current block.
412 * ============================================================
413 * Compute right eigenvectors.
415 * IV is index of column in current block.
416 * For complex right vector, uses IV-1 for real part and IV for complex part.
417 * Non-blocked version always uses IV=2;
418 * blocked version starts with IV=NB, goes down to 1 or 2.
419 * (Note the "0-th" column is used for 1-norms computed above.)
429 * previous iteration (ki+1) was second of conjugate pair,
430 * so this ki is first of conjugate pair; skip to end of loop
433 ELSE IF( KI.EQ.1 ) THEN
434 * last column, so this ki must be real eigenvalue
436 ELSE IF( T( KI, KI-1 ).EQ.ZERO ) THEN
437 * zero on sub-diagonal, so this ki is real eigenvalue
440 * non-zero on sub-diagonal, so this ki is second of conjugate pair
446 IF( .NOT.SELECT( KI ) )
449 IF( .NOT.SELECT( KI-1 ) )
454 * Compute the KI-th eigenvalue (WR,WI).
459 $ WI = SQRT( ABS( T( KI, KI-1 ) ) )*
460 $ SQRT( ABS( T( KI-1, KI ) ) )
461 SMIN = MAX( ULP*( ABS( WR )+ABS( WI ) ), SMLNUM )
465 * --------------------------------------------------------
466 * Real right eigenvector
468 WORK( KI + IV*N ) = ONE
470 * Form right-hand side.
473 WORK( K + IV*N ) = -T( K, KI )
476 * Solve upper quasi-triangular system:
477 * [ T(1:KI-1,1:KI-1) - WR ]*X = SCALE*WORK.
480 DO 60 J = KI - 1, 1, -1
487 IF( T( J, J-1 ).NE.ZERO ) THEN
495 * 1-by-1 diagonal block
497 CALL SLALN2( .FALSE., 1, 1, SMIN, ONE, T( J, J ),
498 $ LDT, ONE, ONE, WORK( J+IV*N ), N, WR,
499 $ ZERO, X, 2, SCALE, XNORM, IERR )
501 * Scale X(1,1) to avoid overflow when updating
502 * the right-hand side.
504 IF( XNORM.GT.ONE ) THEN
505 IF( WORK( J ).GT.BIGNUM / XNORM ) THEN
506 X( 1, 1 ) = X( 1, 1 ) / XNORM
507 SCALE = SCALE / XNORM
514 $ CALL SSCAL( KI, SCALE, WORK( 1+IV*N ), 1 )
515 WORK( J+IV*N ) = X( 1, 1 )
517 * Update right-hand side
519 CALL SAXPY( J-1, -X( 1, 1 ), T( 1, J ), 1,
520 $ WORK( 1+IV*N ), 1 )
524 * 2-by-2 diagonal block
526 CALL SLALN2( .FALSE., 2, 1, SMIN, ONE,
527 $ T( J-1, J-1 ), LDT, ONE, ONE,
528 $ WORK( J-1+IV*N ), N, WR, ZERO, X, 2,
529 $ SCALE, XNORM, IERR )
531 * Scale X(1,1) and X(2,1) to avoid overflow when
532 * updating the right-hand side.
534 IF( XNORM.GT.ONE ) THEN
535 BETA = MAX( WORK( J-1 ), WORK( J ) )
536 IF( BETA.GT.BIGNUM / XNORM ) THEN
537 X( 1, 1 ) = X( 1, 1 ) / XNORM
538 X( 2, 1 ) = X( 2, 1 ) / XNORM
539 SCALE = SCALE / XNORM
546 $ CALL SSCAL( KI, SCALE, WORK( 1+IV*N ), 1 )
547 WORK( J-1+IV*N ) = X( 1, 1 )
548 WORK( J +IV*N ) = X( 2, 1 )
550 * Update right-hand side
552 CALL SAXPY( J-2, -X( 1, 1 ), T( 1, J-1 ), 1,
553 $ WORK( 1+IV*N ), 1 )
554 CALL SAXPY( J-2, -X( 2, 1 ), T( 1, J ), 1,
555 $ WORK( 1+IV*N ), 1 )
559 * Copy the vector x or Q*x to VR and normalize.
562 * ------------------------------
563 * no back-transform: copy x to VR and normalize.
564 CALL SCOPY( KI, WORK( 1 + IV*N ), 1, VR( 1, IS ), 1 )
566 II = ISAMAX( KI, VR( 1, IS ), 1 )
567 REMAX = ONE / ABS( VR( II, IS ) )
568 CALL SSCAL( KI, REMAX, VR( 1, IS ), 1 )
574 ELSE IF( NB.EQ.1 ) THEN
575 * ------------------------------
576 * version 1: back-transform each vector with GEMV, Q*x.
578 $ CALL SGEMV( 'N', N, KI-1, ONE, VR, LDVR,
579 $ WORK( 1 + IV*N ), 1, WORK( KI + IV*N ),
582 II = ISAMAX( N, VR( 1, KI ), 1 )
583 REMAX = ONE / ABS( VR( II, KI ) )
584 CALL SSCAL( N, REMAX, VR( 1, KI ), 1 )
587 * ------------------------------
588 * version 2: back-transform block of vectors with GEMM
589 * zero out below vector
591 WORK( K + IV*N ) = ZERO
594 * back-transform and normalization is done below
598 * --------------------------------------------------------
599 * Complex right eigenvector.
602 * [ ( T(KI-1,KI-1) T(KI-1,KI) ) - (WR + I*WI) ]*X = 0.
603 * [ ( T(KI, KI-1) T(KI, KI) ) ]
605 IF( ABS( T( KI-1, KI ) ).GE.ABS( T( KI, KI-1 ) ) ) THEN
606 WORK( KI-1 + (IV-1)*N ) = ONE
607 WORK( KI + (IV )*N ) = WI / T( KI-1, KI )
609 WORK( KI-1 + (IV-1)*N ) = -WI / T( KI, KI-1 )
610 WORK( KI + (IV )*N ) = ONE
612 WORK( KI + (IV-1)*N ) = ZERO
613 WORK( KI-1 + (IV )*N ) = ZERO
615 * Form right-hand side.
618 WORK( K+(IV-1)*N ) = -WORK( KI-1+(IV-1)*N )*T(K,KI-1)
619 WORK( K+(IV )*N ) = -WORK( KI +(IV )*N )*T(K,KI )
622 * Solve upper quasi-triangular system:
623 * [ T(1:KI-2,1:KI-2) - (WR+i*WI) ]*X = SCALE*(WORK+i*WORK2)
626 DO 90 J = KI - 2, 1, -1
633 IF( T( J, J-1 ).NE.ZERO ) THEN
641 * 1-by-1 diagonal block
643 CALL SLALN2( .FALSE., 1, 2, SMIN, ONE, T( J, J ),
644 $ LDT, ONE, ONE, WORK( J+(IV-1)*N ), N,
645 $ WR, WI, X, 2, SCALE, XNORM, IERR )
647 * Scale X(1,1) and X(1,2) to avoid overflow when
648 * updating the right-hand side.
650 IF( XNORM.GT.ONE ) THEN
651 IF( WORK( J ).GT.BIGNUM / XNORM ) THEN
652 X( 1, 1 ) = X( 1, 1 ) / XNORM
653 X( 1, 2 ) = X( 1, 2 ) / XNORM
654 SCALE = SCALE / XNORM
660 IF( SCALE.NE.ONE ) THEN
661 CALL SSCAL( KI, SCALE, WORK( 1+(IV-1)*N ), 1 )
662 CALL SSCAL( KI, SCALE, WORK( 1+(IV )*N ), 1 )
664 WORK( J+(IV-1)*N ) = X( 1, 1 )
665 WORK( J+(IV )*N ) = X( 1, 2 )
667 * Update the right-hand side
669 CALL SAXPY( J-1, -X( 1, 1 ), T( 1, J ), 1,
670 $ WORK( 1+(IV-1)*N ), 1 )
671 CALL SAXPY( J-1, -X( 1, 2 ), T( 1, J ), 1,
672 $ WORK( 1+(IV )*N ), 1 )
676 * 2-by-2 diagonal block
678 CALL SLALN2( .FALSE., 2, 2, SMIN, ONE,
679 $ T( J-1, J-1 ), LDT, ONE, ONE,
680 $ WORK( J-1+(IV-1)*N ), N, WR, WI, X, 2,
681 $ SCALE, XNORM, IERR )
683 * Scale X to avoid overflow when updating
684 * the right-hand side.
686 IF( XNORM.GT.ONE ) THEN
687 BETA = MAX( WORK( J-1 ), WORK( J ) )
688 IF( BETA.GT.BIGNUM / XNORM ) THEN
690 X( 1, 1 ) = X( 1, 1 )*REC
691 X( 1, 2 ) = X( 1, 2 )*REC
692 X( 2, 1 ) = X( 2, 1 )*REC
693 X( 2, 2 ) = X( 2, 2 )*REC
700 IF( SCALE.NE.ONE ) THEN
701 CALL SSCAL( KI, SCALE, WORK( 1+(IV-1)*N ), 1 )
702 CALL SSCAL( KI, SCALE, WORK( 1+(IV )*N ), 1 )
704 WORK( J-1+(IV-1)*N ) = X( 1, 1 )
705 WORK( J +(IV-1)*N ) = X( 2, 1 )
706 WORK( J-1+(IV )*N ) = X( 1, 2 )
707 WORK( J +(IV )*N ) = X( 2, 2 )
709 * Update the right-hand side
711 CALL SAXPY( J-2, -X( 1, 1 ), T( 1, J-1 ), 1,
712 $ WORK( 1+(IV-1)*N ), 1 )
713 CALL SAXPY( J-2, -X( 2, 1 ), T( 1, J ), 1,
714 $ WORK( 1+(IV-1)*N ), 1 )
715 CALL SAXPY( J-2, -X( 1, 2 ), T( 1, J-1 ), 1,
716 $ WORK( 1+(IV )*N ), 1 )
717 CALL SAXPY( J-2, -X( 2, 2 ), T( 1, J ), 1,
718 $ WORK( 1+(IV )*N ), 1 )
722 * Copy the vector x or Q*x to VR and normalize.
725 * ------------------------------
726 * no back-transform: copy x to VR and normalize.
727 CALL SCOPY( KI, WORK( 1+(IV-1)*N ), 1, VR(1,IS-1), 1 )
728 CALL SCOPY( KI, WORK( 1+(IV )*N ), 1, VR(1,IS ), 1 )
732 EMAX = MAX( EMAX, ABS( VR( K, IS-1 ) )+
733 $ ABS( VR( K, IS ) ) )
736 CALL SSCAL( KI, REMAX, VR( 1, IS-1 ), 1 )
737 CALL SSCAL( KI, REMAX, VR( 1, IS ), 1 )
744 ELSE IF( NB.EQ.1 ) THEN
745 * ------------------------------
746 * version 1: back-transform each vector with GEMV, Q*x.
748 CALL SGEMV( 'N', N, KI-2, ONE, VR, LDVR,
749 $ WORK( 1 + (IV-1)*N ), 1,
750 $ WORK( KI-1 + (IV-1)*N ), VR(1,KI-1), 1)
751 CALL SGEMV( 'N', N, KI-2, ONE, VR, LDVR,
752 $ WORK( 1 + (IV)*N ), 1,
753 $ WORK( KI + (IV)*N ), VR( 1, KI ), 1 )
755 CALL SSCAL( N, WORK(KI-1+(IV-1)*N), VR(1,KI-1), 1)
756 CALL SSCAL( N, WORK(KI +(IV )*N), VR(1,KI ), 1)
761 EMAX = MAX( EMAX, ABS( VR( K, KI-1 ) )+
762 $ ABS( VR( K, KI ) ) )
765 CALL SSCAL( N, REMAX, VR( 1, KI-1 ), 1 )
766 CALL SSCAL( N, REMAX, VR( 1, KI ), 1 )
769 * ------------------------------
770 * version 2: back-transform block of vectors with GEMM
771 * zero out below vector
773 WORK( K + (IV-1)*N ) = ZERO
774 WORK( K + (IV )*N ) = ZERO
776 ISCOMPLEX( IV-1 ) = -IP
779 * back-transform and normalization is done below
784 * --------------------------------------------------------
785 * Blocked version of back-transform
786 * For complex case, KI2 includes both vectors (KI-1 and KI)
793 * Columns IV:NB of work are valid vectors.
794 * When the number of vectors stored reaches NB-1 or NB,
795 * or if this was last vector, do the GEMM
796 IF( (IV.LE.2) .OR. (KI2.EQ.1) ) THEN
797 CALL SGEMM( 'N', 'N', N, NB-IV+1, KI2+NB-IV, ONE,
799 $ WORK( 1 + (IV)*N ), N,
801 $ WORK( 1 + (NB+IV)*N ), N )
804 IF( ISCOMPLEX(K).EQ.0 ) THEN
806 II = ISAMAX( N, WORK( 1 + (NB+K)*N ), 1 )
807 REMAX = ONE / ABS( WORK( II + (NB+K)*N ) )
808 ELSE IF( ISCOMPLEX(K).EQ.1 ) THEN
809 * first eigenvector of conjugate pair
813 $ ABS( WORK( II + (NB+K )*N ) )+
814 $ ABS( WORK( II + (NB+K+1)*N ) ) )
817 * else if ISCOMPLEX(K).EQ.-1
818 * second eigenvector of conjugate pair
819 * reuse same REMAX as previous K
821 CALL SSCAL( N, REMAX, WORK( 1 + (NB+K)*N ), 1 )
823 CALL SLACPY( 'F', N, NB-IV+1,
824 $ WORK( 1 + (NB+IV)*N ), N,
825 $ VR( 1, KI2 ), LDVR )
830 END IF ! blocked back-transform
840 * ============================================================
841 * Compute left eigenvectors.
843 * IV is index of column in current block.
844 * For complex left vector, uses IV for real part and IV+1 for complex part.
845 * Non-blocked version always uses IV=1;
846 * blocked version starts with IV=1, goes up to NB-1 or NB.
847 * (Note the "0-th" column is used for 1-norms computed above.)
853 * previous iteration (ki-1) was first of conjugate pair,
854 * so this ki is second of conjugate pair; skip to end of loop
857 ELSE IF( KI.EQ.N ) THEN
858 * last column, so this ki must be real eigenvalue
860 ELSE IF( T( KI+1, KI ).EQ.ZERO ) THEN
861 * zero on sub-diagonal, so this ki is real eigenvalue
864 * non-zero on sub-diagonal, so this ki is first of conjugate pair
869 IF( .NOT.SELECT( KI ) )
873 * Compute the KI-th eigenvalue (WR,WI).
878 $ WI = SQRT( ABS( T( KI, KI+1 ) ) )*
879 $ SQRT( ABS( T( KI+1, KI ) ) )
880 SMIN = MAX( ULP*( ABS( WR )+ABS( WI ) ), SMLNUM )
884 * --------------------------------------------------------
885 * Real left eigenvector
887 WORK( KI + IV*N ) = ONE
889 * Form right-hand side.
892 WORK( K + IV*N ) = -T( KI, K )
895 * Solve transposed quasi-triangular system:
896 * [ T(KI+1:N,KI+1:N) - WR ]**T * X = SCALE*WORK
909 IF( T( J+1, J ).NE.ZERO ) THEN
917 * 1-by-1 diagonal block
919 * Scale if necessary to avoid overflow when forming
920 * the right-hand side.
922 IF( WORK( J ).GT.VCRIT ) THEN
924 CALL SSCAL( N-KI+1, REC, WORK( KI+IV*N ), 1 )
929 WORK( J+IV*N ) = WORK( J+IV*N ) -
930 $ SDOT( J-KI-1, T( KI+1, J ), 1,
931 $ WORK( KI+1+IV*N ), 1 )
933 * Solve [ T(J,J) - WR ]**T * X = WORK
935 CALL SLALN2( .FALSE., 1, 1, SMIN, ONE, T( J, J ),
936 $ LDT, ONE, ONE, WORK( J+IV*N ), N, WR,
937 $ ZERO, X, 2, SCALE, XNORM, IERR )
942 $ CALL SSCAL( N-KI+1, SCALE, WORK( KI+IV*N ), 1 )
943 WORK( J+IV*N ) = X( 1, 1 )
944 VMAX = MAX( ABS( WORK( J+IV*N ) ), VMAX )
945 VCRIT = BIGNUM / VMAX
949 * 2-by-2 diagonal block
951 * Scale if necessary to avoid overflow when forming
952 * the right-hand side.
954 BETA = MAX( WORK( J ), WORK( J+1 ) )
955 IF( BETA.GT.VCRIT ) THEN
957 CALL SSCAL( N-KI+1, REC, WORK( KI+IV*N ), 1 )
962 WORK( J+IV*N ) = WORK( J+IV*N ) -
963 $ SDOT( J-KI-1, T( KI+1, J ), 1,
964 $ WORK( KI+1+IV*N ), 1 )
966 WORK( J+1+IV*N ) = WORK( J+1+IV*N ) -
967 $ SDOT( J-KI-1, T( KI+1, J+1 ), 1,
968 $ WORK( KI+1+IV*N ), 1 )
971 * [ T(J,J)-WR T(J,J+1) ]**T * X = SCALE*( WORK1 )
972 * [ T(J+1,J) T(J+1,J+1)-WR ] ( WORK2 )
974 CALL SLALN2( .TRUE., 2, 1, SMIN, ONE, T( J, J ),
975 $ LDT, ONE, ONE, WORK( J+IV*N ), N, WR,
976 $ ZERO, X, 2, SCALE, XNORM, IERR )
981 $ CALL SSCAL( N-KI+1, SCALE, WORK( KI+IV*N ), 1 )
982 WORK( J +IV*N ) = X( 1, 1 )
983 WORK( J+1+IV*N ) = X( 2, 1 )
985 VMAX = MAX( ABS( WORK( J +IV*N ) ),
986 $ ABS( WORK( J+1+IV*N ) ), VMAX )
987 VCRIT = BIGNUM / VMAX
992 * Copy the vector x or Q*x to VL and normalize.
995 * ------------------------------
996 * no back-transform: copy x to VL and normalize.
997 CALL SCOPY( N-KI+1, WORK( KI + IV*N ), 1,
1000 II = ISAMAX( N-KI+1, VL( KI, IS ), 1 ) + KI - 1
1001 REMAX = ONE / ABS( VL( II, IS ) )
1002 CALL SSCAL( N-KI+1, REMAX, VL( KI, IS ), 1 )
1004 DO 180 K = 1, KI - 1
1008 ELSE IF( NB.EQ.1 ) THEN
1009 * ------------------------------
1010 * version 1: back-transform each vector with GEMV, Q*x.
1012 $ CALL SGEMV( 'N', N, N-KI, ONE,
1013 $ VL( 1, KI+1 ), LDVL,
1014 $ WORK( KI+1 + IV*N ), 1,
1015 $ WORK( KI + IV*N ), VL( 1, KI ), 1 )
1017 II = ISAMAX( N, VL( 1, KI ), 1 )
1018 REMAX = ONE / ABS( VL( II, KI ) )
1019 CALL SSCAL( N, REMAX, VL( 1, KI ), 1 )
1022 * ------------------------------
1023 * version 2: back-transform block of vectors with GEMM
1024 * zero out above vector
1025 * could go from KI-NV+1 to KI-1
1027 WORK( K + IV*N ) = ZERO
1029 ISCOMPLEX( IV ) = IP
1030 * back-transform and normalization is done below
1034 * --------------------------------------------------------
1035 * Complex left eigenvector.
1038 * [ ( T(KI,KI) T(KI,KI+1) )**T - (WR - I* WI) ]*X = 0.
1039 * [ ( T(KI+1,KI) T(KI+1,KI+1) ) ]
1041 IF( ABS( T( KI, KI+1 ) ).GE.ABS( T( KI+1, KI ) ) ) THEN
1042 WORK( KI + (IV )*N ) = WI / T( KI, KI+1 )
1043 WORK( KI+1 + (IV+1)*N ) = ONE
1045 WORK( KI + (IV )*N ) = ONE
1046 WORK( KI+1 + (IV+1)*N ) = -WI / T( KI+1, KI )
1048 WORK( KI+1 + (IV )*N ) = ZERO
1049 WORK( KI + (IV+1)*N ) = ZERO
1051 * Form right-hand side.
1053 DO 190 K = KI + 2, N
1054 WORK( K+(IV )*N ) = -WORK( KI +(IV )*N )*T(KI, K)
1055 WORK( K+(IV+1)*N ) = -WORK( KI+1+(IV+1)*N )*T(KI+1,K)
1058 * Solve transposed quasi-triangular system:
1059 * [ T(KI+2:N,KI+2:N)**T - (WR-i*WI) ]*X = WORK1+i*WORK2
1065 DO 200 J = KI + 2, N
1072 IF( T( J+1, J ).NE.ZERO ) THEN
1080 * 1-by-1 diagonal block
1082 * Scale if necessary to avoid overflow when
1083 * forming the right-hand side elements.
1085 IF( WORK( J ).GT.VCRIT ) THEN
1087 CALL SSCAL( N-KI+1, REC, WORK(KI+(IV )*N), 1 )
1088 CALL SSCAL( N-KI+1, REC, WORK(KI+(IV+1)*N), 1 )
1093 WORK( J+(IV )*N ) = WORK( J+(IV)*N ) -
1094 $ SDOT( J-KI-2, T( KI+2, J ), 1,
1095 $ WORK( KI+2+(IV)*N ), 1 )
1096 WORK( J+(IV+1)*N ) = WORK( J+(IV+1)*N ) -
1097 $ SDOT( J-KI-2, T( KI+2, J ), 1,
1098 $ WORK( KI+2+(IV+1)*N ), 1 )
1100 * Solve [ T(J,J)-(WR-i*WI) ]*(X11+i*X12)= WK+I*WK2
1102 CALL SLALN2( .FALSE., 1, 2, SMIN, ONE, T( J, J ),
1103 $ LDT, ONE, ONE, WORK( J+IV*N ), N, WR,
1104 $ -WI, X, 2, SCALE, XNORM, IERR )
1106 * Scale if necessary
1108 IF( SCALE.NE.ONE ) THEN
1109 CALL SSCAL( N-KI+1, SCALE, WORK(KI+(IV )*N), 1)
1110 CALL SSCAL( N-KI+1, SCALE, WORK(KI+(IV+1)*N), 1)
1112 WORK( J+(IV )*N ) = X( 1, 1 )
1113 WORK( J+(IV+1)*N ) = X( 1, 2 )
1114 VMAX = MAX( ABS( WORK( J+(IV )*N ) ),
1115 $ ABS( WORK( J+(IV+1)*N ) ), VMAX )
1116 VCRIT = BIGNUM / VMAX
1120 * 2-by-2 diagonal block
1122 * Scale if necessary to avoid overflow when forming
1123 * the right-hand side elements.
1125 BETA = MAX( WORK( J ), WORK( J+1 ) )
1126 IF( BETA.GT.VCRIT ) THEN
1128 CALL SSCAL( N-KI+1, REC, WORK(KI+(IV )*N), 1 )
1129 CALL SSCAL( N-KI+1, REC, WORK(KI+(IV+1)*N), 1 )
1134 WORK( J +(IV )*N ) = WORK( J+(IV)*N ) -
1135 $ SDOT( J-KI-2, T( KI+2, J ), 1,
1136 $ WORK( KI+2+(IV)*N ), 1 )
1138 WORK( J +(IV+1)*N ) = WORK( J+(IV+1)*N ) -
1139 $ SDOT( J-KI-2, T( KI+2, J ), 1,
1140 $ WORK( KI+2+(IV+1)*N ), 1 )
1142 WORK( J+1+(IV )*N ) = WORK( J+1+(IV)*N ) -
1143 $ SDOT( J-KI-2, T( KI+2, J+1 ), 1,
1144 $ WORK( KI+2+(IV)*N ), 1 )
1146 WORK( J+1+(IV+1)*N ) = WORK( J+1+(IV+1)*N ) -
1147 $ SDOT( J-KI-2, T( KI+2, J+1 ), 1,
1148 $ WORK( KI+2+(IV+1)*N ), 1 )
1150 * Solve 2-by-2 complex linear equation
1151 * [ (T(j,j) T(j,j+1) )**T - (wr-i*wi)*I ]*X = SCALE*B
1152 * [ (T(j+1,j) T(j+1,j+1)) ]
1154 CALL SLALN2( .TRUE., 2, 2, SMIN, ONE, T( J, J ),
1155 $ LDT, ONE, ONE, WORK( J+IV*N ), N, WR,
1156 $ -WI, X, 2, SCALE, XNORM, IERR )
1158 * Scale if necessary
1160 IF( SCALE.NE.ONE ) THEN
1161 CALL SSCAL( N-KI+1, SCALE, WORK(KI+(IV )*N), 1)
1162 CALL SSCAL( N-KI+1, SCALE, WORK(KI+(IV+1)*N), 1)
1164 WORK( J +(IV )*N ) = X( 1, 1 )
1165 WORK( J +(IV+1)*N ) = X( 1, 2 )
1166 WORK( J+1+(IV )*N ) = X( 2, 1 )
1167 WORK( J+1+(IV+1)*N ) = X( 2, 2 )
1168 VMAX = MAX( ABS( X( 1, 1 ) ), ABS( X( 1, 2 ) ),
1169 $ ABS( X( 2, 1 ) ), ABS( X( 2, 2 ) ),
1171 VCRIT = BIGNUM / VMAX
1176 * Copy the vector x or Q*x to VL and normalize.
1178 IF( .NOT.OVER ) THEN
1179 * ------------------------------
1180 * no back-transform: copy x to VL and normalize.
1181 CALL SCOPY( N-KI+1, WORK( KI + (IV )*N ), 1,
1183 CALL SCOPY( N-KI+1, WORK( KI + (IV+1)*N ), 1,
1184 $ VL( KI, IS+1 ), 1 )
1188 EMAX = MAX( EMAX, ABS( VL( K, IS ) )+
1189 $ ABS( VL( K, IS+1 ) ) )
1192 CALL SSCAL( N-KI+1, REMAX, VL( KI, IS ), 1 )
1193 CALL SSCAL( N-KI+1, REMAX, VL( KI, IS+1 ), 1 )
1195 DO 230 K = 1, KI - 1
1197 VL( K, IS+1 ) = ZERO
1200 ELSE IF( NB.EQ.1 ) THEN
1201 * ------------------------------
1202 * version 1: back-transform each vector with GEMV, Q*x.
1203 IF( KI.LT.N-1 ) THEN
1204 CALL SGEMV( 'N', N, N-KI-1, ONE,
1205 $ VL( 1, KI+2 ), LDVL,
1206 $ WORK( KI+2 + (IV)*N ), 1,
1207 $ WORK( KI + (IV)*N ),
1209 CALL SGEMV( 'N', N, N-KI-1, ONE,
1210 $ VL( 1, KI+2 ), LDVL,
1211 $ WORK( KI+2 + (IV+1)*N ), 1,
1212 $ WORK( KI+1 + (IV+1)*N ),
1213 $ VL( 1, KI+1 ), 1 )
1215 CALL SSCAL( N, WORK(KI+ (IV )*N), VL(1, KI ), 1)
1216 CALL SSCAL( N, WORK(KI+1+(IV+1)*N), VL(1, KI+1), 1)
1221 EMAX = MAX( EMAX, ABS( VL( K, KI ) )+
1222 $ ABS( VL( K, KI+1 ) ) )
1225 CALL SSCAL( N, REMAX, VL( 1, KI ), 1 )
1226 CALL SSCAL( N, REMAX, VL( 1, KI+1 ), 1 )
1229 * ------------------------------
1230 * version 2: back-transform block of vectors with GEMM
1231 * zero out above vector
1232 * could go from KI-NV+1 to KI-1
1234 WORK( K + (IV )*N ) = ZERO
1235 WORK( K + (IV+1)*N ) = ZERO
1237 ISCOMPLEX( IV ) = IP
1238 ISCOMPLEX( IV+1 ) = -IP
1240 * back-transform and normalization is done below
1245 * --------------------------------------------------------
1246 * Blocked version of back-transform
1247 * For complex case, KI2 includes both vectors (KI and KI+1)
1254 * Columns 1:IV of work are valid vectors.
1255 * When the number of vectors stored reaches NB-1 or NB,
1256 * or if this was last vector, do the GEMM
1257 IF( (IV.GE.NB-1) .OR. (KI2.EQ.N) ) THEN
1258 CALL SGEMM( 'N', 'N', N, IV, N-KI2+IV, ONE,
1259 $ VL( 1, KI2-IV+1 ), LDVL,
1260 $ WORK( KI2-IV+1 + (1)*N ), N,
1262 $ WORK( 1 + (NB+1)*N ), N )
1265 IF( ISCOMPLEX(K).EQ.0) THEN
1267 II = ISAMAX( N, WORK( 1 + (NB+K)*N ), 1 )
1268 REMAX = ONE / ABS( WORK( II + (NB+K)*N ) )
1269 ELSE IF( ISCOMPLEX(K).EQ.1) THEN
1270 * first eigenvector of conjugate pair
1274 $ ABS( WORK( II + (NB+K )*N ) )+
1275 $ ABS( WORK( II + (NB+K+1)*N ) ) )
1278 * else if ISCOMPLEX(K).EQ.-1
1279 * second eigenvector of conjugate pair
1280 * reuse same REMAX as previous K
1282 CALL SSCAL( N, REMAX, WORK( 1 + (NB+K)*N ), 1 )
1284 CALL SLACPY( 'F', N, IV,
1285 $ WORK( 1 + (NB+1)*N ), N,
1286 $ VL( 1, KI2-IV+1 ), LDVL )
1291 END IF ! blocked back-transform