Contents


NAME

     zggesx - compute for a pair of N-by-N  complex  nonsymmetric
     matrices  (A,B),  the  generalized  eigenvalues, the complex
     Schur form (S,T),

SYNOPSIS

     SUBROUTINE ZGGESX(JOBVSL, JOBVSR, SORT, DELCTG, SENSE, N, A, LDA, B,
           LDB, SDIM, ALPHA, BETA, VSL, LDVSL, VSR, LDVSR, RCONDE, RCONDV,
           WORK, LWORK, RWORK, IWORK, LIWORK, BWORK, INFO)

     CHARACTER * 1 JOBVSL, JOBVSR, SORT, SENSE
     DOUBLE  COMPLEX  A(LDA,*),  B(LDB,*),   ALPHA(*),   BETA(*),
     VSL(LDVSL,*), VSR(LDVSR,*), WORK(*)
     INTEGER N, LDA, LDB, SDIM, LDVSL, LDVSR, LWORK, LIWORK, INFO
     INTEGER IWORK(*)
     LOGICAL DELCTG
     LOGICAL BWORK(*)
     DOUBLE PRECISION RCONDE(*), RCONDV(*), RWORK(*)

     SUBROUTINE ZGGESX_64(JOBVSL, JOBVSR, SORT, DELCTG, SENSE, N, A, LDA,
           B, LDB, SDIM, ALPHA, BETA, VSL, LDVSL, VSR, LDVSR, RCONDE,
           RCONDV, WORK, LWORK, RWORK, IWORK, LIWORK, BWORK, INFO)

     CHARACTER * 1 JOBVSL, JOBVSR, SORT, SENSE
     DOUBLE  COMPLEX  A(LDA,*),  B(LDB,*),   ALPHA(*),   BETA(*),
     VSL(LDVSL,*), VSR(LDVSR,*), WORK(*)
     INTEGER*8 N, LDA, LDB, SDIM, LDVSL,  LDVSR,  LWORK,  LIWORK,
     INFO
     INTEGER*8 IWORK(*)
     LOGICAL*8 DELCTG
     LOGICAL*8 BWORK(*)
     DOUBLE PRECISION RCONDE(*), RCONDV(*), RWORK(*)

  F95 INTERFACE
     SUBROUTINE GGESX(JOBVSL, JOBVSR, SORT, [DELCTG], SENSE, [N], A, [LDA],
            B, [LDB], SDIM, ALPHA, BETA, VSL, [LDVSL], VSR, [LDVSR], RCONDE,
            RCONDV, [WORK], [LWORK], [RWORK], [IWORK], [LIWORK], [BWORK],
            [INFO])

     CHARACTER(LEN=1) :: JOBVSL, JOBVSR, SORT, SENSE
     COMPLEX(8), DIMENSION(:) :: ALPHA, BETA, WORK
     COMPLEX(8), DIMENSION(:,:) :: A, B, VSL, VSR
     INTEGER :: N, LDA, LDB, SDIM, LDVSL, LDVSR,  LWORK,  LIWORK,
     INFO
     INTEGER, DIMENSION(:) :: IWORK
     LOGICAL :: DELCTG
     LOGICAL, DIMENSION(:) :: BWORK
     REAL(8), DIMENSION(:) :: RCONDE, RCONDV, RWORK
     SUBROUTINE GGESX_64(JOBVSL, JOBVSR, SORT, [DELCTG], SENSE, [N], A, [LDA],
            B, [LDB], SDIM, ALPHA, BETA, VSL, [LDVSL], VSR, [LDVSR], RCONDE,
            RCONDV, [WORK], [LWORK], [RWORK], [IWORK], [LIWORK], [BWORK],
            [INFO])

     CHARACTER(LEN=1) :: JOBVSL, JOBVSR, SORT, SENSE
     COMPLEX(8), DIMENSION(:) :: ALPHA, BETA, WORK
     COMPLEX(8), DIMENSION(:,:) :: A, B, VSL, VSR
     INTEGER(8) ::  N,  LDA,  LDB,  SDIM,  LDVSL,  LDVSR,  LWORK,
     LIWORK, INFO
     INTEGER(8), DIMENSION(:) :: IWORK
     LOGICAL(8) :: DELCTG
     LOGICAL(8), DIMENSION(:) :: BWORK
     REAL(8), DIMENSION(:) :: RCONDE, RCONDV, RWORK

  C INTERFACE
     #include <sunperf.h>

     void   zggesx(char   jobvsl,   char   jobvsr,   char   sort,
               int(*delctg)(doublecomplex,doublecomplex),    char
               sense, int n, doublecomplex  *a,  int  lda,  doub-
               lecomplex  *b,  int  ldb, int *sdim, doublecomplex
               *alpha, doublecomplex *beta,  doublecomplex  *vsl,
               int  ldvsl,  doublecomplex *vsr, int ldvsr, double
               *rconde, double *rcondv, int *info);

     void  zggesx_64(char  jobvsl,  char   jobvsr,   char   sort,
               long(*delctg)(doublecomplex,doublecomplex),   char
               sense, long n, doublecomplex *a, long  lda,  doub-
               lecomplex  *b, long ldb, long *sdim, doublecomplex
               *alpha, doublecomplex *beta,  doublecomplex  *vsl,
               long ldvsl, doublecomplex *vsr, long ldvsr, double
               *rconde, double *rcondv, long *info);

PURPOSE

     zggesx computes for a pair of  N-by-N  complex  nonsymmetric
     matrices  (A,B),  the  generalized  eigenvalues, the complex
     Schur form (S,T), and, optionally,  the  left  and/or  right
     matrices  of  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 diag-
     onal blocks of the upper triangular matrix S and  the  upper
     triangular  matrix T; computes a reciprocal condition number
     for the average of the selected  eigenvalues  (RCONDE);  and
     computes  a  reciprocal  condition  number for the right and
     left  deflating  subspaces  corresponding  to  the  selected
     eigenvalues  (RCONDV).  The  leading  columns of VSL and VSR
     then form an orthonormal basis for  the  corresponding  left
     and right eigenspaces (deflating subspaces).

     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 or for both being zero.

     A pair of matrices (S,T) is  in  generalized  complex  Schur
     form if T is upper triangular with non-negative diagonal and
     S is upper triangular.

ARGUMENTS

     JOBVSL (input)
               = 'N':  do not compute the left Schur vectors;
               = 'V':  compute the left Schur vectors.

     JOBVSR (input)
               = 'N':  do not compute the right Schur vectors;
               = 'V':  compute the right Schur vectors.

     SORT (input)
               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 DELCTG).

     DELCTG (input)
               LOGICAL FUNCTION of two DOUBLE  COMPLEX  arguments
               DELCTG  must  be  declared EXTERNAL in the calling
               subroutine.  If SORT = 'N', DELCTG is  not  refer-
               enced.   If  SORT  = 'S', DELCTG is used to select
               eigenvalues to sort to the top left of  the  Schur
               form.  Note that a selected complex eigenvalue may
               no  longer  satisfy   DELCTG(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+3 see INFO below).

     SENSE (input)
               Determines which reciprocal condition numbers  are
               computed.  = 'N' : None are computed;
               = 'E' : Computed for average  of  selected  eigen-
               values only;
               = 'V' : Computed for selected deflating  subspaces
               only;
               = 'B' : Computed for both.  If SENSE =  'E',  'V',
               or 'B', SORT must equal 'S'.

     N (input) The order of the matrices A, B, VSL, and  VSR.   N
               >= 0.

     A (input/output)
               DOUBLE COMPLEX 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.

     LDA (input)
               The leading dimension of A.  LDA >= max(1,N).

     B (input/output)
               DOUBLE COMPLEX 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.

     LDB (input)
               The leading dimension of B.  LDB >= max(1,N).

     SDIM (output)
               If SORT = 'N', SDIM = 0.  If SORT =  'S',  SDIM  =
               number  of  eigenvalues  (after sorting) for which
               DELCTG is true.

     ALPHA (output)
               On exit, ALPHA(j)/BETA(j), j=1,...,N, will be  the
               generalized     eigenvalues.      ALPHA(j)     and
               BETA(j),j=1,...,N  are the diagonals of  the  com-
               plex  Schur  form  (S,T).   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).

     BETA (output)
               See description of ALPHA.

     VSL (input)
               If JOBVSL = 'V', VSL will contain the  left  Schur
               vectors.  Not referenced if JOBVSL = 'N'.

     LDVSL (input)
               The leading dimension of  the  matrix  VSL.  LDVSL
               >=1, and if JOBVSL = 'V', LDVSL >= N.

     VSR (input)
               If JOBVSR = 'V', VSR will contain the right  Schur
               vectors.  Not referenced if JOBVSR = 'N'.

     LDVSR (input)
               The leading dimension of the matrix VSR. LDVSR  >=
               1, and if JOBVSR = 'V', LDVSR >= N.

     RCONDE (output)
               If SENSE = 'E' or  'B',  RCONDE(1)  and  RCONDE(2)
               contain  the  reciprocal condition numbers for the
               average of the selected eigenvalues.   Not  refer-
               enced if SENSE = 'N' or 'V'.

     RCONDV (output)
               If SENSE = 'V' or  'B',  RCONDV(1)  and  RCONDV(2)
               contain  the  reciprocal  condition number for the
               selected deflating subspaces.  Not  referenced  if
               SENSE = 'N' or 'E'.

     WORK (workspace)
               On exit, if INFO = 0, WORK(1) returns the  optimal
               LWORK.

     LWORK (input)
               The dimension of the array WORK.   LWORK  >=  2*N.
               If  SENSE  =  'E',  'V', or 'B', LWORK >= MAX(2*N,
               2*SDIM*(N-SDIM)).
     RWORK (workspace)
               dimension(8*N) Real workspace.

     IWORK (workspace/output)
               Not referenced if SENSE = 'N'.  On exit, if INFO =
               0, IWORK(1) returns the optimal LIWORK.

     LIWORK (input)
               The dimension of the array WORK. LIWORK >= N+2.

     BWORK (workspace)
               dimension(N) Not referenced if SORT = 'N'.

     INFO (output)
               = 0:  successful exit
               < 0:  if INFO = -i, the i-th argument had an ille-
               gal 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 CHGEQZ
               =N+2: after reordering, roundoff changed values of
               some  complex  eigenvalues  so that leading eigen-
               values in the Generalized  Schur  form  no  longer
               satisfy  DELCTG=.TRUE.   This could also be caused
               due  to  scaling.   =N+3:  reordering  failed   in
               CTGSEN.