Go to main content
Oracle Developer Studio 12.5 Man Pages

Exit Print View

Updated: June 2017
 
 

dggevx (3p)

Name

dggevx - N real nonsymmetric matrices (A,B), and, optionally, the left and/or right generalized eigenvectors

Synopsis

SUBROUTINE DGGEVX(BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, B, LDB,
ALPHAR, ALPHAI, BETA, VL, LDVL, VR, LDVR, ILO, IHI, LSCALE,
RSCALE, ABNRM, BBNRM, RCONDE, RCONDV, WORK, LWORK, IWORK, BWORK,
INFO)

CHARACTER*1 BALANC, JOBVL, JOBVR, SENSE
INTEGER N, LDA, LDB, LDVL, LDVR, ILO, IHI, LWORK, INFO
INTEGER IWORK(*)
LOGICAL BWORK(*)
DOUBLE PRECISION ABNRM, BBNRM
DOUBLE PRECISION A(LDA,*),  B(LDB,*),  ALPHAR(*),  ALPHAI(*),  BETA(*),
VL(LDVL,*),  VR(LDVR,*),  LSCALE(*),  RSCALE(*),  RCONDE(*), RCONDV(*),
WORK(*)

SUBROUTINE DGGEVX_64(BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, B, LDB,
ALPHAR, ALPHAI, BETA, VL, LDVL, VR, LDVR, ILO, IHI, LSCALE,
RSCALE, ABNRM, BBNRM, RCONDE, RCONDV, WORK, LWORK, IWORK, BWORK,
INFO)

CHARACTER*1 BALANC, JOBVL, JOBVR, SENSE
INTEGER*8 N, LDA, LDB, LDVL, LDVR, ILO, IHI, LWORK, INFO
INTEGER*8 IWORK(*)
LOGICAL*8 BWORK(*)
DOUBLE PRECISION ABNRM, BBNRM
DOUBLE PRECISION A(LDA,*),  B(LDB,*),  ALPHAR(*),  ALPHAI(*),  BETA(*),
VL(LDVL,*),  VR(LDVR,*),  LSCALE(*),  RSCALE(*),  RCONDE(*), RCONDV(*),
WORK(*)




F95 INTERFACE
SUBROUTINE GGEVX(BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, B, LDB,
ALPHAR, ALPHAI, BETA, VL, LDVL, VR, LDVR, ILO, IHI, LSCALE,
RSCALE, ABNRM, BBNRM, RCONDE, RCONDV, WORK, LWORK, IWORK,
BWORK, INFO)

CHARACTER(LEN=1) :: BALANC, JOBVL, JOBVR, SENSE
INTEGER :: N, LDA, LDB, LDVL, LDVR, ILO, IHI, LWORK, INFO
INTEGER, DIMENSION(:) :: IWORK
LOGICAL, DIMENSION(:) :: BWORK
REAL(8) :: ABNRM, BBNRM
REAL(8), DIMENSION(:) :: ALPHAR, ALPHAI, BETA, LSCALE, RSCALE,  RCONDE,
RCONDV, WORK
REAL(8), DIMENSION(:,:) :: A, B, VL, VR

SUBROUTINE GGEVX_64(BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, B,
LDB, ALPHAR, ALPHAI, BETA, VL, LDVL, VR, LDVR, ILO, IHI,
LSCALE, RSCALE, ABNRM, BBNRM, RCONDE, RCONDV, WORK, LWORK,
IWORK, BWORK, INFO)

CHARACTER(LEN=1) :: BALANC, JOBVL, JOBVR, SENSE
INTEGER(8) :: N, LDA, LDB, LDVL, LDVR, ILO, IHI, LWORK, INFO
INTEGER(8), DIMENSION(:) :: IWORK
LOGICAL(8), DIMENSION(:) :: BWORK
REAL(8) :: ABNRM, BBNRM
REAL(8),  DIMENSION(:) :: ALPHAR, ALPHAI, BETA, LSCALE, RSCALE, RCONDE,
RCONDV, WORK
REAL(8), DIMENSION(:,:) :: A, B, VL, VR




C INTERFACE
#include <sunperf.h>

void dggevx(char balanc, char jobvl, char jobvr,  char  sense,  int  n,
double  *a, int lda, double *b, int ldb, double *alphar, dou-
ble *alphai, double *beta, double *vl, int ldvl, double  *vr,
int ldvr, int *ilo, int *ihi, double *lscale, double *rscale,
double *abnrm, double *bbnrm, double *rconde, double *rcondv,
int *info);

void dggevx_64(char balanc, char jobvl, char jobvr, char sense, long n,
double *a, long lda, double *b,  long  ldb,  double  *alphar,
double  *alphai,  double *beta, double *vl, long ldvl, double
*vr, long ldvr, long *ilo, long *ihi, double *lscale,  double
*rscale, double *abnrm, double *bbnrm, double *rconde, double
*rcondv, long *info);

Description

Oracle Solaris Studio Performance Library                           dggevx(3P)



NAME
       dggevx - compute for a pair of N-by-N real nonsymmetric matrices (A,B),
       and, optionally, the left and/or right generalized eigenvectors


SYNOPSIS
       SUBROUTINE DGGEVX(BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, B, LDB,
             ALPHAR, ALPHAI, BETA, VL, LDVL, VR, LDVR, ILO, IHI, LSCALE,
             RSCALE, ABNRM, BBNRM, RCONDE, RCONDV, WORK, LWORK, IWORK, BWORK,
             INFO)

       CHARACTER*1 BALANC, JOBVL, JOBVR, SENSE
       INTEGER N, LDA, LDB, LDVL, LDVR, ILO, IHI, LWORK, INFO
       INTEGER IWORK(*)
       LOGICAL BWORK(*)
       DOUBLE PRECISION ABNRM, BBNRM
       DOUBLE PRECISION A(LDA,*),  B(LDB,*),  ALPHAR(*),  ALPHAI(*),  BETA(*),
       VL(LDVL,*),  VR(LDVR,*),  LSCALE(*),  RSCALE(*),  RCONDE(*), RCONDV(*),
       WORK(*)

       SUBROUTINE DGGEVX_64(BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, B, LDB,
             ALPHAR, ALPHAI, BETA, VL, LDVL, VR, LDVR, ILO, IHI, LSCALE,
             RSCALE, ABNRM, BBNRM, RCONDE, RCONDV, WORK, LWORK, IWORK, BWORK,
             INFO)

       CHARACTER*1 BALANC, JOBVL, JOBVR, SENSE
       INTEGER*8 N, LDA, LDB, LDVL, LDVR, ILO, IHI, LWORK, INFO
       INTEGER*8 IWORK(*)
       LOGICAL*8 BWORK(*)
       DOUBLE PRECISION ABNRM, BBNRM
       DOUBLE PRECISION A(LDA,*),  B(LDB,*),  ALPHAR(*),  ALPHAI(*),  BETA(*),
       VL(LDVL,*),  VR(LDVR,*),  LSCALE(*),  RSCALE(*),  RCONDE(*), RCONDV(*),
       WORK(*)




   F95 INTERFACE
       SUBROUTINE GGEVX(BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, B, LDB,
              ALPHAR, ALPHAI, BETA, VL, LDVL, VR, LDVR, ILO, IHI, LSCALE,
              RSCALE, ABNRM, BBNRM, RCONDE, RCONDV, WORK, LWORK, IWORK,
              BWORK, INFO)

       CHARACTER(LEN=1) :: BALANC, JOBVL, JOBVR, SENSE
       INTEGER :: N, LDA, LDB, LDVL, LDVR, ILO, IHI, LWORK, INFO
       INTEGER, DIMENSION(:) :: IWORK
       LOGICAL, DIMENSION(:) :: BWORK
       REAL(8) :: ABNRM, BBNRM
       REAL(8), DIMENSION(:) :: ALPHAR, ALPHAI, BETA, LSCALE, RSCALE,  RCONDE,
       RCONDV, WORK
       REAL(8), DIMENSION(:,:) :: A, B, VL, VR

       SUBROUTINE GGEVX_64(BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, B,
              LDB, ALPHAR, ALPHAI, BETA, VL, LDVL, VR, LDVR, ILO, IHI,
              LSCALE, RSCALE, ABNRM, BBNRM, RCONDE, RCONDV, WORK, LWORK,
              IWORK, BWORK, INFO)

       CHARACTER(LEN=1) :: BALANC, JOBVL, JOBVR, SENSE
       INTEGER(8) :: N, LDA, LDB, LDVL, LDVR, ILO, IHI, LWORK, INFO
       INTEGER(8), DIMENSION(:) :: IWORK
       LOGICAL(8), DIMENSION(:) :: BWORK
       REAL(8) :: ABNRM, BBNRM
       REAL(8),  DIMENSION(:) :: ALPHAR, ALPHAI, BETA, LSCALE, RSCALE, RCONDE,
       RCONDV, WORK
       REAL(8), DIMENSION(:,:) :: A, B, VL, VR




   C INTERFACE
       #include <sunperf.h>

       void dggevx(char balanc, char jobvl, char jobvr,  char  sense,  int  n,
                 double  *a, int lda, double *b, int ldb, double *alphar, dou-
                 ble *alphai, double *beta, double *vl, int ldvl, double  *vr,
                 int ldvr, int *ilo, int *ihi, double *lscale, double *rscale,
                 double *abnrm, double *bbnrm, double *rconde, double *rcondv,
                 int *info);

       void dggevx_64(char balanc, char jobvl, char jobvr, char sense, long n,
                 double *a, long lda, double *b,  long  ldb,  double  *alphar,
                 double  *alphai,  double *beta, double *vl, long ldvl, double
                 *vr, long ldvr, long *ilo, long *ihi, double *lscale,  double
                 *rscale, double *abnrm, double *bbnrm, double *rconde, double
                 *rcondv, long *info);



PURPOSE
       dggevx computes for a pair of N-by-N real nonsymmetric  matrices  (A,B)
       the generalized eigenvalues, and optionally, the left and/or right gen-
       eralized eigenvectors.

       Optionally also, it computes a balancing transformation to improve  the
       conditioning  of  the  eigenvalues  and eigenvectors (ILO, IHI, LSCALE,
       RSCALE, ABNRM, and BBNRM), reciprocal condition numbers for the  eigen-
       values  (RCONDE), and reciprocal condition numbers for the right eigen-
       vectors (RCONDV).

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

       The right eigenvector v(j) corresponding to the eigenvalue lambda(j) of
       (A,B) satisfies

                        A * v(j) = lambda(j) * B * v(j) .

       The left eigenvector u(j) corresponding to the eigenvalue lambda(j)  of
       (A,B) satisfies

                        u(j)**H * A  = lambda(j) * u(j)**H * B.

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


ARGUMENTS
       BALANC (input)
                 Specifies the balance option to be performed.  = 'N':  do not
                 diagonally scale or permute;
                 = 'P':  permute only;
                 = 'S':  scale only;
                 = 'B':  both permute and scale.  Computed  reciprocal  condi-
                 tion  numbers will be for the matrices after permuting and/or
                 balancing. Permuting does not change  condition  numbers  (in
                 exact arithmetic), but balancing does.


       JOBVL (input)
                 = 'N':  do not compute the left generalized eigenvectors;
                 = 'V':  compute the left generalized eigenvectors.


       JOBVR (input)
                 = 'N':  do not compute the right generalized eigenvectors;
                 = 'V':  compute the right generalized eigenvectors.


       SENSE (input)
                 Determines  which  reciprocal condition numbers are computed.
                 = 'N': none are computed;
                 = 'E': computed for eigenvalues only;
                 = 'V': computed for eigenvectors only;
                 = 'B': computed for eigenvalues and eigenvectors.


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


       A (input/output)
                 On entry, the matrix A in the pair (A,B).   On  exit,  A  has
                 been  overwritten.  If JOBVL='V' or JOBVR='V' or both, then A
                 contains the first part of the real Schur form of  the  "bal-
                 anced" versions of the input A and B.


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


       B (input/output)
                 On  entry,  the  matrix  B in the pair (A,B).  On exit, B has
                 been overwritten. If JOBVL='V' or JOBVR='V' or both,  then  B
                 contains  the second part of the real Schur form of the "bal-
                 anced" versions of the input A and B.


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


       ALPHAR (output)
                 On exit, (ALPHAR(j) + ALPHAI(j)*i)/BETA(j),  j=1,...,N,  will
                 be  the  generalized eigenvalues.  If ALPHAI(j) is zero, then
                 the j-th eigenvalue is real; if positive, then the  j-th  and
                 (j+1)-st  eigenvalues  are  a  complex  conjugate  pair, with
                 ALPHAI(j+1) negative.

                 Note: the quotients ALPHAR(j)/BETA(j)  and  ALPHAI(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,  ALPHAR  and ALPHAI will be always less
                 than and usually comparable with norm(A)  in  magnitude,  and
                 BETA always less than and usually comparable with norm(B).


       ALPHAI (output)
                 See the description of ALPHAR.


       BETA (output)
                 See the description of ALPHAR.


       VL (output)
                 If  JOBVL  =  'V',  the left eigenvectors u(j) are stored one
                 after another in the columns of VL,  in  the  same  order  as
                 their  eigenvalues. If the j-th eigenvalue is real, then u(j)
                 = VL(:,j), the j-th column of VL. If the  j-th  and  (j+1)-th
                 eigenvalues  form  a  complex  conjugate  pair,  then  u(j) =
                 VL(:,j)+i*VL(:,j+1) and u(j+1) =  VL(:,j)-i*VL(:,j+1).   Each
                 eigenvector  will  be  scaled  so  the largest component have
                 abs(real part) + abs(imag. part)  =  1.   Not  referenced  if
                 JOBVL = 'N'.


       LDVL (input)
                 The  leading  dimension  of  the matrix VL. LDVL >= 1, and if
                 JOBVL = 'V', LDVL >= N.


       VR (output)
                 If JOBVR = 'V', the right eigenvectors v(j)  are  stored  one
                 after  another  in  the  columns  of VR, in the same order as
                 their eigenvalues. If the j-th eigenvalue is real, then  v(j)
                 =  VR(:,j),  the  j-th column of VR. If the j-th and (j+1)-th
                 eigenvalues form  a  complex  conjugate  pair,  then  v(j)  =
                 VR(:,j)+i*VR(:,j+1)  and  v(j+1) = VR(:,j)-i*VR(:,j+1).  Each
                 eigenvector will be scaled  so  the  largest  component  have
                 abs(real  part)  +  abs(imag.  part)  = 1.  Not referenced if
                 JOBVR = 'N'.


       LDVR (input)
                 The leading dimension of the matrix VR. LDVR  >=  1,  and  if
                 JOBVR = 'V', LDVR >= N.


       ILO (output)
                 ILO  and  IHI are integer values such that on exit A(i,j) = 0
                 and B(i,j) = 0  if  i  >  j  and  j  =  1,...,ILO-1  or  i  =
                 IHI+1,...,N.  If BALANC = 'N' or 'S', ILO = 1 and IHI = N.


       IHI (output)
                 See the description of ILO.


       LSCALE (output)
                 Details  of  the  permutations and scaling factors applied to
                 the left side of A and B.  If PL(j) is the index of  the  row
                 interchanged  with  row  j,  and  DL(j) is the scaling factor
                 applied to row j, then LSCALE(j) = PL(j)  for j = 1,...,ILO-1
                 =  DL(j)   for  j = ILO,...,IHI = PL(j)  for j = IHI+1,...,N.
                 The order in which the interchanges are made is N  to  IHI+1,
                 then 1 to ILO-1.


       RSCALE (output)
                 Details  of  the  permutations and scaling factors applied to
                 the right side of A and B.  If PR(j) is the index of the col-
                 umn interchanged with column j, and DR(j) is the scaling fac-
                 tor applied to column j, then RSCALE(j)  =  PR(j)   for  j  =
                 1,...,ILO-1  =  DR(j)   for  j = ILO,...,IHI = PR(j)  for j =
                 IHI+1,...,N The order in which the interchanges are made is N
                 to IHI+1, then 1 to ILO-1.


       ABNRM (output)
                 The one-norm of the balanced matrix A.


       BBNRM (output)
                 The one-norm of the balanced matrix B.


       RCONDE (output)
                 If  SENSE  =  'E' or 'B', the reciprocal condition numbers of
                 the selected eigenvalues, stored in consecutive  elements  of
                 the  array.  For  a complex conjugate pair of eigenvalues two
                 consecutive elements of RCONDE are set  to  the  same  value.
                 Thus  RCONDE(j), RCONDV(j), and the j-th columns of VL and VR
                 all correspond to the same eigenpair (but not in general  the
                 j-th  eigenpair,  unless  all  eigenpairs  are selected).  If
                 SENSE = 'V', RCONDE is not referenced.


       RCONDV (output)
                 If SENSE = 'V' or 'B',  the  estimated  reciprocal  condition
                 numbers  of  the selected eigenvectors, stored in consecutive
                 elements of the array. For a complex eigenvector two consecu-
                 tive elements of RCONDV are set to the same value. If the ei-
                 genvalues cannot be reordered to compute RCONDV(j), RCONDV(j)
                 is set to 0; this can only occur when the true value would be
                 very small anyway.  If SENSE = 'E', RCONDV is not referenced.


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


       LWORK (input)
                 The  dimension  of  the  array WORK. LWORK >= max(1,6*N).  If
                 SENSE = 'E', LWORK >= 12*N.  If SENSE = 'V' or 'B', LWORK  >=
                 2*N*N+12*N+16.

                 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.


       IWORK (workspace)
                 dimension(N+6) If SENSE = 'E', IWORK is not referenced.


       BWORK (workspace)
                 dimension(N) If SENSE = 'N', BWORK is not referenced.


       INFO (output)
                 = 0:  successful exit
                 < 0:  if INFO = -i, the i-th argument had an illegal value.
                 = 1,...,N: The QZ iteration  failed.   No  eigenvectors  have
                 been calculated, but ALPHAR(j), ALPHAI(j), and BETA(j) should
                 be correct for j=INFO+1,...,N.  > N:   =N+1:  other  than  QZ
                 iteration failed in DHGEQZ.
                 =N+2: error return from DTGEVC.

FURTHER DETAILS
       Balancing  a matrix pair (A,B) includes, first, permuting rows and col-
       umns to  isolate  eigenvalues,  second,  applying  diagonal  similarity
       transformation  to the rows and columns to make the rows and columns as
       close in norm as possible. The computed  reciprocal  condition  numbers
       correspond  to the balanced matrix. Permuting rows and columns will not
       change the condition numbers (in exact arithmetic) but diagonal scaling
       will.   For  further  explanation of balancing, see section 4.11.1.2 of
       LAPACK Users' Guide.

       An approximate error bound on the chordal  distance  between  the  i-th
       computed generalized eigenvalue w and the corresponding exact eigenval-
       ue lambda is
       hord(w, lambda) <= EPS * norm(ABNRM, BBNRM) / RCONDE(I)

       An approximate error bound for the  angle  between  the  i-th  computed
       eigenvector VL(i) or VR(i) is given by
       PS * norm(ABNRM, BBNRM) / DIF(i).

       For  further explanation of the reciprocal condition numbers RCONDE and
       RCONDV, see section 4.11 of LAPACK User's Guide.




                                  7 Nov 2015                        dggevx(3P)