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

◆ zgees()

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

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

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

Purpose:
!>
!> ZGEES 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.
!> The leading columns of Z then form an orthonormal basis for the
!> invariant subspace corresponding to the selected eigenvalues.
!>
!> 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.
!>          The eigenvalue W(j) is selected if SELECT(W(j)) is true.
!> 
[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 has been 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; if
!>          JOBVS = 'V', LDVS >= N.
!> 
[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).
!>          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]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 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 193 of file zgees.f.

195*
196* -- LAPACK driver routine --
197* -- LAPACK is a software package provided by Univ. of Tennessee, --
198* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
199*
200* .. Scalar Arguments ..
201 CHARACTER JOBVS, SORT
202 INTEGER INFO, LDA, LDVS, LWORK, N, SDIM
203* ..
204* .. Array Arguments ..
205 LOGICAL BWORK( * )
206 DOUBLE PRECISION RWORK( * )
207 COMPLEX*16 A( LDA, * ), VS( LDVS, * ), W( * ), WORK( * )
208* ..
209* .. Function Arguments ..
210 LOGICAL SELECT
211 EXTERNAL SELECT
212* ..
213*
214* =====================================================================
215*
216* .. Parameters ..
217 DOUBLE PRECISION ZERO, ONE
218 parameter( zero = 0.0d0, one = 1.0d0 )
219* ..
220* .. Local Scalars ..
221 LOGICAL LQUERY, SCALEA, WANTST, WANTVS
222 INTEGER HSWORK, I, IBAL, ICOND, IERR, IEVAL, IHI, ILO,
223 $ ITAU, IWRK, MAXWRK, MINWRK
224 DOUBLE PRECISION ANRM, BIGNUM, CSCALE, EPS, S, SEP, SMLNUM
225* ..
226* .. Local Arrays ..
227 DOUBLE PRECISION DUM( 1 )
228* ..
229* .. External Subroutines ..
230 EXTERNAL xerbla, zcopy, zgebak, zgebal,
231 $ zgehrd,
233* ..
234* .. External Functions ..
235 LOGICAL LSAME
236 INTEGER ILAENV
237 DOUBLE PRECISION DLAMCH, ZLANGE
238 EXTERNAL lsame, ilaenv, dlamch, zlange
239* ..
240* .. Intrinsic Functions ..
241 INTRINSIC max, sqrt
242* ..
243* .. Executable Statements ..
244*
245* Test the input arguments
246*
247 info = 0
248 lquery = ( lwork.EQ.-1 )
249 wantvs = lsame( jobvs, 'V' )
250 wantst = lsame( sort, 'S' )
251 IF( ( .NOT.wantvs ) .AND. ( .NOT.lsame( jobvs, 'N' ) ) ) THEN
252 info = -1
253 ELSE IF( ( .NOT.wantst ) .AND.
254 $ ( .NOT.lsame( sort, 'N' ) ) ) THEN
255 info = -2
256 ELSE IF( n.LT.0 ) THEN
257 info = -4
258 ELSE IF( lda.LT.max( 1, n ) ) THEN
259 info = -6
260 ELSE IF( ldvs.LT.1 .OR. ( wantvs .AND. ldvs.LT.n ) ) THEN
261 info = -10
262 END IF
263*
264* Compute workspace
265* (Note: Comments in the code beginning "Workspace:" describe the
266* minimal amount of workspace needed at that point in the code,
267* as well as the preferred amount for good performance.
268* CWorkspace refers to complex workspace, and RWorkspace to real
269* workspace. NB refers to the optimal block size for the
270* immediately following subroutine, as returned by ILAENV.
271* HSWORK refers to the workspace preferred by ZHSEQR, as
272* calculated below. HSWORK is computed assuming ILO=1 and IHI=N,
273* the worst case.)
274*
275 IF( info.EQ.0 ) THEN
276 IF( n.EQ.0 ) THEN
277 minwrk = 1
278 maxwrk = 1
279 ELSE
280 maxwrk = n + n*ilaenv( 1, 'ZGEHRD', ' ', n, 1, n, 0 )
281 minwrk = 2*n
282*
283 CALL zhseqr( 'S', jobvs, n, 1, n, a, lda, w, vs, ldvs,
284 $ work, -1, ieval )
285 hswork = int( work( 1 ) )
286*
287 IF( .NOT.wantvs ) THEN
288 maxwrk = max( maxwrk, hswork )
289 ELSE
290 maxwrk = max( maxwrk, n + ( n - 1 )*ilaenv( 1,
291 $ 'ZUNGHR',
292 $ ' ', n, 1, n, -1 ) )
293 maxwrk = max( maxwrk, hswork )
294 END IF
295 END IF
296 work( 1 ) = maxwrk
297*
298 IF( lwork.LT.minwrk .AND. .NOT.lquery ) THEN
299 info = -12
300 END IF
301 END IF
302*
303 IF( info.NE.0 ) THEN
304 CALL xerbla( 'ZGEES ', -info )
305 RETURN
306 ELSE IF( lquery ) THEN
307 RETURN
308 END IF
309*
310* Quick return if possible
311*
312 IF( n.EQ.0 ) THEN
313 sdim = 0
314 RETURN
315 END IF
316*
317* Get machine constants
318*
319 eps = dlamch( 'P' )
320 smlnum = dlamch( 'S' )
321 bignum = one / smlnum
322 smlnum = sqrt( smlnum ) / eps
323 bignum = one / smlnum
324*
325* Scale A if max element outside range [SMLNUM,BIGNUM]
326*
327 anrm = zlange( 'M', n, n, a, lda, dum )
328 scalea = .false.
329 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
330 scalea = .true.
331 cscale = smlnum
332 ELSE IF( anrm.GT.bignum ) THEN
333 scalea = .true.
334 cscale = bignum
335 END IF
336 IF( scalea )
337 $ CALL zlascl( 'G', 0, 0, anrm, cscale, n, n, a, lda, ierr )
338*
339* Permute the matrix to make it more nearly triangular
340* (CWorkspace: none)
341* (RWorkspace: need N)
342*
343 ibal = 1
344 CALL zgebal( 'P', n, a, lda, ilo, ihi, rwork( ibal ), ierr )
345*
346* Reduce to upper Hessenberg form
347* (CWorkspace: need 2*N, prefer N+N*NB)
348* (RWorkspace: none)
349*
350 itau = 1
351 iwrk = n + itau
352 CALL zgehrd( n, ilo, ihi, a, lda, work( itau ), work( iwrk ),
353 $ lwork-iwrk+1, ierr )
354*
355 IF( wantvs ) THEN
356*
357* Copy Householder vectors to VS
358*
359 CALL zlacpy( 'L', n, n, a, lda, vs, ldvs )
360*
361* Generate unitary matrix in VS
362* (CWorkspace: need 2*N-1, prefer N+(N-1)*NB)
363* (RWorkspace: none)
364*
365 CALL zunghr( n, ilo, ihi, vs, ldvs, work( itau ),
366 $ work( iwrk ),
367 $ lwork-iwrk+1, ierr )
368 END IF
369*
370 sdim = 0
371*
372* Perform QR iteration, accumulating Schur vectors in VS if desired
373* (CWorkspace: need 1, prefer HSWORK (see comments) )
374* (RWorkspace: none)
375*
376 iwrk = itau
377 CALL zhseqr( 'S', jobvs, n, ilo, ihi, a, lda, w, vs, ldvs,
378 $ work( iwrk ), lwork-iwrk+1, ieval )
379 IF( ieval.GT.0 )
380 $ info = ieval
381*
382* Sort eigenvalues if desired
383*
384 IF( wantst .AND. info.EQ.0 ) THEN
385 IF( scalea )
386 $ CALL zlascl( 'G', 0, 0, cscale, anrm, n, 1, w, n, ierr )
387 DO 10 i = 1, n
388 bwork( i ) = SELECT( w( i ) )
389 10 CONTINUE
390*
391* Reorder eigenvalues and transform Schur vectors
392* (CWorkspace: none)
393* (RWorkspace: none)
394*
395 CALL ztrsen( 'N', jobvs, bwork, n, a, lda, vs, ldvs, w,
396 $ sdim,
397 $ s, sep, work( iwrk ), lwork-iwrk+1, icond )
398 END IF
399*
400 IF( wantvs ) THEN
401*
402* Undo balancing
403* (CWorkspace: none)
404* (RWorkspace: need N)
405*
406 CALL zgebak( 'P', 'R', n, ilo, ihi, rwork( ibal ), n, vs,
407 $ ldvs,
408 $ ierr )
409 END IF
410*
411 IF( scalea ) THEN
412*
413* Undo scaling for the Schur form of A
414*
415 CALL zlascl( 'U', 0, 0, cscale, anrm, n, n, a, lda, ierr )
416 CALL zcopy( n, a, lda+1, w, 1 )
417 END IF
418*
419 work( 1 ) = maxwrk
420 RETURN
421*
422* End of ZGEES
423*
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
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: