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

◆ zgges()

subroutine zgges ( character  jobvsl,
character  jobvsr,
character  sort,
external  selctg,
integer  n,
complex*16, dimension( lda, * )  a,
integer  lda,
complex*16, dimension( ldb, * )  b,
integer  ldb,
integer  sdim,
complex*16, dimension( * )  alpha,
complex*16, dimension( * )  beta,
complex*16, dimension( ldvsl, * )  vsl,
integer  ldvsl,
complex*16, dimension( ldvsr, * )  vsr,
integer  ldvsr,
complex*16, dimension( * )  work,
integer  lwork,
double precision, dimension( * )  rwork,
logical, dimension( * )  bwork,
integer  info 
)

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

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

Purpose:
 ZGGES computes for a pair of N-by-N complex nonsymmetric matrices
 (A,B), the generalized eigenvalues, the generalized complex Schur
 form (S, T), and optionally left and/or right Schur vectors (VSL
 and VSR). This gives the generalized Schur factorization

         (A,B) = ( (VSL)*S*(VSR)**H, (VSL)*T*(VSR)**H )

 where (VSR)**H is the conjugate-transpose of VSR.

 Optionally, it also orders the eigenvalues so that a selected cluster
 of eigenvalues appears in the leading diagonal blocks of the upper
 triangular matrix S and the upper triangular matrix T. The leading
 columns of VSL and VSR then form an unitary basis for the
 corresponding left and right eigenspaces (deflating subspaces).

 (If only the generalized eigenvalues are needed, use the driver
 ZGGEV instead, which is faster.)

 A generalized eigenvalue for a pair of matrices (A,B) is a scalar w
 or a ratio alpha/beta = w, such that  A - w*B is singular.  It is
 usually represented as the pair (alpha,beta), as there is a
 reasonable interpretation for beta=0, and even for both being zero.

 A pair of matrices (S,T) is in generalized complex Schur form if S
 and T are upper triangular and, in addition, the diagonal elements
 of T are non-negative real numbers.
Parameters
[in]JOBVSL
          JOBVSL is CHARACTER*1
          = 'N':  do not compute the left Schur vectors;
          = 'V':  compute the left Schur vectors.
[in]JOBVSR
          JOBVSR is CHARACTER*1
          = 'N':  do not compute the right Schur vectors;
          = 'V':  compute the right Schur vectors.
[in]SORT
          SORT is CHARACTER*1
          Specifies whether or not to order the eigenvalues on the
          diagonal of the generalized Schur form.
          = 'N':  Eigenvalues are not ordered;
          = 'S':  Eigenvalues are ordered (see SELCTG).
[in]SELCTG
          SELCTG is a LOGICAL FUNCTION of two COMPLEX*16 arguments
          SELCTG must be declared EXTERNAL in the calling subroutine.
          If SORT = 'N', SELCTG is not referenced.
          If SORT = 'S', SELCTG is used to select eigenvalues to sort
          to the top left of the Schur form.
          An eigenvalue ALPHA(j)/BETA(j) is selected if
          SELCTG(ALPHA(j),BETA(j)) is true.

          Note that a selected complex eigenvalue may no longer satisfy
          SELCTG(ALPHA(j),BETA(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 matrices A, B, VSL, and VSR.  N >= 0.
[in,out]A
          A is COMPLEX*16 array, dimension (LDA, N)
          On entry, the first of the pair of matrices.
          On exit, A has been overwritten by its generalized Schur
          form S.
[in]LDA
          LDA is INTEGER
          The leading dimension of A.  LDA >= max(1,N).
[in,out]B
          B is COMPLEX*16 array, dimension (LDB, N)
          On entry, the second of the pair of matrices.
          On exit, B has been overwritten by its generalized Schur
          form T.
[in]LDB
          LDB is INTEGER
          The leading dimension of B.  LDB >= max(1,N).
[out]SDIM
          SDIM is INTEGER
          If SORT = 'N', SDIM = 0.
          If SORT = 'S', SDIM = number of eigenvalues (after sorting)
          for which SELCTG is true.
[out]ALPHA
          ALPHA is COMPLEX*16 array, dimension (N)
[out]BETA
          BETA is COMPLEX*16 array, dimension (N)
          On exit,  ALPHA(j)/BETA(j), j=1,...,N, will be the
          generalized eigenvalues.  ALPHA(j), j=1,...,N  and  BETA(j),
          j=1,...,N  are the diagonals of the complex Schur form (A,B)
          output by ZGGES. The  BETA(j) will be non-negative real.

          Note: the quotients ALPHA(j)/BETA(j) may easily over- or
          underflow, and BETA(j) may even be zero.  Thus, the user
          should avoid naively computing the ratio alpha/beta.
          However, ALPHA will be always less than and usually
          comparable with norm(A) in magnitude, and BETA always less
          than and usually comparable with norm(B).
[out]VSL
          VSL is COMPLEX*16 array, dimension (LDVSL,N)
          If JOBVSL = 'V', VSL will contain the left Schur vectors.
          Not referenced if JOBVSL = 'N'.
[in]LDVSL
          LDVSL is INTEGER
          The leading dimension of the matrix VSL. LDVSL >= 1, and
          if JOBVSL = 'V', LDVSL >= N.
[out]VSR
          VSR is COMPLEX*16 array, dimension (LDVSR,N)
          If JOBVSR = 'V', VSR will contain the right Schur vectors.
          Not referenced if JOBVSR = 'N'.
[in]LDVSR
          LDVSR is INTEGER
          The leading dimension of the matrix VSR. LDVSR >= 1, and
          if JOBVSR = 'V', LDVSR >= 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 (8*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.
          =1,...,N:
                The QZ iteration failed.  (A,B) are not in Schur
                form, but ALPHA(j) and BETA(j) should be correct for
                j=INFO+1,...,N.
          > N:  =N+1: other than QZ iteration failed in ZHGEQZ
                =N+2: after reordering, roundoff changed values of
                      some complex eigenvalues so that leading
                      eigenvalues in the Generalized Schur form no
                      longer satisfy SELCTG=.TRUE.  This could also
                      be caused due to scaling.
                =N+3: reordering failed in ZTGSEN.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 267 of file zgges.f.

270*
271* -- LAPACK driver routine --
272* -- LAPACK is a software package provided by Univ. of Tennessee, --
273* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
274*
275* .. Scalar Arguments ..
276 CHARACTER JOBVSL, JOBVSR, SORT
277 INTEGER INFO, LDA, LDB, LDVSL, LDVSR, LWORK, N, SDIM
278* ..
279* .. Array Arguments ..
280 LOGICAL BWORK( * )
281 DOUBLE PRECISION RWORK( * )
282 COMPLEX*16 A( LDA, * ), ALPHA( * ), B( LDB, * ),
283 $ BETA( * ), VSL( LDVSL, * ), VSR( LDVSR, * ),
284 $ WORK( * )
285* ..
286* .. Function Arguments ..
287 LOGICAL SELCTG
288 EXTERNAL selctg
289* ..
290*
291* =====================================================================
292*
293* .. Parameters ..
294 DOUBLE PRECISION ZERO, ONE
295 parameter( zero = 0.0d0, one = 1.0d0 )
296 COMPLEX*16 CZERO, CONE
297 parameter( czero = ( 0.0d0, 0.0d0 ),
298 $ cone = ( 1.0d0, 0.0d0 ) )
299* ..
300* .. Local Scalars ..
301 LOGICAL CURSL, ILASCL, ILBSCL, ILVSL, ILVSR, LASTSL,
302 $ LQUERY, WANTST
303 INTEGER I, ICOLS, IERR, IHI, IJOBVL, IJOBVR, ILEFT,
304 $ ILO, IRIGHT, IROWS, IRWRK, ITAU, IWRK, LWKMIN,
305 $ LWKOPT
306 DOUBLE PRECISION ANRM, ANRMTO, BIGNUM, BNRM, BNRMTO, EPS, PVSL,
307 $ PVSR, SMLNUM
308* ..
309* .. Local Arrays ..
310 INTEGER IDUM( 1 )
311 DOUBLE PRECISION DIF( 2 )
312* ..
313* .. External Subroutines ..
314 EXTERNAL xerbla, zgeqrf, zggbak, zggbal, zgghrd, zhgeqz,
316* ..
317* .. External Functions ..
318 LOGICAL LSAME
319 INTEGER ILAENV
320 DOUBLE PRECISION DLAMCH, ZLANGE
321 EXTERNAL lsame, ilaenv, dlamch, zlange
322* ..
323* .. Intrinsic Functions ..
324 INTRINSIC max, sqrt
325* ..
326* .. Executable Statements ..
327*
328* Decode the input arguments
329*
330 IF( lsame( jobvsl, 'N' ) ) THEN
331 ijobvl = 1
332 ilvsl = .false.
333 ELSE IF( lsame( jobvsl, 'V' ) ) THEN
334 ijobvl = 2
335 ilvsl = .true.
336 ELSE
337 ijobvl = -1
338 ilvsl = .false.
339 END IF
340*
341 IF( lsame( jobvsr, 'N' ) ) THEN
342 ijobvr = 1
343 ilvsr = .false.
344 ELSE IF( lsame( jobvsr, 'V' ) ) THEN
345 ijobvr = 2
346 ilvsr = .true.
347 ELSE
348 ijobvr = -1
349 ilvsr = .false.
350 END IF
351*
352 wantst = lsame( sort, 'S' )
353*
354* Test the input arguments
355*
356 info = 0
357 lquery = ( lwork.EQ.-1 )
358 IF( ijobvl.LE.0 ) THEN
359 info = -1
360 ELSE IF( ijobvr.LE.0 ) THEN
361 info = -2
362 ELSE IF( ( .NOT.wantst ) .AND. ( .NOT.lsame( sort, 'N' ) ) ) THEN
363 info = -3
364 ELSE IF( n.LT.0 ) THEN
365 info = -5
366 ELSE IF( lda.LT.max( 1, n ) ) THEN
367 info = -7
368 ELSE IF( ldb.LT.max( 1, n ) ) THEN
369 info = -9
370 ELSE IF( ldvsl.LT.1 .OR. ( ilvsl .AND. ldvsl.LT.n ) ) THEN
371 info = -14
372 ELSE IF( ldvsr.LT.1 .OR. ( ilvsr .AND. ldvsr.LT.n ) ) THEN
373 info = -16
374 END IF
375*
376* Compute workspace
377* (Note: Comments in the code beginning "Workspace:" describe the
378* minimal amount of workspace needed at that point in the code,
379* as well as the preferred amount for good performance.
380* NB refers to the optimal block size for the immediately
381* following subroutine, as returned by ILAENV.)
382*
383 IF( info.EQ.0 ) THEN
384 lwkmin = max( 1, 2*n )
385 lwkopt = max( 1, n + n*ilaenv( 1, 'ZGEQRF', ' ', n, 1, n, 0 ) )
386 lwkopt = max( lwkopt, n +
387 $ n*ilaenv( 1, 'ZUNMQR', ' ', n, 1, n, -1 ) )
388 IF( ilvsl ) THEN
389 lwkopt = max( lwkopt, n +
390 $ n*ilaenv( 1, 'ZUNGQR', ' ', n, 1, n, -1 ) )
391 END IF
392 work( 1 ) = lwkopt
393*
394 IF( lwork.LT.lwkmin .AND. .NOT.lquery )
395 $ info = -18
396 END IF
397*
398 IF( info.NE.0 ) THEN
399 CALL xerbla( 'ZGGES ', -info )
400 RETURN
401 ELSE IF( lquery ) THEN
402 RETURN
403 END IF
404*
405* Quick return if possible
406*
407 IF( n.EQ.0 ) THEN
408 sdim = 0
409 RETURN
410 END IF
411*
412* Get machine constants
413*
414 eps = dlamch( 'P' )
415 smlnum = dlamch( 'S' )
416 bignum = one / smlnum
417 smlnum = sqrt( smlnum ) / eps
418 bignum = one / smlnum
419*
420* Scale A if max element outside range [SMLNUM,BIGNUM]
421*
422 anrm = zlange( 'M', n, n, a, lda, rwork )
423 ilascl = .false.
424 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
425 anrmto = smlnum
426 ilascl = .true.
427 ELSE IF( anrm.GT.bignum ) THEN
428 anrmto = bignum
429 ilascl = .true.
430 END IF
431*
432 IF( ilascl )
433 $ CALL zlascl( 'G', 0, 0, anrm, anrmto, n, n, a, lda, ierr )
434*
435* Scale B if max element outside range [SMLNUM,BIGNUM]
436*
437 bnrm = zlange( 'M', n, n, b, ldb, rwork )
438 ilbscl = .false.
439 IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
440 bnrmto = smlnum
441 ilbscl = .true.
442 ELSE IF( bnrm.GT.bignum ) THEN
443 bnrmto = bignum
444 ilbscl = .true.
445 END IF
446*
447 IF( ilbscl )
448 $ CALL zlascl( 'G', 0, 0, bnrm, bnrmto, n, n, b, ldb, ierr )
449*
450* Permute the matrix to make it more nearly triangular
451* (Real Workspace: need 6*N)
452*
453 ileft = 1
454 iright = n + 1
455 irwrk = iright + n
456 CALL zggbal( 'P', n, a, lda, b, ldb, ilo, ihi, rwork( ileft ),
457 $ rwork( iright ), rwork( irwrk ), ierr )
458*
459* Reduce B to triangular form (QR decomposition of B)
460* (Complex Workspace: need N, prefer N*NB)
461*
462 irows = ihi + 1 - ilo
463 icols = n + 1 - ilo
464 itau = 1
465 iwrk = itau + irows
466 CALL zgeqrf( irows, icols, b( ilo, ilo ), ldb, work( itau ),
467 $ work( iwrk ), lwork+1-iwrk, ierr )
468*
469* Apply the orthogonal transformation to matrix A
470* (Complex Workspace: need N, prefer N*NB)
471*
472 CALL zunmqr( 'L', 'C', irows, icols, irows, b( ilo, ilo ), ldb,
473 $ work( itau ), a( ilo, ilo ), lda, work( iwrk ),
474 $ lwork+1-iwrk, ierr )
475*
476* Initialize VSL
477* (Complex Workspace: need N, prefer N*NB)
478*
479 IF( ilvsl ) THEN
480 CALL zlaset( 'Full', n, n, czero, cone, vsl, ldvsl )
481 IF( irows.GT.1 ) THEN
482 CALL zlacpy( 'L', irows-1, irows-1, b( ilo+1, ilo ), ldb,
483 $ vsl( ilo+1, ilo ), ldvsl )
484 END IF
485 CALL zungqr( irows, irows, irows, vsl( ilo, ilo ), ldvsl,
486 $ work( itau ), work( iwrk ), lwork+1-iwrk, ierr )
487 END IF
488*
489* Initialize VSR
490*
491 IF( ilvsr )
492 $ CALL zlaset( 'Full', n, n, czero, cone, vsr, ldvsr )
493*
494* Reduce to generalized Hessenberg form
495* (Workspace: none needed)
496*
497 CALL zgghrd( jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb, vsl,
498 $ ldvsl, vsr, ldvsr, ierr )
499*
500 sdim = 0
501*
502* Perform QZ algorithm, computing Schur vectors if desired
503* (Complex Workspace: need N)
504* (Real Workspace: need N)
505*
506 iwrk = itau
507 CALL zhgeqz( 'S', jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb,
508 $ alpha, beta, vsl, ldvsl, vsr, ldvsr, work( iwrk ),
509 $ lwork+1-iwrk, rwork( irwrk ), ierr )
510 IF( ierr.NE.0 ) THEN
511 IF( ierr.GT.0 .AND. ierr.LE.n ) THEN
512 info = ierr
513 ELSE IF( ierr.GT.n .AND. ierr.LE.2*n ) THEN
514 info = ierr - n
515 ELSE
516 info = n + 1
517 END IF
518 GO TO 30
519 END IF
520*
521* Sort eigenvalues ALPHA/BETA if desired
522* (Workspace: none needed)
523*
524 IF( wantst ) THEN
525*
526* Undo scaling on eigenvalues before selecting
527*
528 IF( ilascl )
529 $ CALL zlascl( 'G', 0, 0, anrm, anrmto, n, 1, alpha, n, ierr )
530 IF( ilbscl )
531 $ CALL zlascl( 'G', 0, 0, bnrm, bnrmto, n, 1, beta, n, ierr )
532*
533* Select eigenvalues
534*
535 DO 10 i = 1, n
536 bwork( i ) = selctg( alpha( i ), beta( i ) )
537 10 CONTINUE
538*
539 CALL ztgsen( 0, ilvsl, ilvsr, bwork, n, a, lda, b, ldb, alpha,
540 $ beta, vsl, ldvsl, vsr, ldvsr, sdim, pvsl, pvsr,
541 $ dif, work( iwrk ), lwork-iwrk+1, idum, 1, ierr )
542 IF( ierr.EQ.1 )
543 $ info = n + 3
544*
545 END IF
546*
547* Apply back-permutation to VSL and VSR
548* (Workspace: none needed)
549*
550 IF( ilvsl )
551 $ CALL zggbak( 'P', 'L', n, ilo, ihi, rwork( ileft ),
552 $ rwork( iright ), n, vsl, ldvsl, ierr )
553 IF( ilvsr )
554 $ CALL zggbak( 'P', 'R', n, ilo, ihi, rwork( ileft ),
555 $ rwork( iright ), n, vsr, ldvsr, ierr )
556*
557* Undo scaling
558*
559 IF( ilascl ) THEN
560 CALL zlascl( 'U', 0, 0, anrmto, anrm, n, n, a, lda, ierr )
561 CALL zlascl( 'G', 0, 0, anrmto, anrm, n, 1, alpha, n, ierr )
562 END IF
563*
564 IF( ilbscl ) THEN
565 CALL zlascl( 'U', 0, 0, bnrmto, bnrm, n, n, b, ldb, ierr )
566 CALL zlascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
567 END IF
568*
569 IF( wantst ) THEN
570*
571* Check if reordering is correct
572*
573 lastsl = .true.
574 sdim = 0
575 DO 20 i = 1, n
576 cursl = selctg( alpha( i ), beta( i ) )
577 IF( cursl )
578 $ sdim = sdim + 1
579 IF( cursl .AND. .NOT.lastsl )
580 $ info = n + 2
581 lastsl = cursl
582 20 CONTINUE
583*
584 END IF
585*
586 30 CONTINUE
587*
588 work( 1 ) = lwkopt
589*
590 RETURN
591*
592* End of ZGGES
593*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgeqrf(m, n, a, lda, tau, work, lwork, info)
ZGEQRF
Definition zgeqrf.f:146
subroutine zggbak(job, side, n, ilo, ihi, lscale, rscale, m, v, ldv, info)
ZGGBAK
Definition zggbak.f:148
subroutine zggbal(job, n, a, lda, b, ldb, ilo, ihi, lscale, rscale, work, info)
ZGGBAL
Definition zggbal.f:177
subroutine zgghrd(compq, compz, n, ilo, ihi, a, lda, b, ldb, q, ldq, z, ldz, info)
ZGGHRD
Definition zgghrd.f:204
subroutine zhgeqz(job, compq, compz, n, ilo, ihi, h, ldh, t, ldt, alpha, beta, q, ldq, z, ldz, work, lwork, rwork, info)
ZHGEQZ
Definition zhgeqz.f:284
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
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:115
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:143
subroutine zlaset(uplo, m, n, alpha, beta, a, lda)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition zlaset.f:106
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine ztgsen(ijob, wantq, wantz, select, n, a, lda, b, ldb, alpha, beta, q, ldq, z, ldz, m, pl, pr, dif, work, lwork, iwork, liwork, info)
ZTGSEN
Definition ztgsen.f:433
subroutine zungqr(m, n, k, a, lda, tau, work, lwork, info)
ZUNGQR
Definition zungqr.f:128
subroutine zunmqr(side, trans, m, n, k, a, lda, tau, c, ldc, work, lwork, info)
ZUNMQR
Definition zunmqr.f:167
Here is the call graph for this function:
Here is the caller graph for this function: