LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
sgeesx.f
Go to the documentation of this file.
1*> \brief <b> SGEESX computes the eigenvalues, the Schur form, and, optionally, the matrix of Schur vectors for GE matrices</b>
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download SGEESX + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sgeesx.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sgeesx.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sgeesx.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE SGEESX( JOBVS, SORT, SELECT, SENSE, N, A, LDA, SDIM,
20* WR, WI, VS, LDVS, RCONDE, RCONDV, WORK, LWORK,
21* IWORK, LIWORK, BWORK, INFO )
22*
23* .. Scalar Arguments ..
24* CHARACTER JOBVS, SENSE, SORT
25* INTEGER INFO, LDA, LDVS, LIWORK, LWORK, N, SDIM
26* REAL RCONDE, RCONDV
27* ..
28* .. Array Arguments ..
29* LOGICAL BWORK( * )
30* INTEGER IWORK( * )
31* REAL A( LDA, * ), VS( LDVS, * ), WI( * ), WORK( * ),
32* $ WR( * )
33* ..
34* .. Function Arguments ..
35* LOGICAL SELECT
36* EXTERNAL SELECT
37* ..
38*
39*
40*> \par Purpose:
41* =============
42*>
43*> \verbatim
44*>
45*> SGEESX computes for an N-by-N real nonsymmetric matrix A, the
46*> eigenvalues, the real Schur form T, and, optionally, the matrix of
47*> Schur vectors Z. This gives the Schur factorization A = Z*T*(Z**T).
48*>
49*> Optionally, it also orders the eigenvalues on the diagonal of the
50*> real Schur form so that selected eigenvalues are at the top left;
51*> computes a reciprocal condition number for the average of the
52*> selected eigenvalues (RCONDE); and computes a reciprocal condition
53*> number for the right invariant subspace corresponding to the
54*> selected eigenvalues (RCONDV). The leading columns of Z form an
55*> orthonormal basis for this invariant subspace.
56*>
57*> For further explanation of the reciprocal condition numbers RCONDE
58*> and RCONDV, see Section 4.10 of the LAPACK Users' Guide (where
59*> these quantities are called s and sep respectively).
60*>
61*> A real matrix is in real Schur form if it is upper quasi-triangular
62*> with 1-by-1 and 2-by-2 blocks. 2-by-2 blocks will be standardized in
63*> the form
64*> [ a b ]
65*> [ c a ]
66*>
67*> where b*c < 0. The eigenvalues of such a block are a +- sqrt(bc).
68*> \endverbatim
69*
70* Arguments:
71* ==========
72*
73*> \param[in] JOBVS
74*> \verbatim
75*> JOBVS is CHARACTER*1
76*> = 'N': Schur vectors are not computed;
77*> = 'V': Schur vectors are computed.
78*> \endverbatim
79*>
80*> \param[in] SORT
81*> \verbatim
82*> SORT is CHARACTER*1
83*> Specifies whether or not to order the eigenvalues on the
84*> diagonal of the Schur form.
85*> = 'N': Eigenvalues are not ordered;
86*> = 'S': Eigenvalues are ordered (see SELECT).
87*> \endverbatim
88*>
89*> \param[in] SELECT
90*> \verbatim
91*> SELECT is a LOGICAL FUNCTION of two REAL arguments
92*> SELECT must be declared EXTERNAL in the calling subroutine.
93*> If SORT = 'S', SELECT is used to select eigenvalues to sort
94*> to the top left of the Schur form.
95*> If SORT = 'N', SELECT is not referenced.
96*> An eigenvalue WR(j)+sqrt(-1)*WI(j) is selected if
97*> SELECT(WR(j),WI(j)) is true; i.e., if either one of a
98*> complex conjugate pair of eigenvalues is selected, then both
99*> are. Note that a selected complex eigenvalue may no longer
100*> satisfy SELECT(WR(j),WI(j)) = .TRUE. after ordering, since
101*> ordering may change the value of complex eigenvalues
102*> (especially if the eigenvalue is ill-conditioned); in this
103*> case INFO may be set to N+3 (see INFO below).
104*> \endverbatim
105*>
106*> \param[in] SENSE
107*> \verbatim
108*> SENSE is CHARACTER*1
109*> Determines which reciprocal condition numbers are computed.
110*> = 'N': None are computed;
111*> = 'E': Computed for average of selected eigenvalues only;
112*> = 'V': Computed for selected right invariant subspace only;
113*> = 'B': Computed for both.
114*> If SENSE = 'E', 'V' or 'B', SORT must equal 'S'.
115*> \endverbatim
116*>
117*> \param[in] N
118*> \verbatim
119*> N is INTEGER
120*> The order of the matrix A. N >= 0.
121*> \endverbatim
122*>
123*> \param[in,out] A
124*> \verbatim
125*> A is REAL array, dimension (LDA, N)
126*> On entry, the N-by-N matrix A.
127*> On exit, A is overwritten by its real Schur form T.
128*> \endverbatim
129*>
130*> \param[in] LDA
131*> \verbatim
132*> LDA is INTEGER
133*> The leading dimension of the array A. LDA >= max(1,N).
134*> \endverbatim
135*>
136*> \param[out] SDIM
137*> \verbatim
138*> SDIM is INTEGER
139*> If SORT = 'N', SDIM = 0.
140*> If SORT = 'S', SDIM = number of eigenvalues (after sorting)
141*> for which SELECT is true. (Complex conjugate
142*> pairs for which SELECT is true for either
143*> eigenvalue count as 2.)
144*> \endverbatim
145*>
146*> \param[out] WR
147*> \verbatim
148*> WR is REAL array, dimension (N)
149*> \endverbatim
150*>
151*> \param[out] WI
152*> \verbatim
153*> WI is REAL array, dimension (N)
154*> WR and WI contain the real and imaginary parts, respectively,
155*> of the computed eigenvalues, in the same order that they
156*> appear on the diagonal of the output Schur form T. Complex
157*> conjugate pairs of eigenvalues appear consecutively with the
158*> eigenvalue having the positive imaginary part first.
159*> \endverbatim
160*>
161*> \param[out] VS
162*> \verbatim
163*> VS is REAL array, dimension (LDVS,N)
164*> If JOBVS = 'V', VS contains the orthogonal matrix Z of Schur
165*> vectors.
166*> If JOBVS = 'N', VS is not referenced.
167*> \endverbatim
168*>
169*> \param[in] LDVS
170*> \verbatim
171*> LDVS is INTEGER
172*> The leading dimension of the array VS. LDVS >= 1, and if
173*> JOBVS = 'V', LDVS >= N.
174*> \endverbatim
175*>
176*> \param[out] RCONDE
177*> \verbatim
178*> RCONDE is REAL
179*> If SENSE = 'E' or 'B', RCONDE contains the reciprocal
180*> condition number for the average of the selected eigenvalues.
181*> Not referenced if SENSE = 'N' or 'V'.
182*> \endverbatim
183*>
184*> \param[out] RCONDV
185*> \verbatim
186*> RCONDV is REAL
187*> If SENSE = 'V' or 'B', RCONDV contains the reciprocal
188*> condition number for the selected right invariant subspace.
189*> Not referenced if SENSE = 'N' or 'E'.
190*> \endverbatim
191*>
192*> \param[out] WORK
193*> \verbatim
194*> WORK is REAL array, dimension (MAX(1,LWORK))
195*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
196*> \endverbatim
197*>
198*> \param[in] LWORK
199*> \verbatim
200*> LWORK is INTEGER
201*> The dimension of the array WORK. LWORK >= max(1,3*N).
202*> Also, if SENSE = 'E' or 'V' or 'B',
203*> LWORK >= N+2*SDIM*(N-SDIM), where SDIM is the number of
204*> selected eigenvalues computed by this routine. Note that
205*> N+2*SDIM*(N-SDIM) <= N+N*N/2. Note also that an error is only
206*> returned if LWORK < max(1,3*N), but if SENSE = 'E' or 'V' or
207*> 'B' this may not be large enough.
208*> For good performance, LWORK must generally be larger.
209*>
210*> If LWORK = -1, then a workspace query is assumed; the routine
211*> only calculates upper bounds on the optimal sizes of the
212*> arrays WORK and IWORK, returns these values as the first
213*> entries of the WORK and IWORK arrays, and no error messages
214*> related to LWORK or LIWORK are issued by XERBLA.
215*> \endverbatim
216*>
217*> \param[out] IWORK
218*> \verbatim
219*> IWORK is INTEGER array, dimension (MAX(1,LIWORK))
220*> On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK.
221*> \endverbatim
222*>
223*> \param[in] LIWORK
224*> \verbatim
225*> LIWORK is INTEGER
226*> The dimension of the array IWORK.
227*> LIWORK >= 1; if SENSE = 'V' or 'B', LIWORK >= SDIM*(N-SDIM).
228*> Note that SDIM*(N-SDIM) <= N*N/4. Note also that an error is
229*> only returned if LIWORK < 1, but if SENSE = 'V' or 'B' this
230*> may not be large enough.
231*>
232*> If LIWORK = -1, then a workspace query is assumed; the
233*> routine only calculates upper bounds on the optimal sizes of
234*> the arrays WORK and IWORK, returns these values as the first
235*> entries of the WORK and IWORK arrays, and no error messages
236*> related to LWORK or LIWORK are issued by XERBLA.
237*> \endverbatim
238*>
239*> \param[out] BWORK
240*> \verbatim
241*> BWORK is LOGICAL array, dimension (N)
242*> Not referenced if SORT = 'N'.
243*> \endverbatim
244*>
245*> \param[out] INFO
246*> \verbatim
247*> INFO is INTEGER
248*> = 0: successful exit
249*> < 0: if INFO = -i, the i-th argument had an illegal value.
250*> > 0: if INFO = i, and i is
251*> <= N: the QR algorithm failed to compute all the
252*> eigenvalues; elements 1:ILO-1 and i+1:N of WR and WI
253*> contain those eigenvalues which have converged; if
254*> JOBVS = 'V', VS contains the transformation which
255*> reduces A to its partially converged Schur form.
256*> = N+1: the eigenvalues could not be reordered because some
257*> eigenvalues were too close to separate (the problem
258*> is very ill-conditioned);
259*> = N+2: after reordering, roundoff changed values of some
260*> complex eigenvalues so that leading eigenvalues in
261*> the Schur form no longer satisfy SELECT=.TRUE. This
262*> could also be caused by underflow due to scaling.
263*> \endverbatim
264*
265* Authors:
266* ========
267*
268*> \author Univ. of Tennessee
269*> \author Univ. of California Berkeley
270*> \author Univ. of Colorado Denver
271*> \author NAG Ltd.
272*
273*> \ingroup geesx
274*
275* =====================================================================
276 SUBROUTINE sgeesx( JOBVS, SORT, SELECT, SENSE, N, A, LDA, SDIM,
277 $ WR, WI, VS, LDVS, RCONDE, RCONDV, WORK, LWORK,
278 $ IWORK, LIWORK, BWORK, INFO )
279*
280* -- LAPACK driver routine --
281* -- LAPACK is a software package provided by Univ. of Tennessee, --
282* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
283*
284* .. Scalar Arguments ..
285 CHARACTER JOBVS, SENSE, SORT
286 INTEGER INFO, LDA, LDVS, LIWORK, LWORK, N, SDIM
287 REAL RCONDE, RCONDV
288* ..
289* .. Array Arguments ..
290 LOGICAL BWORK( * )
291 INTEGER IWORK( * )
292 REAL A( LDA, * ), VS( LDVS, * ), WI( * ), WORK( * ),
293 $ wr( * )
294* ..
295* .. Function Arguments ..
296 LOGICAL SELECT
297 EXTERNAL SELECT
298* ..
299*
300* =====================================================================
301*
302* .. Parameters ..
303 REAL ZERO, ONE
304 PARAMETER ( ZERO = 0.0e0, one = 1.0e0 )
305* ..
306* .. Local Scalars ..
307 LOGICAL CURSL, LASTSL, LQUERY, LST2SL, SCALEA, WANTSB,
308 $ WANTSE, WANTSN, WANTST, WANTSV, WANTVS
309 INTEGER HSWORK, I, I1, I2, IBAL, ICOND, IERR, IEVAL,
310 $ IHI, ILO, INXT, IP, ITAU, IWRK, LWRK, LIWRK,
311 $ maxwrk, minwrk
312 REAL ANRM, BIGNUM, CSCALE, EPS, SMLNUM
313* ..
314* .. Local Arrays ..
315 REAL DUM( 1 )
316* ..
317* .. External Subroutines ..
318 EXTERNAL scopy, sgebak, sgebal, sgehrd,
319 $ shseqr,
321* ..
322* .. External Functions ..
323 LOGICAL LSAME
324 INTEGER ILAENV
325 REAL SLAMCH, SLANGE, SROUNDUP_LWORK
326 EXTERNAL lsame, ilaenv, slamch,
327 $ slange, sroundup_lwork
328* ..
329* .. Intrinsic Functions ..
330 INTRINSIC max, sqrt
331* ..
332* .. Executable Statements ..
333*
334* Test the input arguments
335*
336 info = 0
337 wantvs = lsame( jobvs, 'V' )
338 wantst = lsame( sort, 'S' )
339 wantsn = lsame( sense, 'N' )
340 wantse = lsame( sense, 'E' )
341 wantsv = lsame( sense, 'V' )
342 wantsb = lsame( sense, 'B' )
343 lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
344*
345 IF( ( .NOT.wantvs ) .AND. ( .NOT.lsame( jobvs, 'N' ) ) ) THEN
346 info = -1
347 ELSE IF( ( .NOT.wantst ) .AND.
348 $ ( .NOT.lsame( sort, 'N' ) ) ) THEN
349 info = -2
350 ELSE IF( .NOT.( wantsn .OR. wantse .OR. wantsv .OR. wantsb ) .OR.
351 $ ( .NOT.wantst .AND. .NOT.wantsn ) ) THEN
352 info = -4
353 ELSE IF( n.LT.0 ) THEN
354 info = -5
355 ELSE IF( lda.LT.max( 1, n ) ) THEN
356 info = -7
357 ELSE IF( ldvs.LT.1 .OR. ( wantvs .AND. ldvs.LT.n ) ) THEN
358 info = -12
359 END IF
360*
361* Compute workspace
362* (Note: Comments in the code beginning "RWorkspace:" describe the
363* minimal amount of real workspace needed at that point in the
364* code, as well as the preferred amount for good performance.
365* IWorkspace refers to integer workspace.
366* NB refers to the optimal block size for the immediately
367* following subroutine, as returned by ILAENV.
368* HSWORK refers to the workspace preferred by SHSEQR, as
369* calculated below. HSWORK is computed assuming ILO=1 and IHI=N,
370* the worst case.
371* If SENSE = 'E', 'V' or 'B', then the amount of workspace needed
372* depends on SDIM, which is computed by the routine STRSEN later
373* in the code.)
374*
375 IF( info.EQ.0 ) THEN
376 liwrk = 1
377 IF( n.EQ.0 ) THEN
378 minwrk = 1
379 lwrk = 1
380 ELSE
381 maxwrk = 2*n + n*ilaenv( 1, 'SGEHRD', ' ', n, 1, n, 0 )
382 minwrk = 3*n
383*
384 CALL shseqr( 'S', jobvs, n, 1, n, a, lda, wr, wi, vs,
385 $ ldvs,
386 $ work, -1, ieval )
387 hswork = int( work( 1 ) )
388*
389 IF( .NOT.wantvs ) THEN
390 maxwrk = max( maxwrk, n + hswork )
391 ELSE
392 maxwrk = max( maxwrk, 2*n + ( n - 1 )*ilaenv( 1,
393 $ 'SORGHR', ' ', n, 1, n, -1 ) )
394 maxwrk = max( maxwrk, n + hswork )
395 END IF
396 lwrk = maxwrk
397 IF( .NOT.wantsn )
398 $ lwrk = max( lwrk, n + ( n*n )/2 )
399 IF( wantsv .OR. wantsb )
400 $ liwrk = ( n*n )/4
401 END IF
402 iwork( 1 ) = liwrk
403 work( 1 ) = sroundup_lwork(lwrk)
404*
405 IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
406 info = -16
407 ELSE IF( liwork.LT.1 .AND. .NOT.lquery ) THEN
408 info = -18
409 END IF
410 END IF
411*
412 IF( info.NE.0 ) THEN
413 CALL xerbla( 'SGEESX', -info )
414 RETURN
415 ELSE IF( lquery ) THEN
416 RETURN
417 END IF
418*
419* Quick return if possible
420*
421 IF( n.EQ.0 ) THEN
422 sdim = 0
423 RETURN
424 END IF
425*
426* Get machine constants
427*
428 eps = slamch( 'P' )
429 smlnum = slamch( 'S' )
430 bignum = one / smlnum
431 smlnum = sqrt( smlnum ) / eps
432 bignum = one / smlnum
433*
434* Scale A if max element outside range [SMLNUM,BIGNUM]
435*
436 anrm = slange( 'M', n, n, a, lda, dum )
437 scalea = .false.
438 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
439 scalea = .true.
440 cscale = smlnum
441 ELSE IF( anrm.GT.bignum ) THEN
442 scalea = .true.
443 cscale = bignum
444 END IF
445 IF( scalea )
446 $ CALL slascl( 'G', 0, 0, anrm, cscale, n, n, a, lda, ierr )
447*
448* Permute the matrix to make it more nearly triangular
449* (RWorkspace: need N)
450*
451 ibal = 1
452 CALL sgebal( 'P', n, a, lda, ilo, ihi, work( ibal ), ierr )
453*
454* Reduce to upper Hessenberg form
455* (RWorkspace: need 3*N, prefer 2*N+N*NB)
456*
457 itau = n + ibal
458 iwrk = n + itau
459 CALL sgehrd( n, ilo, ihi, a, lda, work( itau ), work( iwrk ),
460 $ lwork-iwrk+1, ierr )
461*
462 IF( wantvs ) THEN
463*
464* Copy Householder vectors to VS
465*
466 CALL slacpy( 'L', n, n, a, lda, vs, ldvs )
467*
468* Generate orthogonal matrix in VS
469* (RWorkspace: need 3*N-1, prefer 2*N+(N-1)*NB)
470*
471 CALL sorghr( n, ilo, ihi, vs, ldvs, work( itau ),
472 $ work( iwrk ),
473 $ lwork-iwrk+1, ierr )
474 END IF
475*
476 sdim = 0
477*
478* Perform QR iteration, accumulating Schur vectors in VS if desired
479* (RWorkspace: need N+1, prefer N+HSWORK (see comments) )
480*
481 iwrk = itau
482 CALL shseqr( 'S', jobvs, n, ilo, ihi, a, lda, wr, wi, vs, ldvs,
483 $ work( iwrk ), lwork-iwrk+1, ieval )
484 IF( ieval.GT.0 )
485 $ info = ieval
486*
487* Sort eigenvalues if desired
488*
489 IF( wantst .AND. info.EQ.0 ) THEN
490 IF( scalea ) THEN
491 CALL slascl( 'G', 0, 0, cscale, anrm, n, 1, wr, n, ierr )
492 CALL slascl( 'G', 0, 0, cscale, anrm, n, 1, wi, n, ierr )
493 END IF
494 DO 10 i = 1, n
495 bwork( i ) = SELECT( wr( i ), wi( i ) )
496 10 CONTINUE
497*
498* Reorder eigenvalues, transform Schur vectors, and compute
499* reciprocal condition numbers
500* (RWorkspace: if SENSE is not 'N', need N+2*SDIM*(N-SDIM)
501* otherwise, need N )
502* (IWorkspace: if SENSE is 'V' or 'B', need SDIM*(N-SDIM)
503* otherwise, need 0 )
504*
505 CALL strsen( sense, jobvs, bwork, n, a, lda, vs, ldvs, wr,
506 $ wi,
507 $ sdim, rconde, rcondv, work( iwrk ), lwork-iwrk+1,
508 $ iwork, liwork, icond )
509 IF( .NOT.wantsn )
510 $ maxwrk = max( maxwrk, n+2*sdim*( n-sdim ) )
511 IF( icond.EQ.-15 ) THEN
512*
513* Not enough real workspace
514*
515 info = -16
516 ELSE IF( icond.EQ.-17 ) THEN
517*
518* Not enough integer workspace
519*
520 info = -18
521 ELSE IF( icond.GT.0 ) THEN
522*
523* STRSEN failed to reorder or to restore standard Schur form
524*
525 info = icond + n
526 END IF
527 END IF
528*
529 IF( wantvs ) THEN
530*
531* Undo balancing
532* (RWorkspace: need N)
533*
534 CALL sgebak( 'P', 'R', n, ilo, ihi, work( ibal ), n, vs,
535 $ ldvs,
536 $ ierr )
537 END IF
538*
539 IF( scalea ) THEN
540*
541* Undo scaling for the Schur form of A
542*
543 CALL slascl( 'H', 0, 0, cscale, anrm, n, n, a, lda, ierr )
544 CALL scopy( n, a, lda+1, wr, 1 )
545 IF( ( wantsv .OR. wantsb ) .AND. info.EQ.0 ) THEN
546 dum( 1 ) = rcondv
547 CALL slascl( 'G', 0, 0, cscale, anrm, 1, 1, dum, 1,
548 $ ierr )
549 rcondv = dum( 1 )
550 END IF
551 IF( cscale.EQ.smlnum ) THEN
552*
553* If scaling back towards underflow, adjust WI if an
554* offdiagonal element of a 2-by-2 block in the Schur form
555* underflows.
556*
557 IF( ieval.GT.0 ) THEN
558 i1 = ieval + 1
559 i2 = ihi - 1
560 CALL slascl( 'G', 0, 0, cscale, anrm, ilo-1, 1, wi, n,
561 $ ierr )
562 ELSE IF( wantst ) THEN
563 i1 = 1
564 i2 = n - 1
565 ELSE
566 i1 = ilo
567 i2 = ihi - 1
568 END IF
569 inxt = i1 - 1
570 DO 20 i = i1, i2
571 IF( i.LT.inxt )
572 $ GO TO 20
573 IF( wi( i ).EQ.zero ) THEN
574 inxt = i + 1
575 ELSE
576 IF( a( i+1, i ).EQ.zero ) THEN
577 wi( i ) = zero
578 wi( i+1 ) = zero
579 ELSE IF( a( i+1, i ).NE.zero .AND. a( i, i+1 ).EQ.
580 $ zero ) THEN
581 wi( i ) = zero
582 wi( i+1 ) = zero
583 IF( i.GT.1 )
584 $ CALL sswap( i-1, a( 1, i ), 1, a( 1, i+1 ),
585 $ 1 )
586 IF( n.GT.i+1 )
587 $ CALL sswap( n-i-1, a( i, i+2 ), lda,
588 $ a( i+1, i+2 ), lda )
589 IF( wantvs ) THEN
590 CALL sswap( n, vs( 1, i ), 1, vs( 1, i+1 ),
591 $ 1 )
592 END IF
593 a( i, i+1 ) = a( i+1, i )
594 a( i+1, i ) = zero
595 END IF
596 inxt = i + 2
597 END IF
598 20 CONTINUE
599 END IF
600 CALL slascl( 'G', 0, 0, cscale, anrm, n-ieval, 1,
601 $ wi( ieval+1 ), max( n-ieval, 1 ), ierr )
602 END IF
603*
604 IF( wantst .AND. info.EQ.0 ) THEN
605*
606* Check if reordering successful
607*
608 lastsl = .true.
609 lst2sl = .true.
610 sdim = 0
611 ip = 0
612 DO 30 i = 1, n
613 cursl = SELECT( wr( i ), wi( i ) )
614 IF( wi( i ).EQ.zero ) THEN
615 IF( cursl )
616 $ sdim = sdim + 1
617 ip = 0
618 IF( cursl .AND. .NOT.lastsl )
619 $ info = n + 2
620 ELSE
621 IF( ip.EQ.1 ) THEN
622*
623* Last eigenvalue of conjugate pair
624*
625 cursl = cursl .OR. lastsl
626 lastsl = cursl
627 IF( cursl )
628 $ sdim = sdim + 2
629 ip = -1
630 IF( cursl .AND. .NOT.lst2sl )
631 $ info = n + 2
632 ELSE
633*
634* First eigenvalue of conjugate pair
635*
636 ip = 1
637 END IF
638 END IF
639 lst2sl = lastsl
640 lastsl = cursl
641 30 CONTINUE
642 END IF
643*
644 work( 1 ) = sroundup_lwork(maxwrk)
645 IF( wantsv .OR. wantsb ) THEN
646 iwork( 1 ) = sdim*(n-sdim)
647 ELSE
648 iwork( 1 ) = 1
649 END IF
650*
651 RETURN
652*
653* End of SGEESX
654*
655 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
subroutine sgebak(job, side, n, ilo, ihi, scale, m, v, ldv, info)
SGEBAK
Definition sgebak.f:128
subroutine sgebal(job, n, a, lda, ilo, ihi, scale, info)
SGEBAL
Definition sgebal.f:161
subroutine sgeesx(jobvs, sort, select, sense, n, a, lda, sdim, wr, wi, vs, ldvs, rconde, rcondv, work, lwork, iwork, liwork, bwork, info)
SGEESX computes the eigenvalues, the Schur form, and, optionally, the matrix of Schur vectors for GE ...
Definition sgeesx.f:279
subroutine sgehrd(n, ilo, ihi, a, lda, tau, work, lwork, info)
SGEHRD
Definition sgehrd.f:166
subroutine shseqr(job, compz, n, ilo, ihi, h, ldh, wr, wi, z, ldz, work, lwork, info)
SHSEQR
Definition shseqr.f:314
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:101
subroutine slascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
SLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition slascl.f:142
subroutine sswap(n, sx, incx, sy, incy)
SSWAP
Definition sswap.f:82
subroutine strsen(job, compq, select, n, t, ldt, q, ldq, wr, wi, m, s, sep, work, lwork, iwork, liwork, info)
STRSEN
Definition strsen.f:313
subroutine sorghr(n, ilo, ihi, a, lda, tau, work, lwork, info)
SORGHR
Definition sorghr.f:125