?hpgvx
?hpgvx
Computes selected eigenvalues and, optionally, eigenvectors of a generalized Hermitian positive-definite eigenproblem with matrices in packed storage.
Syntax
lapack_int LAPACKE_chpgvx
(
int
matrix_layout
,
lapack_int
itype
,
char
jobz
,
char
range
,
char
uplo
,
lapack_int
n
,
lapack_complex_float*
ap
,
lapack_complex_float*
bp
,
float
vl
,
float
vu
,
lapack_int
il
,
lapack_int
iu
,
float
abstol
,
lapack_int*
m
,
float*
w
,
lapack_complex_float*
z
,
lapack_int
ldz
,
lapack_int*
ifail
);
lapack_int LAPACKE_zhpgvx
(
int
matrix_layout
,
lapack_int
itype
,
char
jobz
,
char
range
,
char
uplo
,
lapack_int
n
,
lapack_complex_double*
ap
,
lapack_complex_double*
bp
,
double
vl
,
double
vu
,
lapack_int
il
,
lapack_int
iu
,
double
abstol
,
lapack_int*
m
,
double*
w
,
lapack_complex_double*
z
,
lapack_int
ldz
,
lapack_int*
ifail
);
Include Files
- mkl.h
Description
The routine computes selected eigenvalues, and optionally, the eigenvectors of a complex generalized Hermitian positive-definite eigenproblem, of the form
A
*x
= λ
*B
*x
A
*B
*x
= λ
*x
B
*A
*x
= λ
*x
Here
A
and B
are assumed to be Hermitian, stored in packed format, 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. Input Parameters
- matrix_layout
- Specifies whether matrix storage layout is row major (LAPACK_ROW_MAJOR) or column major (LAPACK_COL_MAJOR).
- itype
- Must be 1 or 2 or 3. Specifies the problem type to be solved:if, the problem type isitype= 1;A*x=lambda*B*xif, the problem type isitype= 2;A*B*x=lambda*xif, the problem type isitype= 3.B*A*x=lambda*x
- jobz
- Must be'N'or'V'.If, then compute eigenvalues only.jobz='N'If, then compute eigenvalues and eigenvectors.jobz='V'
- range
- Must be'A'or'V'or'I'.If, the routine computes all eigenvalues.range='A'If, the routine computes eigenvaluesrange='V'w[in the half-open interval:i]vl<w[i].≤vuIf, the routine computes eigenvalues with indicesrange='I'iltoiu.
- uplo
- Must be'U'or'L'.If, arraysuplo='U'apandbpstore the upper triangles ofAandB;If, arraysuplo='L'apandbpstore the lower triangles ofAandB.
- n
- The order of the matricesAandB().n≥0
- ap,bp
- Arrays:apcontains the packed upper or lower triangle of the Hermitian matrixA, as specified byuplo.The dimension ofapmust be at least max(1,n*(n+1)/2).bpcontains the packed upper or lower triangle of the Hermitian matrixB, as specified byuplo.The dimension ofbpmust be at least max(1,n*(n+1)/2).
- vl,vu
- If, the lower and upper bounds of the interval to be searched for eigenvalues.range='V'Constraint:.vl<vuIforrange='A''I',vlandvuare not referenced.
- il,iu
- If, the indices in ascending order of the smallest and largest eigenvalues to be returned.range='I'Constraint:1, if≤il≤iu≤n;n> 0andil=1iu=0if.n= 0Iforrange='A''V',ilandiuare not referenced.
- abstol
- The absolute error tolerance for the eigenvalues.SeeApplication Notesfor more information.
- ldz
- The leading dimension of the output arrayz;. Ifldz≥1,jobz='V'ldz≥max(1,n)for column major layout and.ldz≥max(1,m) for row major layout
Output Parameters
- ap
- On exit, the contents ofapare overwritten.
- bp
- On exit, contains the triangular factorUorLfrom the Cholesky factorizationorB=U*HU, in the same storage format asB=L*LHB.
- m
- The total number of eigenvalues found,0. If≤m≤n,range='A', and ifm=n,range='I'.m=iu-il+1
- w
- Array, size at least max(1,n).If, contains the eigenvalues in ascending order.info= 0
- z
- Arrayz(size at least max(1,.ldz*m) for column major layout and max(1,ldz*n) for row major layout)If, then ifjobz='V', the firstinfo= 0mcolumns ofzcontain the orthonormal eigenvectors of the matrixAcorresponding to the selected eigenvalues, with thei-th column ofzholding the eigenvector associated withw(i). The eigenvectors are normalized as follows:iforitype= 12,;Z*HB*Z= Iif,itype= 3;Z*inv(HB)*Z= IIf, thenjobz='N'zis not referenced.If an eigenvector fails to converge, then that column ofzcontains the latest approximation to the eigenvector, and the index of the eigenvector is returned inifail.Note: you must ensure that at least max(1,m) columns are supplied in the arrayz; if, the exact value ofrange='V'mis not known in advance and an upper bound must be used.
- ifail
- Array, size at least max(1,n).If, then ifjobz='V', the firstinfo= 0melements ofifailare zero; if, theinfo> 0ifailcontains the indices of the eigenvectors that failed to converge.If, thenjobz='N'ifailis not referenced.
Return Values
This function returns a value
info
.If , the execution is successful.
info
=0If , the
info
= -i
i
-th parameter had an illegal value.If ,
info
> 0cpptrf
/zpptrf
and chpevx
/zhpevx
returned an error code:If ,
info
= i
≤
n
chpevx
/zhpevx
failed to converge, and i
eigenvectors failed to converge. Their indices are stored in the array ifail
;If , for
info
= n
+ i
1
, then the leading minor of order ≤
i
≤
n
i
of B
is not positive-definite. The factorization of B
could not be completed and no eigenvalues or eigenvectors were computed.Application Notes
An approximate eigenvalue is accepted as converged when it is determined to lie in an interval [a,b] of width less than or equal to , where
abstol
+ε
*max(|a|,|b|)ε
is the machine precision. If is used as tolerance, where
abstol
is less than or equal to zero, then ε
*||T
||1
T
is 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*?lamch
('S'), not zero. If this routine returns with , indicating that some eigenvectors did not converge, try setting
info
> 0abstol
to 2*?lamch
('S').