LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dgeesx.f
Go to the documentation of this file.
1*> \brief <b> DGEESX 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 DGEESX + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgeesx.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgeesx.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgeesx.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE DGEESX( 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* DOUBLE PRECISION RCONDE, RCONDV
27* ..
28* .. Array Arguments ..
29* LOGICAL BWORK( * )
30* INTEGER IWORK( * )
31* DOUBLE PRECISION 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*> DGEESX 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION array, dimension (N)
149*> \endverbatim
150*>
151*> \param[out] WI
152*> \verbatim
153*> WI is DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION
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 DOUBLE PRECISION
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 DOUBLE PRECISION 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 dgeesx( 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 DOUBLE PRECISION RCONDE, RCONDV
288* ..
289* .. Array Arguments ..
290 LOGICAL BWORK( * )
291 INTEGER IWORK( * )
292 DOUBLE PRECISION 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 DOUBLE PRECISION ZERO, ONE
304 PARAMETER ( ZERO = 0.0d0, one = 1.0d0 )
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, LIWRK, LWRK,
311 $ maxwrk, minwrk
312 DOUBLE PRECISION ANRM, BIGNUM, CSCALE, EPS, SMLNUM
313* ..
314* .. Local Arrays ..
315 DOUBLE PRECISION DUM( 1 )
316* ..
317* .. External Subroutines ..
318 EXTERNAL dcopy, dgebak, dgebal, dgehrd, dhseqr,
319 $ dlacpy,
321* ..
322* .. External Functions ..
323 LOGICAL LSAME
324 INTEGER ILAENV
325 DOUBLE PRECISION DLAMCH, DLANGE
326 EXTERNAL lsame, ilaenv, dlamch, dlange
327* ..
328* .. Intrinsic Functions ..
329 INTRINSIC max, sqrt
330* ..
331* .. Executable Statements ..
332*
333* Test the input arguments
334*
335 info = 0
336 wantvs = lsame( jobvs, 'V' )
337 wantst = lsame( sort, 'S' )
338 wantsn = lsame( sense, 'N' )
339 wantse = lsame( sense, 'E' )
340 wantsv = lsame( sense, 'V' )
341 wantsb = lsame( sense, 'B' )
342 lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
343*
344 IF( ( .NOT.wantvs ) .AND. ( .NOT.lsame( jobvs, 'N' ) ) ) THEN
345 info = -1
346 ELSE IF( ( .NOT.wantst ) .AND.
347 $ ( .NOT.lsame( sort, 'N' ) ) ) THEN
348 info = -2
349 ELSE IF( .NOT.( wantsn .OR. wantse .OR. wantsv .OR. wantsb ) .OR.
350 $ ( .NOT.wantst .AND. .NOT.wantsn ) ) THEN
351 info = -4
352 ELSE IF( n.LT.0 ) THEN
353 info = -5
354 ELSE IF( lda.LT.max( 1, n ) ) THEN
355 info = -7
356 ELSE IF( ldvs.LT.1 .OR. ( wantvs .AND. ldvs.LT.n ) ) THEN
357 info = -12
358 END IF
359*
360* Compute workspace
361* (Note: Comments in the code beginning "RWorkspace:" describe the
362* minimal amount of real workspace needed at that point in the
363* code, as well as the preferred amount for good performance.
364* IWorkspace refers to integer workspace.
365* NB refers to the optimal block size for the immediately
366* following subroutine, as returned by ILAENV.
367* HSWORK refers to the workspace preferred by DHSEQR, as
368* calculated below. HSWORK is computed assuming ILO=1 and IHI=N,
369* the worst case.
370* If SENSE = 'E', 'V' or 'B', then the amount of workspace needed
371* depends on SDIM, which is computed by the routine DTRSEN later
372* in the code.)
373*
374 IF( info.EQ.0 ) THEN
375 liwrk = 1
376 IF( n.EQ.0 ) THEN
377 minwrk = 1
378 lwrk = 1
379 ELSE
380 maxwrk = 2*n + n*ilaenv( 1, 'DGEHRD', ' ', n, 1, n, 0 )
381 minwrk = 3*n
382*
383 CALL dhseqr( 'S', jobvs, n, 1, n, a, lda, wr, wi, vs,
384 $ ldvs,
385 $ work, -1, ieval )
386 hswork = int( work( 1 ) )
387*
388 IF( .NOT.wantvs ) THEN
389 maxwrk = max( maxwrk, n + hswork )
390 ELSE
391 maxwrk = max( maxwrk, 2*n + ( n - 1 )*ilaenv( 1,
392 $ 'DORGHR', ' ', n, 1, n, -1 ) )
393 maxwrk = max( maxwrk, n + hswork )
394 END IF
395 lwrk = maxwrk
396 IF( .NOT.wantsn )
397 $ lwrk = max( lwrk, n + ( n*n )/2 )
398 IF( wantsv .OR. wantsb )
399 $ liwrk = ( n*n )/4
400 END IF
401 iwork( 1 ) = liwrk
402 work( 1 ) = lwrk
403*
404 IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
405 info = -16
406 ELSE IF( liwork.LT.1 .AND. .NOT.lquery ) THEN
407 info = -18
408 END IF
409 END IF
410*
411 IF( info.NE.0 ) THEN
412 CALL xerbla( 'DGEESX', -info )
413 RETURN
414 ELSE IF( lquery ) THEN
415 RETURN
416 END IF
417*
418* Quick return if possible
419*
420 IF( n.EQ.0 ) THEN
421 sdim = 0
422 RETURN
423 END IF
424*
425* Get machine constants
426*
427 eps = dlamch( 'P' )
428 smlnum = dlamch( 'S' )
429 bignum = one / smlnum
430 smlnum = sqrt( smlnum ) / eps
431 bignum = one / smlnum
432*
433* Scale A if max element outside range [SMLNUM,BIGNUM]
434*
435 anrm = dlange( 'M', n, n, a, lda, dum )
436 scalea = .false.
437 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
438 scalea = .true.
439 cscale = smlnum
440 ELSE IF( anrm.GT.bignum ) THEN
441 scalea = .true.
442 cscale = bignum
443 END IF
444 IF( scalea )
445 $ CALL dlascl( 'G', 0, 0, anrm, cscale, n, n, a, lda, ierr )
446*
447* Permute the matrix to make it more nearly triangular
448* (RWorkspace: need N)
449*
450 ibal = 1
451 CALL dgebal( 'P', n, a, lda, ilo, ihi, work( ibal ), ierr )
452*
453* Reduce to upper Hessenberg form
454* (RWorkspace: need 3*N, prefer 2*N+N*NB)
455*
456 itau = n + ibal
457 iwrk = n + itau
458 CALL dgehrd( n, ilo, ihi, a, lda, work( itau ), work( iwrk ),
459 $ lwork-iwrk+1, ierr )
460*
461 IF( wantvs ) THEN
462*
463* Copy Householder vectors to VS
464*
465 CALL dlacpy( 'L', n, n, a, lda, vs, ldvs )
466*
467* Generate orthogonal matrix in VS
468* (RWorkspace: need 3*N-1, prefer 2*N+(N-1)*NB)
469*
470 CALL dorghr( n, ilo, ihi, vs, ldvs, work( itau ),
471 $ work( iwrk ),
472 $ lwork-iwrk+1, ierr )
473 END IF
474*
475 sdim = 0
476*
477* Perform QR iteration, accumulating Schur vectors in VS if desired
478* (RWorkspace: need N+1, prefer N+HSWORK (see comments) )
479*
480 iwrk = itau
481 CALL dhseqr( 'S', jobvs, n, ilo, ihi, a, lda, wr, wi, vs, ldvs,
482 $ work( iwrk ), lwork-iwrk+1, ieval )
483 IF( ieval.GT.0 )
484 $ info = ieval
485*
486* Sort eigenvalues if desired
487*
488 IF( wantst .AND. info.EQ.0 ) THEN
489 IF( scalea ) THEN
490 CALL dlascl( 'G', 0, 0, cscale, anrm, n, 1, wr, n, ierr )
491 CALL dlascl( 'G', 0, 0, cscale, anrm, n, 1, wi, n, ierr )
492 END IF
493 DO 10 i = 1, n
494 bwork( i ) = SELECT( wr( i ), wi( i ) )
495 10 CONTINUE
496*
497* Reorder eigenvalues, transform Schur vectors, and compute
498* reciprocal condition numbers
499* (RWorkspace: if SENSE is not 'N', need N+2*SDIM*(N-SDIM)
500* otherwise, need N )
501* (IWorkspace: if SENSE is 'V' or 'B', need SDIM*(N-SDIM)
502* otherwise, need 0 )
503*
504 CALL dtrsen( sense, jobvs, bwork, n, a, lda, vs, ldvs, wr,
505 $ wi,
506 $ sdim, rconde, rcondv, work( iwrk ), lwork-iwrk+1,
507 $ iwork, liwork, icond )
508 IF( .NOT.wantsn )
509 $ maxwrk = max( maxwrk, n+2*sdim*( n-sdim ) )
510 IF( icond.EQ.-15 ) THEN
511*
512* Not enough real workspace
513*
514 info = -16
515 ELSE IF( icond.EQ.-17 ) THEN
516*
517* Not enough integer workspace
518*
519 info = -18
520 ELSE IF( icond.GT.0 ) THEN
521*
522* DTRSEN failed to reorder or to restore standard Schur form
523*
524 info = icond + n
525 END IF
526 END IF
527*
528 IF( wantvs ) THEN
529*
530* Undo balancing
531* (RWorkspace: need N)
532*
533 CALL dgebak( 'P', 'R', n, ilo, ihi, work( ibal ), n, vs,
534 $ ldvs,
535 $ ierr )
536 END IF
537*
538 IF( scalea ) THEN
539*
540* Undo scaling for the Schur form of A
541*
542 CALL dlascl( 'H', 0, 0, cscale, anrm, n, n, a, lda, ierr )
543 CALL dcopy( n, a, lda+1, wr, 1 )
544 IF( ( wantsv .OR. wantsb ) .AND. info.EQ.0 ) THEN
545 dum( 1 ) = rcondv
546 CALL dlascl( 'G', 0, 0, cscale, anrm, 1, 1, dum, 1,
547 $ ierr )
548 rcondv = dum( 1 )
549 END IF
550 IF( cscale.EQ.smlnum ) THEN
551*
552* If scaling back towards underflow, adjust WI if an
553* offdiagonal element of a 2-by-2 block in the Schur form
554* underflows.
555*
556 IF( ieval.GT.0 ) THEN
557 i1 = ieval + 1
558 i2 = ihi - 1
559 CALL dlascl( 'G', 0, 0, cscale, anrm, ilo-1, 1, wi, n,
560 $ ierr )
561 ELSE IF( wantst ) THEN
562 i1 = 1
563 i2 = n - 1
564 ELSE
565 i1 = ilo
566 i2 = ihi - 1
567 END IF
568 inxt = i1 - 1
569 DO 20 i = i1, i2
570 IF( i.LT.inxt )
571 $ GO TO 20
572 IF( wi( i ).EQ.zero ) THEN
573 inxt = i + 1
574 ELSE
575 IF( a( i+1, i ).EQ.zero ) THEN
576 wi( i ) = zero
577 wi( i+1 ) = zero
578 ELSE IF( a( i+1, i ).NE.zero .AND. a( i, i+1 ).EQ.
579 $ zero ) THEN
580 wi( i ) = zero
581 wi( i+1 ) = zero
582 IF( i.GT.1 )
583 $ CALL dswap( i-1, a( 1, i ), 1, a( 1, i+1 ),
584 $ 1 )
585 IF( n.GT.i+1 )
586 $ CALL dswap( n-i-1, a( i, i+2 ), lda,
587 $ a( i+1, i+2 ), lda )
588 IF( wantvs ) THEN
589 CALL dswap( n, vs( 1, i ), 1, vs( 1, i+1 ),
590 $ 1 )
591 END IF
592 a( i, i+1 ) = a( i+1, i )
593 a( i+1, i ) = zero
594 END IF
595 inxt = i + 2
596 END IF
597 20 CONTINUE
598 END IF
599 CALL dlascl( 'G', 0, 0, cscale, anrm, n-ieval, 1,
600 $ wi( ieval+1 ), max( n-ieval, 1 ), ierr )
601 END IF
602*
603 IF( wantst .AND. info.EQ.0 ) THEN
604*
605* Check if reordering successful
606*
607 lastsl = .true.
608 lst2sl = .true.
609 sdim = 0
610 ip = 0
611 DO 30 i = 1, n
612 cursl = SELECT( wr( i ), wi( i ) )
613 IF( wi( i ).EQ.zero ) THEN
614 IF( cursl )
615 $ sdim = sdim + 1
616 ip = 0
617 IF( cursl .AND. .NOT.lastsl )
618 $ info = n + 2
619 ELSE
620 IF( ip.EQ.1 ) THEN
621*
622* Last eigenvalue of conjugate pair
623*
624 cursl = cursl .OR. lastsl
625 lastsl = cursl
626 IF( cursl )
627 $ sdim = sdim + 2
628 ip = -1
629 IF( cursl .AND. .NOT.lst2sl )
630 $ info = n + 2
631 ELSE
632*
633* First eigenvalue of conjugate pair
634*
635 ip = 1
636 END IF
637 END IF
638 lst2sl = lastsl
639 lastsl = cursl
640 30 CONTINUE
641 END IF
642*
643 work( 1 ) = maxwrk
644 IF( wantsv .OR. wantsb ) THEN
645 iwork( 1 ) = max( 1, sdim*( n-sdim ) )
646 ELSE
647 iwork( 1 ) = 1
648 END IF
649*
650 RETURN
651*
652* End of DGEESX
653*
654 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dgebak(job, side, n, ilo, ihi, scale, m, v, ldv, info)
DGEBAK
Definition dgebak.f:128
subroutine dgebal(job, n, a, lda, ilo, ihi, scale, info)
DGEBAL
Definition dgebal.f:161
subroutine dgeesx(jobvs, sort, select, sense, n, a, lda, sdim, wr, wi, vs, ldvs, rconde, rcondv, work, lwork, iwork, liwork, bwork, info)
DGEESX computes the eigenvalues, the Schur form, and, optionally, the matrix of Schur vectors for GE ...
Definition dgeesx.f:279
subroutine dgehrd(n, ilo, ihi, a, lda, tau, work, lwork, info)
DGEHRD
Definition dgehrd.f:166
subroutine dhseqr(job, compz, n, ilo, ihi, h, ldh, wr, wi, z, ldz, work, lwork, info)
DHSEQR
Definition dhseqr.f:314
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:101
subroutine dlascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
DLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition dlascl.f:142
subroutine dswap(n, dx, incx, dy, incy)
DSWAP
Definition dswap.f:82
subroutine dtrsen(job, compq, select, n, t, ldt, q, ldq, wr, wi, m, s, sep, work, lwork, iwork, liwork, info)
DTRSEN
Definition dtrsen.f:312
subroutine dorghr(n, ilo, ihi, a, lda, tau, work, lwork, info)
DORGHR
Definition dorghr.f:125