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

◆ cgees()

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

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

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

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