dspgvx - compute selected eigenvalues, and optionally, eigenvectors of a real generalized symmetric-definite eigenproblem, of the form A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x
SUBROUTINE DSPGVX(ITYPE, JOBZ, RANGE, UPLO, N, AP, BP, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, IWORK, IFAIL, INFO) CHARACTER*1 JOBZ, RANGE, UPLO INTEGER ITYPE, N, IL, IU, M, LDZ, INFO INTEGER IWORK(*), IFAIL(*) DOUBLE PRECISION VL, VU, ABSTOL DOUBLE PRECISION AP(*), BP(*), W(*), Z(LDZ,*), WORK(*) SUBROUTINE DSPGVX_64(ITYPE, JOBZ, RANGE, UPLO, N, AP, BP, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, IWORK, IFAIL, INFO) CHARACTER*1 JOBZ, RANGE, UPLO INTEGER*8 ITYPE, N, IL, IU, M, LDZ, INFO INTEGER*8 IWORK(*), IFAIL(*) DOUBLE PRECISION VL, VU, ABSTOL DOUBLE PRECISION AP(*), BP(*), W(*), Z(LDZ,*), WORK(*) F95 INTERFACE SUBROUTINE SPGVX(ITYPE, JOBZ, RANGE, UPLO, N, AP, BP, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, IWORK, IFAIL, INFO) CHARACTER(LEN=1) :: JOBZ, RANGE, UPLO INTEGER :: ITYPE, N, IL, IU, M, LDZ, INFO INTEGER, DIMENSION(:) :: IWORK, IFAIL REAL(8) :: VL, VU, ABSTOL REAL(8), DIMENSION(:) :: AP, BP, W, WORK REAL(8), DIMENSION(:,:) :: Z SUBROUTINE SPGVX_64(ITYPE, JOBZ, RANGE, UPLO, N, AP, BP, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, IWORK, IFAIL, INFO) CHARACTER(LEN=1) :: JOBZ, RANGE, UPLO INTEGER(8) :: ITYPE, N, IL, IU, M, LDZ, INFO INTEGER(8), DIMENSION(:) :: IWORK, IFAIL REAL(8) :: VL, VU, ABSTOL REAL(8), DIMENSION(:) :: AP, BP, W, WORK REAL(8), DIMENSION(:,:) :: Z C INTERFACE #include <sunperf.h> void dspgvx(int itype, char jobz, char range, char uplo, int n, double *ap, double *bp, double vl, double vu, int il, int iu, double abstol, int *m, double *w, double *z, int ldz, int *ifail, int *info); void dspgvx_64(long itype, char jobz, char range, char uplo, long n, double *ap, double *bp, double vl, double vu, long il, long iu, double abstol, long *m, double *w, double *z, long ldz, long *ifail, long *info);
Oracle Solaris Studio Performance Library dspgvx(3P)
NAME
dspgvx - compute selected eigenvalues, and optionally, eigenvectors of
a real generalized symmetric-definite eigenproblem, of the form
A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x
SYNOPSIS
SUBROUTINE DSPGVX(ITYPE, JOBZ, RANGE, UPLO, N, AP, BP, VL, VU, IL,
IU, ABSTOL, M, W, Z, LDZ, WORK, IWORK, IFAIL, INFO)
CHARACTER*1 JOBZ, RANGE, UPLO
INTEGER ITYPE, N, IL, IU, M, LDZ, INFO
INTEGER IWORK(*), IFAIL(*)
DOUBLE PRECISION VL, VU, ABSTOL
DOUBLE PRECISION AP(*), BP(*), W(*), Z(LDZ,*), WORK(*)
SUBROUTINE DSPGVX_64(ITYPE, JOBZ, RANGE, UPLO, N, AP, BP, VL, VU, IL,
IU, ABSTOL, M, W, Z, LDZ, WORK, IWORK, IFAIL, INFO)
CHARACTER*1 JOBZ, RANGE, UPLO
INTEGER*8 ITYPE, N, IL, IU, M, LDZ, INFO
INTEGER*8 IWORK(*), IFAIL(*)
DOUBLE PRECISION VL, VU, ABSTOL
DOUBLE PRECISION AP(*), BP(*), W(*), Z(LDZ,*), WORK(*)
F95 INTERFACE
SUBROUTINE SPGVX(ITYPE, JOBZ, RANGE, UPLO, N, AP, BP, VL, VU, IL,
IU, ABSTOL, M, W, Z, LDZ, WORK, IWORK, IFAIL, INFO)
CHARACTER(LEN=1) :: JOBZ, RANGE, UPLO
INTEGER :: ITYPE, N, IL, IU, M, LDZ, INFO
INTEGER, DIMENSION(:) :: IWORK, IFAIL
REAL(8) :: VL, VU, ABSTOL
REAL(8), DIMENSION(:) :: AP, BP, W, WORK
REAL(8), DIMENSION(:,:) :: Z
SUBROUTINE SPGVX_64(ITYPE, JOBZ, RANGE, UPLO, N, AP, BP, VL, VU,
IL, IU, ABSTOL, M, W, Z, LDZ, WORK, IWORK, IFAIL, INFO)
CHARACTER(LEN=1) :: JOBZ, RANGE, UPLO
INTEGER(8) :: ITYPE, N, IL, IU, M, LDZ, INFO
INTEGER(8), DIMENSION(:) :: IWORK, IFAIL
REAL(8) :: VL, VU, ABSTOL
REAL(8), DIMENSION(:) :: AP, BP, W, WORK
REAL(8), DIMENSION(:,:) :: Z
C INTERFACE
#include <sunperf.h>
void dspgvx(int itype, char jobz, char range, char uplo, int n, double
*ap, double *bp, double vl, double vu, int il, int iu, double
abstol, int *m, double *w, double *z, int ldz, int *ifail,
int *info);
void dspgvx_64(long itype, char jobz, char range, char uplo, long n,
double *ap, double *bp, double vl, double vu, long il, long
iu, double abstol, long *m, double *w, double *z, long ldz,
long *ifail, long *info);
PURPOSE
dspgvx computes selected eigenvalues, and optionally, eigenvectors of a
real generalized symmetric-definite eigenproblem, of the form
A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. Here A and B
are assumed to be symmetric, stored in packed storage, and B is also
positive definite. Eigenvalues and eigenvectors can be selected by
specifying either a range of values or a range of indices for the
desired eigenvalues.
ARGUMENTS
ITYPE (input)
Specifies the problem type to be solved:
= 1: A*x = (lambda)*B*x
= 2: A*B*x = (lambda)*x
= 3: B*A*x = (lambda)*x
JOBZ (input)
= 'N': Compute eigenvalues only;
= 'V': Compute eigenvalues and eigenvectors.
RANGE (input)
= 'A': all eigenvalues will be found.
= 'V': all eigenvalues in the half-open interval (VL,VU] will
be found. = 'I': the IL-th through IU-th eigenvalues will be
found.
UPLO (input)
= 'U': Upper triangle of A and B are stored;
= 'L': Lower triangle of A and B are stored.
N (input) The order of the matrix pencil (A,B). N >= 0.
AP (input/output)
Double precision array, dimension (N*(N+1)/2) On entry, the
upper or lower triangle of the symmetric matrix A, packed
columnwise in a linear array. The j-th column of A is stored
in the array AP as follows: if UPLO = 'U', AP(i + (j-1)*j/2)
= A(i,j) for 1<=i<=j; if UPLO = 'L', AP(i + (j-1)*(2*n-j)/2)
= A(i,j) for j<=i<=n.
On exit, the contents of AP are destroyed.
BP (input/output)
Double precision array, dimension (N*(N+1)/2) On entry, the
upper or lower triangle of the symmetric matrix B, packed
columnwise in a linear array. The j-th column of B is stored
in the array BP as follows: if UPLO = 'U', BP(i + (j-1)*j/2)
= B(i,j) for 1<=i<=j; if UPLO = 'L', BP(i + (j-1)*(2*n-j)/2)
= B(i,j) for j<=i<=n.
On exit, the triangular factor U or L from the Cholesky fac-
torization B = U**T*U or B = L*L**T, in the same storage for-
mat as B.
VL (input)
If RANGE='V', the lower and upper bounds of the interval to
be searched for eigenvalues. VL < VU. Not referenced if
RANGE = 'A' or 'I'.
VU (input)
See the description of VL.
IL (input)
If RANGE='I', the indices (in ascending order) of the small-
est and largest eigenvalues to be returned. 1 <= IL <= IU <=
N, if N > 0; IL = 1 and IU = 0 if N = 0. Not referenced if
RANGE = 'A' or 'V'.
IU (input)
See the description of IL.
ABSTOL (input)
The absolute error tolerance for the eigenvalues. An approx-
imate eigenvalue is accepted as converged when it is deter-
mined to lie in an interval [a,b] of width less than or equal
to
ABSTOL + EPS * max( |a|,|b| ) ,
where EPS is the machine precision. If ABSTOL is less than
or equal to zero, then EPS*|T| will be used in its place,
where |T| is the 1-norm of the tridiagonal matrix obtained by
reducing A to tridiagonal form.
Eigenvalues will be computed most accurately when ABSTOL is
set to twice the underflow threshold 2*DLAMCH('S'), not zero.
If this routine returns with INFO>0, indicating that some
eigenvectors did not converge, try setting ABSTOL to
2*DLAMCH('S').
M (output)
The total number of eigenvalues found. 0 <= M <= N. If
RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
W (output)
Double precision array, dimension (N) On normal exit, the
first M elements contain the selected eigenvalues in ascend-
ing order.
Z (output)
Double precision array, dimension (LDZ, max(1,M)) If JOBZ =
'N', then Z is not referenced. If JOBZ = 'V', then if INFO =
0, the first M columns of Z contain the orthonormal eigenvec-
tors of the matrix A corresponding to the selected eigenval-
ues, with the i-th column of Z holding the eigenvector asso-
ciated with W(i). The eigenvectors are normalized as fol-
lows: if ITYPE = 1 or 2, Z**T*B*Z = I; if ITYPE = 3,
Z**T*inv(B)*Z = I.
If an eigenvector fails to converge, then that column of Z
contains the latest approximation to the eigenvector, and the
index of the eigenvector is returned in IFAIL. Note: the
user must ensure that at least max(1,M) columns are supplied
in the array Z; if RANGE = 'V', the exact value of M is not
known in advance and an upper bound must be used.
LDZ (input)
The leading dimension of the array Z. LDZ >= 1, and if JOBZ
= 'V', LDZ >= max(1,N).
WORK (workspace)
Double precision array, dimension(8*N)
IWORK (workspace)
Integer array, dimension(5*N)
IFAIL (output)
Integer array, dimension (N) If JOBZ = 'V', then if INFO = 0,
the first M elements of IFAIL are zero. If INFO > 0, then
IFAIL contains the indices of the eigenvectors that failed to
converge. If JOBZ = 'N', then IFAIL is not referenced.
INFO (output)
= 0: successful exit
< 0: if INFO = -i, the i-th argument had an illegal value
> 0: DPPTRF or DSPEVX returned an error code:
<= N: if INFO = i, DSPEVX failed to converge; i eigenvectors
failed to converge. Their indices are stored in array IFAIL.
> N: if INFO = N + i, for 1 <= i <= N, then the leading
minor of order i of B is not positive definite. The factor-
ization of B could not be completed and no eigenvalues or
eigenvectors were computed.
FURTHER DETAILS
Based on contributions by
Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
7 Nov 2015 dspgvx(3P)