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

◆ zgeesx()

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

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

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

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