LAPACK 3.11.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 dlabad, xerbla, zgeqrf, zggbak, zggbal, zgghrd,
316 $ zunmqr
317* ..
318* .. External Functions ..
319 LOGICAL LSAME
320 INTEGER ILAENV
321 DOUBLE PRECISION DLAMCH, ZLANGE
322 EXTERNAL lsame, ilaenv, dlamch, zlange
323* ..
324* .. Intrinsic Functions ..
325 INTRINSIC max, sqrt
326* ..
327* .. Executable Statements ..
328*
329* Decode the input arguments
330*
331 IF( lsame( jobvsl, 'N' ) ) THEN
332 ijobvl = 1
333 ilvsl = .false.
334 ELSE IF( lsame( jobvsl, 'V' ) ) THEN
335 ijobvl = 2
336 ilvsl = .true.
337 ELSE
338 ijobvl = -1
339 ilvsl = .false.
340 END IF
341*
342 IF( lsame( jobvsr, 'N' ) ) THEN
343 ijobvr = 1
344 ilvsr = .false.
345 ELSE IF( lsame( jobvsr, 'V' ) ) THEN
346 ijobvr = 2
347 ilvsr = .true.
348 ELSE
349 ijobvr = -1
350 ilvsr = .false.
351 END IF
352*
353 wantst = lsame( sort, 'S' )
354*
355* Test the input arguments
356*
357 info = 0
358 lquery = ( lwork.EQ.-1 )
359 IF( ijobvl.LE.0 ) THEN
360 info = -1
361 ELSE IF( ijobvr.LE.0 ) THEN
362 info = -2
363 ELSE IF( ( .NOT.wantst ) .AND. ( .NOT.lsame( sort, 'N' ) ) ) THEN
364 info = -3
365 ELSE IF( n.LT.0 ) THEN
366 info = -5
367 ELSE IF( lda.LT.max( 1, n ) ) THEN
368 info = -7
369 ELSE IF( ldb.LT.max( 1, n ) ) THEN
370 info = -9
371 ELSE IF( ldvsl.LT.1 .OR. ( ilvsl .AND. ldvsl.LT.n ) ) THEN
372 info = -14
373 ELSE IF( ldvsr.LT.1 .OR. ( ilvsr .AND. ldvsr.LT.n ) ) THEN
374 info = -16
375 END IF
376*
377* Compute workspace
378* (Note: Comments in the code beginning "Workspace:" describe the
379* minimal amount of workspace needed at that point in the code,
380* as well as the preferred amount for good performance.
381* NB refers to the optimal block size for the immediately
382* following subroutine, as returned by ILAENV.)
383*
384 IF( info.EQ.0 ) THEN
385 lwkmin = max( 1, 2*n )
386 lwkopt = max( 1, n + n*ilaenv( 1, 'ZGEQRF', ' ', n, 1, n, 0 ) )
387 lwkopt = max( lwkopt, n +
388 $ n*ilaenv( 1, 'ZUNMQR', ' ', n, 1, n, -1 ) )
389 IF( ilvsl ) THEN
390 lwkopt = max( lwkopt, n +
391 $ n*ilaenv( 1, 'ZUNGQR', ' ', n, 1, n, -1 ) )
392 END IF
393 work( 1 ) = lwkopt
394*
395 IF( lwork.LT.lwkmin .AND. .NOT.lquery )
396 $ info = -18
397 END IF
398*
399 IF( info.NE.0 ) THEN
400 CALL xerbla( 'ZGGES ', -info )
401 RETURN
402 ELSE IF( lquery ) THEN
403 RETURN
404 END IF
405*
406* Quick return if possible
407*
408 IF( n.EQ.0 ) THEN
409 sdim = 0
410 RETURN
411 END IF
412*
413* Get machine constants
414*
415 eps = dlamch( 'P' )
416 smlnum = dlamch( 'S' )
417 bignum = one / smlnum
418 CALL dlabad( smlnum, bignum )
419 smlnum = sqrt( smlnum ) / eps
420 bignum = one / smlnum
421*
422* Scale A if max element outside range [SMLNUM,BIGNUM]
423*
424 anrm = zlange( 'M', n, n, a, lda, rwork )
425 ilascl = .false.
426 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
427 anrmto = smlnum
428 ilascl = .true.
429 ELSE IF( anrm.GT.bignum ) THEN
430 anrmto = bignum
431 ilascl = .true.
432 END IF
433*
434 IF( ilascl )
435 $ CALL zlascl( 'G', 0, 0, anrm, anrmto, n, n, a, lda, ierr )
436*
437* Scale B if max element outside range [SMLNUM,BIGNUM]
438*
439 bnrm = zlange( 'M', n, n, b, ldb, rwork )
440 ilbscl = .false.
441 IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
442 bnrmto = smlnum
443 ilbscl = .true.
444 ELSE IF( bnrm.GT.bignum ) THEN
445 bnrmto = bignum
446 ilbscl = .true.
447 END IF
448*
449 IF( ilbscl )
450 $ CALL zlascl( 'G', 0, 0, bnrm, bnrmto, n, n, b, ldb, ierr )
451*
452* Permute the matrix to make it more nearly triangular
453* (Real Workspace: need 6*N)
454*
455 ileft = 1
456 iright = n + 1
457 irwrk = iright + n
458 CALL zggbal( 'P', n, a, lda, b, ldb, ilo, ihi, rwork( ileft ),
459 $ rwork( iright ), rwork( irwrk ), ierr )
460*
461* Reduce B to triangular form (QR decomposition of B)
462* (Complex Workspace: need N, prefer N*NB)
463*
464 irows = ihi + 1 - ilo
465 icols = n + 1 - ilo
466 itau = 1
467 iwrk = itau + irows
468 CALL zgeqrf( irows, icols, b( ilo, ilo ), ldb, work( itau ),
469 $ work( iwrk ), lwork+1-iwrk, ierr )
470*
471* Apply the orthogonal transformation to matrix A
472* (Complex Workspace: need N, prefer N*NB)
473*
474 CALL zunmqr( 'L', 'C', irows, icols, irows, b( ilo, ilo ), ldb,
475 $ work( itau ), a( ilo, ilo ), lda, work( iwrk ),
476 $ lwork+1-iwrk, ierr )
477*
478* Initialize VSL
479* (Complex Workspace: need N, prefer N*NB)
480*
481 IF( ilvsl ) THEN
482 CALL zlaset( 'Full', n, n, czero, cone, vsl, ldvsl )
483 IF( irows.GT.1 ) THEN
484 CALL zlacpy( 'L', irows-1, irows-1, b( ilo+1, ilo ), ldb,
485 $ vsl( ilo+1, ilo ), ldvsl )
486 END IF
487 CALL zungqr( irows, irows, irows, vsl( ilo, ilo ), ldvsl,
488 $ work( itau ), work( iwrk ), lwork+1-iwrk, ierr )
489 END IF
490*
491* Initialize VSR
492*
493 IF( ilvsr )
494 $ CALL zlaset( 'Full', n, n, czero, cone, vsr, ldvsr )
495*
496* Reduce to generalized Hessenberg form
497* (Workspace: none needed)
498*
499 CALL zgghrd( jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb, vsl,
500 $ ldvsl, vsr, ldvsr, ierr )
501*
502 sdim = 0
503*
504* Perform QZ algorithm, computing Schur vectors if desired
505* (Complex Workspace: need N)
506* (Real Workspace: need N)
507*
508 iwrk = itau
509 CALL zhgeqz( 'S', jobvsl, jobvsr, n, ilo, ihi, a, lda, b, ldb,
510 $ alpha, beta, vsl, ldvsl, vsr, ldvsr, work( iwrk ),
511 $ lwork+1-iwrk, rwork( irwrk ), ierr )
512 IF( ierr.NE.0 ) THEN
513 IF( ierr.GT.0 .AND. ierr.LE.n ) THEN
514 info = ierr
515 ELSE IF( ierr.GT.n .AND. ierr.LE.2*n ) THEN
516 info = ierr - n
517 ELSE
518 info = n + 1
519 END IF
520 GO TO 30
521 END IF
522*
523* Sort eigenvalues ALPHA/BETA if desired
524* (Workspace: none needed)
525*
526 IF( wantst ) THEN
527*
528* Undo scaling on eigenvalues before selecting
529*
530 IF( ilascl )
531 $ CALL zlascl( 'G', 0, 0, anrm, anrmto, n, 1, alpha, n, ierr )
532 IF( ilbscl )
533 $ CALL zlascl( 'G', 0, 0, bnrm, bnrmto, n, 1, beta, n, ierr )
534*
535* Select eigenvalues
536*
537 DO 10 i = 1, n
538 bwork( i ) = selctg( alpha( i ), beta( i ) )
539 10 CONTINUE
540*
541 CALL ztgsen( 0, ilvsl, ilvsr, bwork, n, a, lda, b, ldb, alpha,
542 $ beta, vsl, ldvsl, vsr, ldvsr, sdim, pvsl, pvsr,
543 $ dif, work( iwrk ), lwork-iwrk+1, idum, 1, ierr )
544 IF( ierr.EQ.1 )
545 $ info = n + 3
546*
547 END IF
548*
549* Apply back-permutation to VSL and VSR
550* (Workspace: none needed)
551*
552 IF( ilvsl )
553 $ CALL zggbak( 'P', 'L', n, ilo, ihi, rwork( ileft ),
554 $ rwork( iright ), n, vsl, ldvsl, ierr )
555 IF( ilvsr )
556 $ CALL zggbak( 'P', 'R', n, ilo, ihi, rwork( ileft ),
557 $ rwork( iright ), n, vsr, ldvsr, ierr )
558*
559* Undo scaling
560*
561 IF( ilascl ) THEN
562 CALL zlascl( 'U', 0, 0, anrmto, anrm, n, n, a, lda, ierr )
563 CALL zlascl( 'G', 0, 0, anrmto, anrm, n, 1, alpha, n, ierr )
564 END IF
565*
566 IF( ilbscl ) THEN
567 CALL zlascl( 'U', 0, 0, bnrmto, bnrm, n, n, b, ldb, ierr )
568 CALL zlascl( 'G', 0, 0, bnrmto, bnrm, n, 1, beta, n, ierr )
569 END IF
570*
571 IF( wantst ) THEN
572*
573* Check if reordering is correct
574*
575 lastsl = .true.
576 sdim = 0
577 DO 20 i = 1, n
578 cursl = selctg( alpha( i ), beta( i ) )
579 IF( cursl )
580 $ sdim = sdim + 1
581 IF( cursl .AND. .NOT.lastsl )
582 $ info = n + 2
583 lastsl = cursl
584 20 CONTINUE
585*
586 END IF
587*
588 30 CONTINUE
589*
590 work( 1 ) = lwkopt
591*
592 RETURN
593*
594* End of ZGGES
595*
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:69
subroutine dlabad(SMALL, LARGE)
DLABAD
Definition: dlabad.f:74
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV
Definition: ilaenv.f:162
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:53
subroutine zggbal(JOB, N, A, LDA, B, LDB, ILO, IHI, LSCALE, RSCALE, WORK, INFO)
ZGGBAL
Definition: zggbal.f:177
subroutine zggbak(JOB, SIDE, N, ILO, IHI, LSCALE, RSCALE, M, V, LDV, INFO)
ZGGBAK
Definition: zggbak.f:148
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 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
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 zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:103
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
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 zgghrd(COMPQ, COMPZ, N, ILO, IHI, A, LDA, B, LDB, Q, LDQ, Z, LDZ, INFO)
ZGGHRD
Definition: zgghrd.f:204
subroutine zunmqr(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
ZUNMQR
Definition: zunmqr.f:167
subroutine zgeqrf(M, N, A, LDA, TAU, WORK, LWORK, INFO)
ZGEQRF VARIANT: left-looking Level 3 BLAS of the algorithm.
Definition: zgeqrf.f:152
Here is the call graph for this function:
Here is the caller graph for this function: