GVCRG

   more...
Computes all of the eigenvalues and eigenvectors of a generalized real eigensystem Az = λBz.
Required Arguments
A — Real matrix of order N. (Input)
B — Real matrix of order N. (Input)
ALPHA — Complex vector of size N containing scalars αi. If βi  0, λi = αi / βii = 1, , n are the eigenvalues of the system.
BETAV — Vector of size N containing scalars βi. (Output)
EVEC — Complex matrix of order N. (Output)
The J-th eigenvector, corresponding to λJ, is stored in the J-th column. Each vector is normalized to have Euclidean length equal to the value one.
Optional Arguments
N — Order of the matrices A and B. (Input)
Default: N = SIZE (A,2).
LDA — Leading dimension of A exactly as specified in the dimension statement in the calling program. (Input)
Default: LDA = SIZE (A,1).
LDB — Leading dimension of B exactly as specified in the dimension statement in the calling program. (Input)
Default: LDB = SIZE (B,1).
LDEVEC — Leading dimension of EVEC exactly as specified in the dimension statement in the calling program. (Input)
Default: LDEVEC = SIZE (EVEC,1).
FORTRAN 90 Interface
Generic: CALL GVCRG (A, B, ALPHA, BETAV, EVEC [,])
Specific: The specific interface names are S_GVCRG and D_GVCRG.
FORTRAN 77 Interface
Single: CALL GVCRG (N, A, LDA, B, LDB, ALPHA, BETAV, EVEC, LDEVEC)
Double: The double precision name is DGVCRG.
Description
Routine GVCRG computes the complex eigenvalues and eigenvectors of the generalized eigensystem Ax = λBx where A and B are real matrices of order N. The eigenvalues for this problem can be infinite; so instead of returning λ, GVCRG returns complex numbers α and real numbers β. If β is nonzero, then λ = α/β. For problems with small ∣β∣ users can choose to solve the mathematically equivalent problem Bx = μAx where μ = λ-1.
The first step of the QZ algorithm is to simultaneously reduce A to upper Hessenberg form and B to upper triangular form. Then, orthogonal transformations are used to reduce A to quasi-upper-triangular form while keeping B upper triangular. The generalized eigenvalues and eigenvectors for the reduced problem are then computed.
The underlying code is based on either EISPACK or LAPACK code depending upon which supporting libraries are used during linking. For a detailed explanation, see “Using ScaLAPACK, LAPACK, LINPACK, and EISPACK” in the Introduction section of this manual.
Comments
1. Workspace may be explicitly provided, if desired, by use of G8CRG/DG8CRG. The reference is:
CALL G8CRG (N, A, LDA, B, LDB, ALPHA, BETAV, EVEC, LDEVEC, ACOPY, BCOPY, ECOPY, RWK, CWK, IWK)
The additional arguments are as follows:
ACOPY — Work array of size N2. The arrays A and ACOPY may be the same, in which case the first N2 elements of A will be destroyed.
BCOPY — Work array of size N2. The arrays B and BCOPY may be the same, in which case the first N2 elements of B will be destroyed.
ECOPY — Work array of size N2.
RWK — Work array of size N.
CWK — Complex work array of size N.
IWK — Integer work array of size N.
2. Integer Options with Chapter 11 Options Manager
1 This option uses eight values to solve memory bank conflict (access inefficiency) problems. In routine G8CRG, the internal or working leading dimensions of ACOPY and ECOPY are both increased by IVAL(3) when N is a multiple of IVAL(4). The values IVAL(3) and IVAL(4) are temporarily replaced by IVAL(1) and IVAL(2), respectively, in routine GVCRG. Analogous comments hold for the array BCOPY and the option values IVAL(5) – IVAL(8). Additional memory allocation and option value restoration are automatically done in GVCRG. There is no requirement that users change existing applications that use GVCRG or G8CRG. Default values for the option are IVAL(*) = 1, 16, 0, 1, 1, 16, 0, 1. Items 5–8 in IVAL(*) are for the generalized eigenvalue problem and are not used in GVCRG.
Example
In this example, DATA statements are used to set A and B. The eigenvalues, eigenvectors and performance index are computed and printed for the systems Ax = λBx and Bx = μAx where μ = λ-1. For more details about the performance index, see routine GPIRG.
 
USE IMSL_LIBRARIES
 
IMPLICIT NONE
INTEGER LDA, LDB, LDEVEC, N
PARAMETER (N=3, LDA=N, LDB=N, LDEVEC=N)
!
INTEGER I, NOUT
REAL A(LDA,N), B(LDB,N), BETAV(N), PI
COMPLEX ALPHA(N), EVAL(N), EVEC(LDEVEC,N)
!
! Define values of A and B:
! A = ( 1.0 0.5 0.0 )
! (-10.0 2.0 0.0 )
! ( 5.0 1.0 0.5 )
!
! B = ( 0.5 0.0 0.0 )
! ( 3.0 3.0 0.0 )
! ( 4.0 0.5 1.0 )
!
! Declare variables
DATA A/1.0, -10.0, 5.0, 0.5, 2.0, 1.0, 0.0, 0.0, 0.5/
DATA B/0.5, 3.0, 4.0, 0.0, 3.0, 0.5, 0.0, 0.0, 1.0/
!
CALL GVCRG (A, B, ALPHA, BETAV, EVEC)
! Compute eigenvalues
DO 10 I=1, N
EVAL(I) = ALPHA(I)/BETAV(I)
10 CONTINUE
! Compute performance index
PI = GPIRG(N,A,B,ALPHA,BETAV,EVEC)
! Print results
CALL UMACH (2, NOUT)
CALL WRCRN ('EVAL', EVAL, 1, N, 1)
CALL WRCRN ('EVEC', EVEC)
WRITE (NOUT,'(/,A,F6.3)') ' Performance index = ', PI
! Solve for reciprocals of values
CALL GVCRG (B, A, ALPHA, BETAV, EVEC)
 
! Compute reciprocals
DO 20 I=1, N
EVAL(I) = ALPHA(I)/BETAV(I)
20 CONTINUE
! Compute performance index
PI = GPIRG(N,B,A,ALPHA,BETAV,EVEC)
! Print results
CALL WRCRN ('EVAL reciprocals', EVAL, 1, N, 1)
CALL WRCRN ('EVEC', EVEC)
WRITE (NOUT,'(/,A,F6.3)') ' Performance index = ', PI
END
Output
EVAL
1 2 3
( 0.833, 1.993) ( 0.833,-1.993) ( 0.500, 0.000)
EVEC
1 2 3
1 (-0.197, 0.150) (-0.197,-0.150) (-0.000, 0.000)
2 (-0.069,-0.568) (-0.069, 0.568) (-0.000, 0.000)
3 ( 0.782, 0.000) ( 0.782, 0.000) ( 1.000, 0.000)
 
Performance index = 0.384
EVAL reciprocals
1 2 3
( 2.000, 0.000) ( 0.179, 0.427) ( 0.179,-0.427)
EVEC
1 2 3
1 ( 0.000, 0.000) (-0.197,-0.150) (-0.197, 0.150)
2 ( 0.000, 0.000) (-0.069, 0.568) (-0.069,-0.568)
3 ( 1.000, 0.000) ( 0.782, 0.000) ( 0.782, 0.000)
 
Performance index = 0.283