LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ sgees()

subroutine sgees ( character jobvs,
character sort,
external select,
integer n,
real, dimension( lda, * ) a,
integer lda,
integer sdim,
real, dimension( * ) wr,
real, dimension( * ) wi,
real, dimension( ldvs, * ) vs,
integer ldvs,
real, dimension( * ) work,
integer lwork,
logical, dimension( * ) bwork,
integer info )

SGEES computes the eigenvalues, the Schur form, and, optionally, the matrix of Schur vectors for GE matrices

Download SGEES + dependencies [TGZ] [ZIP] [TXT]

Purpose:
!>
!> SGEES computes for an N-by-N real nonsymmetric matrix A, the
!> eigenvalues, the real Schur form T, and, optionally, the matrix of
!> Schur vectors Z.  This gives the Schur factorization A = Z*T*(Z**T).
!>
!> Optionally, it also orders the eigenvalues on the diagonal of the
!> real Schur form so that selected eigenvalues are at the top left.
!> The leading columns of Z then form an orthonormal basis for the
!> invariant subspace corresponding to the selected eigenvalues.
!>
!> A matrix is in real Schur form if it is upper quasi-triangular with
!> 1-by-1 and 2-by-2 blocks. 2-by-2 blocks will be standardized in the
!> form
!>         [  a  b  ]
!>         [  c  a  ]
!>
!> where b*c < 0. The eigenvalues of such a block are a +- sqrt(bc).
!> 
Parameters
[in]JOBVS
!>          JOBVS is CHARACTER*1
!>          = 'N': Schur vectors are not computed;
!>          = 'V': Schur vectors are computed.
!> 
[in]SORT
!>          SORT is CHARACTER*1
!>          Specifies whether or not to order the eigenvalues on the
!>          diagonal of the Schur form.
!>          = 'N': Eigenvalues are not ordered;
!>          = 'S': Eigenvalues are ordered (see SELECT).
!> 
[in]SELECT
!>          SELECT is a LOGICAL FUNCTION of two REAL arguments
!>          SELECT must be declared EXTERNAL in the calling subroutine.
!>          If SORT = 'S', SELECT is used to select eigenvalues to sort
!>          to the top left of the Schur form.
!>          If SORT = 'N', SELECT is not referenced.
!>          An eigenvalue WR(j)+sqrt(-1)*WI(j) is selected if
!>          SELECT(WR(j),WI(j)) is true; i.e., if either one of a complex
!>          conjugate pair of eigenvalues is selected, then both complex
!>          eigenvalues are selected.
!>          Note that a selected complex eigenvalue may no longer
!>          satisfy SELECT(WR(j),WI(j)) = .TRUE. after ordering, since
!>          ordering may change the value of complex eigenvalues
!>          (especially if the eigenvalue is ill-conditioned); in this
!>          case INFO is set to N+2 (see INFO below).
!> 
[in]N
!>          N is INTEGER
!>          The order of the matrix A. N >= 0.
!> 
[in,out]A
!>          A is REAL array, dimension (LDA,N)
!>          On entry, the N-by-N matrix A.
!>          On exit, A has been overwritten by its real Schur form T.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[out]SDIM
!>          SDIM is INTEGER
!>          If SORT = 'N', SDIM = 0.
!>          If SORT = 'S', SDIM = number of eigenvalues (after sorting)
!>                         for which SELECT is true. (Complex conjugate
!>                         pairs for which SELECT is true for either
!>                         eigenvalue count as 2.)
!> 
[out]WR
!>          WR is REAL array, dimension (N)
!> 
[out]WI
!>          WI is REAL array, dimension (N)
!>          WR and WI contain the real and imaginary parts,
!>          respectively, of the computed eigenvalues in the same order
!>          that they appear on the diagonal of the output Schur form T.
!>          Complex conjugate pairs of eigenvalues will appear
!>          consecutively with the eigenvalue having the positive
!>          imaginary part first.
!> 
[out]VS
!>          VS is REAL array, dimension (LDVS,N)
!>          If JOBVS = 'V', VS contains the orthogonal matrix Z of Schur
!>          vectors.
!>          If JOBVS = 'N', VS is not referenced.
!> 
[in]LDVS
!>          LDVS is INTEGER
!>          The leading dimension of the array VS.  LDVS >= 1; if
!>          JOBVS = 'V', LDVS >= N.
!> 
[out]WORK
!>          WORK is REAL array, dimension (MAX(1,LWORK))
!>          On exit, if INFO = 0, WORK(1) contains the optimal LWORK.
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The dimension of the array WORK.  LWORK >= max(1,3*N).
!>          For good performance, LWORK must generally be larger.
!>
!>          If LWORK = -1, then a workspace query is assumed; the routine
!>          only calculates the optimal size of the WORK array, returns
!>          this value as the first entry of the WORK array, and no error
!>          message related to LWORK is issued by XERBLA.
!> 
[out]BWORK
!>          BWORK is LOGICAL array, dimension (N)
!>          Not referenced if SORT = 'N'.
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0: successful exit
!>          < 0: if INFO = -i, the i-th argument had an illegal value.
!>          > 0: if INFO = i, and i is
!>             <= N: the QR algorithm failed to compute all the
!>                   eigenvalues; elements 1:ILO-1 and i+1:N of WR and WI
!>                   contain those eigenvalues which have converged; if
!>                   JOBVS = 'V', VS contains the matrix which reduces A
!>                   to its partially converged Schur form.
!>             = N+1: the eigenvalues could not be reordered because some
!>                   eigenvalues were too close to separate (the problem
!>                   is very ill-conditioned);
!>             = N+2: after reordering, roundoff changed values of some
!>                   complex eigenvalues so that leading eigenvalues in
!>                   the Schur form no longer satisfy SELECT=.TRUE.  This
!>                   could also be caused by underflow due to scaling.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 212 of file sgees.f.

214*
215* -- LAPACK driver routine --
216* -- LAPACK is a software package provided by Univ. of Tennessee, --
217* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
218*
219* .. Scalar Arguments ..
220 CHARACTER JOBVS, SORT
221 INTEGER INFO, LDA, LDVS, LWORK, N, SDIM
222* ..
223* .. Array Arguments ..
224 LOGICAL BWORK( * )
225 REAL A( LDA, * ), VS( LDVS, * ), WI( * ), WORK( * ),
226 $ WR( * )
227* ..
228* .. Function Arguments ..
229 LOGICAL SELECT
230 EXTERNAL SELECT
231* ..
232*
233* =====================================================================
234*
235* .. Parameters ..
236 REAL ZERO, ONE
237 parameter( zero = 0.0e0, one = 1.0e0 )
238* ..
239* .. Local Scalars ..
240 LOGICAL CURSL, LASTSL, LQUERY, LST2SL, SCALEA, WANTST,
241 $ WANTVS
242 INTEGER HSWORK, I, I1, I2, IBAL, ICOND, IERR, IEVAL,
243 $ IHI, ILO, INXT, IP, ITAU, IWRK, MAXWRK, MINWRK
244 REAL ANRM, BIGNUM, CSCALE, EPS, S, SEP, SMLNUM
245* ..
246* .. Local Arrays ..
247 INTEGER IDUM( 1 )
248 REAL DUM( 1 )
249* ..
250* .. External Subroutines ..
251 EXTERNAL scopy, sgebak, sgebal, sgehrd, shseqr,
252 $ slacpy,
254* ..
255* .. External Functions ..
256 LOGICAL LSAME
257 INTEGER ILAENV
258 REAL SLAMCH, SLANGE, SROUNDUP_LWORK
259 EXTERNAL lsame, ilaenv, slamch, slange,
261* ..
262* .. Intrinsic Functions ..
263 INTRINSIC max, sqrt
264* ..
265* .. Executable Statements ..
266*
267* Test the input arguments
268*
269 info = 0
270 lquery = ( lwork.EQ.-1 )
271 wantvs = lsame( jobvs, 'V' )
272 wantst = lsame( sort, 'S' )
273 IF( ( .NOT.wantvs ) .AND. ( .NOT.lsame( jobvs, 'N' ) ) ) THEN
274 info = -1
275 ELSE IF( ( .NOT.wantst ) .AND.
276 $ ( .NOT.lsame( sort, 'N' ) ) ) THEN
277 info = -2
278 ELSE IF( n.LT.0 ) THEN
279 info = -4
280 ELSE IF( lda.LT.max( 1, n ) ) THEN
281 info = -6
282 ELSE IF( ldvs.LT.1 .OR. ( wantvs .AND. ldvs.LT.n ) ) THEN
283 info = -11
284 END IF
285*
286* Compute workspace
287* (Note: Comments in the code beginning "Workspace:" describe the
288* minimal amount of workspace needed at that point in the code,
289* as well as the preferred amount for good performance.
290* NB refers to the optimal block size for the immediately
291* following subroutine, as returned by ILAENV.
292* HSWORK refers to the workspace preferred by SHSEQR, as
293* calculated below. HSWORK is computed assuming ILO=1 and IHI=N,
294* the worst case.)
295*
296 IF( info.EQ.0 ) THEN
297 IF( n.EQ.0 ) THEN
298 minwrk = 1
299 maxwrk = 1
300 ELSE
301 maxwrk = 2*n + n*ilaenv( 1, 'SGEHRD', ' ', n, 1, n, 0 )
302 minwrk = 3*n
303*
304 CALL shseqr( 'S', jobvs, n, 1, n, a, lda, wr, wi, vs,
305 $ ldvs,
306 $ work, -1, ieval )
307 hswork = int( work( 1 ) )
308*
309 IF( .NOT.wantvs ) THEN
310 maxwrk = max( maxwrk, n + hswork )
311 ELSE
312 maxwrk = max( maxwrk, 2*n + ( n - 1 )*ilaenv( 1,
313 $ 'SORGHR', ' ', n, 1, n, -1 ) )
314 maxwrk = max( maxwrk, n + hswork )
315 END IF
316 END IF
317 work( 1 ) = sroundup_lwork(maxwrk)
318*
319 IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
320 info = -13
321 END IF
322 END IF
323*
324 IF( info.NE.0 ) THEN
325 CALL xerbla( 'SGEES ', -info )
326 RETURN
327 ELSE IF( lquery ) THEN
328 RETURN
329 END IF
330*
331* Quick return if possible
332*
333 IF( n.EQ.0 ) THEN
334 sdim = 0
335 RETURN
336 END IF
337*
338* Get machine constants
339*
340 eps = slamch( 'P' )
341 smlnum = slamch( 'S' )
342 bignum = one / smlnum
343 smlnum = sqrt( smlnum ) / eps
344 bignum = one / smlnum
345*
346* Scale A if max element outside range [SMLNUM,BIGNUM]
347*
348 anrm = slange( 'M', n, n, a, lda, dum )
349 scalea = .false.
350 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
351 scalea = .true.
352 cscale = smlnum
353 ELSE IF( anrm.GT.bignum ) THEN
354 scalea = .true.
355 cscale = bignum
356 END IF
357 IF( scalea )
358 $ CALL slascl( 'G', 0, 0, anrm, cscale, n, n, a, lda, ierr )
359*
360* Permute the matrix to make it more nearly triangular
361* (Workspace: need N)
362*
363 ibal = 1
364 CALL sgebal( 'P', n, a, lda, ilo, ihi, work( ibal ), ierr )
365*
366* Reduce to upper Hessenberg form
367* (Workspace: need 3*N, prefer 2*N+N*NB)
368*
369 itau = n + ibal
370 iwrk = n + itau
371 CALL sgehrd( n, ilo, ihi, a, lda, work( itau ), work( iwrk ),
372 $ lwork-iwrk+1, ierr )
373*
374 IF( wantvs ) THEN
375*
376* Copy Householder vectors to VS
377*
378 CALL slacpy( 'L', n, n, a, lda, vs, ldvs )
379*
380* Generate orthogonal matrix in VS
381* (Workspace: need 3*N-1, prefer 2*N+(N-1)*NB)
382*
383 CALL sorghr( n, ilo, ihi, vs, ldvs, work( itau ),
384 $ work( iwrk ),
385 $ lwork-iwrk+1, ierr )
386 END IF
387*
388 sdim = 0
389*
390* Perform QR iteration, accumulating Schur vectors in VS if desired
391* (Workspace: need N+1, prefer N+HSWORK (see comments) )
392*
393 iwrk = itau
394 CALL shseqr( 'S', jobvs, n, ilo, ihi, a, lda, wr, wi, vs, ldvs,
395 $ work( iwrk ), lwork-iwrk+1, ieval )
396 IF( ieval.GT.0 )
397 $ info = ieval
398*
399* Sort eigenvalues if desired
400*
401 IF( wantst .AND. info.EQ.0 ) THEN
402 IF( scalea ) THEN
403 CALL slascl( 'G', 0, 0, cscale, anrm, n, 1, wr, n, ierr )
404 CALL slascl( 'G', 0, 0, cscale, anrm, n, 1, wi, n, ierr )
405 END IF
406 DO 10 i = 1, n
407 bwork( i ) = SELECT( wr( i ), wi( i ) )
408 10 CONTINUE
409*
410* Reorder eigenvalues and transform Schur vectors
411* (Workspace: none needed)
412*
413 CALL strsen( 'N', jobvs, bwork, n, a, lda, vs, ldvs, wr, wi,
414 $ sdim, s, sep, work( iwrk ), lwork-iwrk+1, idum, 1,
415 $ icond )
416 IF( icond.GT.0 )
417 $ info = n + icond
418 END IF
419*
420 IF( wantvs ) THEN
421*
422* Undo balancing
423* (Workspace: need N)
424*
425 CALL sgebak( 'P', 'R', n, ilo, ihi, work( ibal ), n, vs,
426 $ ldvs,
427 $ ierr )
428 END IF
429*
430 IF( scalea ) THEN
431*
432* Undo scaling for the Schur form of A
433*
434 CALL slascl( 'H', 0, 0, cscale, anrm, n, n, a, lda, ierr )
435 CALL scopy( n, a, lda+1, wr, 1 )
436 IF( cscale.EQ.smlnum ) THEN
437*
438* If scaling back towards underflow, adjust WI if an
439* offdiagonal element of a 2-by-2 block in the Schur form
440* underflows.
441*
442 IF( ieval.GT.0 ) THEN
443 i1 = ieval + 1
444 i2 = ihi - 1
445 CALL slascl( 'G', 0, 0, cscale, anrm, ilo-1, 1, wi,
446 $ max( ilo-1, 1 ), ierr )
447 ELSE IF( wantst ) THEN
448 i1 = 1
449 i2 = n - 1
450 ELSE
451 i1 = ilo
452 i2 = ihi - 1
453 END IF
454 inxt = i1 - 1
455 DO 20 i = i1, i2
456 IF( i.LT.inxt )
457 $ GO TO 20
458 IF( wi( i ).EQ.zero ) THEN
459 inxt = i + 1
460 ELSE
461 IF( a( i+1, i ).EQ.zero ) THEN
462 wi( i ) = zero
463 wi( i+1 ) = zero
464 ELSE IF( a( i+1, i ).NE.zero .AND. a( i, i+1 ).EQ.
465 $ zero ) THEN
466 wi( i ) = zero
467 wi( i+1 ) = zero
468 IF( i.GT.1 )
469 $ CALL sswap( i-1, a( 1, i ), 1, a( 1, i+1 ),
470 $ 1 )
471 IF( n.GT.i+1 )
472 $ CALL sswap( n-i-1, a( i, i+2 ), lda,
473 $ a( i+1, i+2 ), lda )
474 IF( wantvs ) THEN
475 CALL sswap( n, vs( 1, i ), 1, vs( 1, i+1 ),
476 $ 1 )
477 END IF
478 a( i, i+1 ) = a( i+1, i )
479 a( i+1, i ) = zero
480 END IF
481 inxt = i + 2
482 END IF
483 20 CONTINUE
484 END IF
485*
486* Undo scaling for the imaginary part of the eigenvalues
487*
488 CALL slascl( 'G', 0, 0, cscale, anrm, n-ieval, 1,
489 $ wi( ieval+1 ), max( n-ieval, 1 ), ierr )
490 END IF
491*
492 IF( wantst .AND. info.EQ.0 ) THEN
493*
494* Check if reordering successful
495*
496 lastsl = .true.
497 lst2sl = .true.
498 sdim = 0
499 ip = 0
500 DO 30 i = 1, n
501 cursl = SELECT( wr( i ), wi( i ) )
502 IF( wi( i ).EQ.zero ) THEN
503 IF( cursl )
504 $ sdim = sdim + 1
505 ip = 0
506 IF( cursl .AND. .NOT.lastsl )
507 $ info = n + 2
508 ELSE
509 IF( ip.EQ.1 ) THEN
510*
511* Last eigenvalue of conjugate pair
512*
513 cursl = cursl .OR. lastsl
514 lastsl = cursl
515 IF( cursl )
516 $ sdim = sdim + 2
517 ip = -1
518 IF( cursl .AND. .NOT.lst2sl )
519 $ info = n + 2
520 ELSE
521*
522* First eigenvalue of conjugate pair
523*
524 ip = 1
525 END IF
526 END IF
527 lst2sl = lastsl
528 lastsl = cursl
529 30 CONTINUE
530 END IF
531*
532 work( 1 ) = sroundup_lwork(maxwrk)
533 RETURN
534*
535* End of SGEES
536*
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 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
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
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
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
real function slange(norm, m, n, a, lda, work)
SLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition slange.f:112
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
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
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
Here is the call graph for this function:
Here is the caller graph for this function: