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

◆ cgeesx()

subroutine cgeesx ( character jobvs,
character sort,
external select,
character sense,
integer n,
complex, dimension( lda, * ) a,
integer lda,
integer sdim,
complex, dimension( * ) w,
complex, dimension( ldvs, * ) vs,
integer ldvs,
real rconde,
real rcondv,
complex, dimension( * ) work,
integer lwork,
real, dimension( * ) rwork,
logical, dimension( * ) bwork,
integer info )

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

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

Purpose:
!>
!> CGEESX computes for an N-by-N complex nonsymmetric matrix A, the
!> eigenvalues, the Schur form T, and, optionally, the matrix of Schur
!> vectors Z.  This gives the Schur factorization A = Z*T*(Z**H).
!>
!> Optionally, it also orders the eigenvalues on the diagonal of the
!> Schur form so that selected eigenvalues are at the top left;
!> computes a reciprocal condition number for the average of the
!> selected eigenvalues (RCONDE); and computes a reciprocal condition
!> number for the right invariant subspace corresponding to the
!> selected eigenvalues (RCONDV).  The leading columns of Z form an
!> orthonormal basis for this invariant subspace.
!>
!> For further explanation of the reciprocal condition numbers RCONDE
!> and RCONDV, see Section 4.10 of the LAPACK Users' Guide (where
!> these quantities are called s and sep respectively).
!>
!> A complex matrix is in Schur form if it is upper triangular.
!> 
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 one COMPLEX argument
!>          SELECT must be declared EXTERNAL in the calling subroutine.
!>          If SORT = 'S', SELECT is used to select eigenvalues to order
!>          to the top left of the Schur form.
!>          If SORT = 'N', SELECT is not referenced.
!>          An eigenvalue W(j) is selected if SELECT(W(j)) is true.
!> 
[in]SENSE
!>          SENSE is CHARACTER*1
!>          Determines which reciprocal condition numbers are computed.
!>          = 'N': None are computed;
!>          = 'E': Computed for average of selected eigenvalues only;
!>          = 'V': Computed for selected right invariant subspace only;
!>          = 'B': Computed for both.
!>          If SENSE = 'E', 'V' or 'B', SORT must equal 'S'.
!> 
[in]N
!>          N is INTEGER
!>          The order of the matrix A. N >= 0.
!> 
[in,out]A
!>          A is COMPLEX array, dimension (LDA, N)
!>          On entry, the N-by-N matrix A.
!>          On exit, A is overwritten by its 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 for which
!>                         SELECT is true.
!> 
[out]W
!>          W is COMPLEX array, dimension (N)
!>          W contains the computed eigenvalues, in the same order
!>          that they appear on the diagonal of the output Schur form T.
!> 
[out]VS
!>          VS is COMPLEX array, dimension (LDVS,N)
!>          If JOBVS = 'V', VS contains the unitary 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, and if
!>          JOBVS = 'V', LDVS >= N.
!> 
[out]RCONDE
!>          RCONDE is REAL
!>          If SENSE = 'E' or 'B', RCONDE contains the reciprocal
!>          condition number for the average of the selected eigenvalues.
!>          Not referenced if SENSE = 'N' or 'V'.
!> 
[out]RCONDV
!>          RCONDV is REAL
!>          If SENSE = 'V' or 'B', RCONDV contains the reciprocal
!>          condition number for the selected right invariant subspace.
!>          Not referenced if SENSE = 'N' or 'E'.
!> 
[out]WORK
!>          WORK is COMPLEX array, dimension (MAX(1,LWORK))
!>          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The dimension of the array WORK.  LWORK >= max(1,2*N).
!>          Also, if SENSE = 'E' or 'V' or 'B', LWORK >= 2*SDIM*(N-SDIM),
!>          where SDIM is the number of selected eigenvalues computed by
!>          this routine.  Note that 2*SDIM*(N-SDIM) <= N*N/2. Note also
!>          that an error is only returned if LWORK < max(1,2*N), but if
!>          SENSE = 'E' or 'V' or 'B' this may not be large enough.
!>          For good performance, LWORK must generally be larger.
!>
!>          If LWORK = -1, then a workspace query is assumed; the routine
!>          only calculates upper bound on the optimal size of the
!>          array WORK, returns this value as the first entry of the WORK
!>          array, and no error message related to LWORK is issued by
!>          XERBLA.
!> 
[out]RWORK
!>          RWORK is REAL array, dimension (N)
!> 
[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 W
!>                   contain those eigenvalues which have converged; if
!>                   JOBVS = 'V', VS contains the transformation 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 234 of file cgeesx.f.

238*
239* -- LAPACK driver routine --
240* -- LAPACK is a software package provided by Univ. of Tennessee, --
241* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
242*
243* .. Scalar Arguments ..
244 CHARACTER JOBVS, SENSE, SORT
245 INTEGER INFO, LDA, LDVS, LWORK, N, SDIM
246 REAL RCONDE, RCONDV
247* ..
248* .. Array Arguments ..
249 LOGICAL BWORK( * )
250 REAL RWORK( * )
251 COMPLEX A( LDA, * ), VS( LDVS, * ), W( * ), WORK( * )
252* ..
253* .. Function Arguments ..
254 LOGICAL SELECT
255 EXTERNAL SELECT
256* ..
257*
258* =====================================================================
259*
260* .. Parameters ..
261 REAL ZERO, ONE
262 parameter( zero = 0.0e0, one = 1.0e0 )
263* ..
264* .. Local Scalars ..
265 LOGICAL LQUERY, SCALEA, WANTSB, WANTSE, WANTSN, WANTST,
266 $ WANTSV, WANTVS
267 INTEGER HSWORK, I, IBAL, ICOND, IERR, IEVAL, IHI, ILO,
268 $ ITAU, IWRK, LWRK, MAXWRK, MINWRK
269 REAL ANRM, BIGNUM, CSCALE, EPS, SMLNUM
270* ..
271* .. Local Arrays ..
272 REAL DUM( 1 )
273* ..
274* .. External Subroutines ..
275 EXTERNAL ccopy, cgebak, cgebal, cgehrd, chseqr,
276 $ clacpy,
278* ..
279* .. External Functions ..
280 LOGICAL LSAME
281 INTEGER ILAENV
282 REAL CLANGE, SLAMCH, SROUNDUP_LWORK
283 EXTERNAL lsame, ilaenv, clange, slamch,
285* ..
286* .. Intrinsic Functions ..
287 INTRINSIC max, sqrt
288* ..
289* .. Executable Statements ..
290*
291* Test the input arguments
292*
293 info = 0
294 wantvs = lsame( jobvs, 'V' )
295 wantst = lsame( sort, 'S' )
296 wantsn = lsame( sense, 'N' )
297 wantse = lsame( sense, 'E' )
298 wantsv = lsame( sense, 'V' )
299 wantsb = lsame( sense, 'B' )
300 lquery = ( lwork.EQ.-1 )
301*
302 IF( ( .NOT.wantvs ) .AND. ( .NOT.lsame( jobvs, 'N' ) ) ) THEN
303 info = -1
304 ELSE IF( ( .NOT.wantst ) .AND.
305 $ ( .NOT.lsame( sort, 'N' ) ) ) THEN
306 info = -2
307 ELSE IF( .NOT.( wantsn .OR. wantse .OR. wantsv .OR. wantsb ) .OR.
308 $ ( .NOT.wantst .AND. .NOT.wantsn ) ) THEN
309 info = -4
310 ELSE IF( n.LT.0 ) THEN
311 info = -5
312 ELSE IF( lda.LT.max( 1, n ) ) THEN
313 info = -7
314 ELSE IF( ldvs.LT.1 .OR. ( wantvs .AND. ldvs.LT.n ) ) THEN
315 info = -11
316 END IF
317*
318* Compute workspace
319* (Note: Comments in the code beginning "Workspace:" describe the
320* minimal amount of real workspace needed at that point in the
321* code, as well as the preferred amount for good performance.
322* CWorkspace refers to complex workspace, and RWorkspace to real
323* workspace. NB refers to the optimal block size for the
324* immediately following subroutine, as returned by ILAENV.
325* HSWORK refers to the workspace preferred by CHSEQR, as
326* calculated below. HSWORK is computed assuming ILO=1 and IHI=N,
327* the worst case.
328* If SENSE = 'E', 'V' or 'B', then the amount of workspace needed
329* depends on SDIM, which is computed by the routine CTRSEN later
330* in the code.)
331*
332 IF( info.EQ.0 ) THEN
333 IF( n.EQ.0 ) THEN
334 minwrk = 1
335 lwrk = 1
336 ELSE
337 maxwrk = n + n*ilaenv( 1, 'CGEHRD', ' ', n, 1, n, 0 )
338 minwrk = 2*n
339*
340 CALL chseqr( 'S', jobvs, n, 1, n, a, lda, w, vs, ldvs,
341 $ work, -1, ieval )
342 hswork = int( work( 1 ) )
343*
344 IF( .NOT.wantvs ) THEN
345 maxwrk = max( maxwrk, hswork )
346 ELSE
347 maxwrk = max( maxwrk, n + ( n - 1 )*ilaenv( 1,
348 $ 'CUNGHR',
349 $ ' ', n, 1, n, -1 ) )
350 maxwrk = max( maxwrk, hswork )
351 END IF
352 lwrk = maxwrk
353 IF( .NOT.wantsn )
354 $ lwrk = max( lwrk, ( n*n )/2 )
355 END IF
356 work( 1 ) = sroundup_lwork(lwrk)
357*
358 IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
359 info = -15
360 END IF
361 END IF
362*
363 IF( info.NE.0 ) THEN
364 CALL xerbla( 'CGEESX', -info )
365 RETURN
366 ELSE IF( lquery ) THEN
367 RETURN
368 END IF
369*
370* Quick return if possible
371*
372 IF( n.EQ.0 ) THEN
373 sdim = 0
374 RETURN
375 END IF
376*
377* Get machine constants
378*
379 eps = slamch( 'P' )
380 smlnum = slamch( 'S' )
381 bignum = one / smlnum
382 smlnum = sqrt( smlnum ) / eps
383 bignum = one / smlnum
384*
385* Scale A if max element outside range [SMLNUM,BIGNUM]
386*
387 anrm = clange( 'M', n, n, a, lda, dum )
388 scalea = .false.
389 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
390 scalea = .true.
391 cscale = smlnum
392 ELSE IF( anrm.GT.bignum ) THEN
393 scalea = .true.
394 cscale = bignum
395 END IF
396 IF( scalea )
397 $ CALL clascl( 'G', 0, 0, anrm, cscale, n, n, a, lda, ierr )
398*
399*
400* Permute the matrix to make it more nearly triangular
401* (CWorkspace: none)
402* (RWorkspace: need N)
403*
404 ibal = 1
405 CALL cgebal( 'P', n, a, lda, ilo, ihi, rwork( ibal ), ierr )
406*
407* Reduce to upper Hessenberg form
408* (CWorkspace: need 2*N, prefer N+N*NB)
409* (RWorkspace: none)
410*
411 itau = 1
412 iwrk = n + itau
413 CALL cgehrd( n, ilo, ihi, a, lda, work( itau ), work( iwrk ),
414 $ lwork-iwrk+1, ierr )
415*
416 IF( wantvs ) THEN
417*
418* Copy Householder vectors to VS
419*
420 CALL clacpy( 'L', n, n, a, lda, vs, ldvs )
421*
422* Generate unitary matrix in VS
423* (CWorkspace: need 2*N-1, prefer N+(N-1)*NB)
424* (RWorkspace: none)
425*
426 CALL cunghr( n, ilo, ihi, vs, ldvs, work( itau ),
427 $ work( iwrk ),
428 $ lwork-iwrk+1, ierr )
429 END IF
430*
431 sdim = 0
432*
433* Perform QR iteration, accumulating Schur vectors in VS if desired
434* (CWorkspace: need 1, prefer HSWORK (see comments) )
435* (RWorkspace: none)
436*
437 iwrk = itau
438 CALL chseqr( 'S', jobvs, n, ilo, ihi, a, lda, w, vs, ldvs,
439 $ work( iwrk ), lwork-iwrk+1, ieval )
440 IF( ieval.GT.0 )
441 $ info = ieval
442*
443* Sort eigenvalues if desired
444*
445 IF( wantst .AND. info.EQ.0 ) THEN
446 IF( scalea )
447 $ CALL clascl( 'G', 0, 0, cscale, anrm, n, 1, w, n, ierr )
448 DO 10 i = 1, n
449 bwork( i ) = SELECT( w( i ) )
450 10 CONTINUE
451*
452* Reorder eigenvalues, transform Schur vectors, and compute
453* reciprocal condition numbers
454* (CWorkspace: if SENSE is not 'N', need 2*SDIM*(N-SDIM)
455* otherwise, need none )
456* (RWorkspace: none)
457*
458 CALL ctrsen( sense, jobvs, bwork, n, a, lda, vs, ldvs, w,
459 $ sdim,
460 $ rconde, rcondv, work( iwrk ), lwork-iwrk+1,
461 $ icond )
462 IF( .NOT.wantsn )
463 $ maxwrk = max( maxwrk, 2*sdim*( n-sdim ) )
464 IF( icond.EQ.-14 ) THEN
465*
466* Not enough complex workspace
467*
468 info = -15
469 END IF
470 END IF
471*
472 IF( wantvs ) THEN
473*
474* Undo balancing
475* (CWorkspace: none)
476* (RWorkspace: need N)
477*
478 CALL cgebak( 'P', 'R', n, ilo, ihi, rwork( ibal ), n, vs,
479 $ ldvs,
480 $ ierr )
481 END IF
482*
483 IF( scalea ) THEN
484*
485* Undo scaling for the Schur form of A
486*
487 CALL clascl( 'U', 0, 0, cscale, anrm, n, n, a, lda, ierr )
488 CALL ccopy( n, a, lda+1, w, 1 )
489 IF( ( wantsv .OR. wantsb ) .AND. info.EQ.0 ) THEN
490 dum( 1 ) = rcondv
491 CALL slascl( 'G', 0, 0, cscale, anrm, 1, 1, dum, 1,
492 $ ierr )
493 rcondv = dum( 1 )
494 END IF
495 END IF
496*
497 work( 1 ) = sroundup_lwork(maxwrk)
498 RETURN
499*
500* End of CGEESX
501*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ccopy(n, cx, incx, cy, incy)
CCOPY
Definition ccopy.f:81
subroutine cgebak(job, side, n, ilo, ihi, scale, m, v, ldv, info)
CGEBAK
Definition cgebak.f:129
subroutine cgebal(job, n, a, lda, ilo, ihi, scale, info)
CGEBAL
Definition cgebal.f:163
subroutine cgehrd(n, ilo, ihi, a, lda, tau, work, lwork, info)
CGEHRD
Definition cgehrd.f:166
subroutine chseqr(job, compz, n, ilo, ihi, h, ldh, w, z, ldz, work, lwork, info)
CHSEQR
Definition chseqr.f:297
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:101
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
real function clange(norm, m, n, a, lda, work)
CLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition clange.f:113
subroutine clascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
CLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition clascl.f:142
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 ctrsen(job, compq, select, n, t, ldt, q, ldq, w, m, s, sep, work, lwork, info)
CTRSEN
Definition ctrsen.f:263
subroutine cunghr(n, ilo, ihi, a, lda, tau, work, lwork, info)
CUNGHR
Definition cunghr.f:125
Here is the call graph for this function:
Here is the caller graph for this function: