?hpgv
?hpgv
Computes all eigenvalues and, optionally, eigenvectors of a complex generalized Hermitian positive-definite eigenproblem with matrices in packed storage.
Syntax
lapack_int LAPACKE_chpgv
(
int
matrix_layout
,
lapack_int
itype
,
char
jobz
,
char
uplo
,
lapack_int
n
,
lapack_complex_float*
ap
,
lapack_complex_float*
bp
,
float*
w
,
lapack_complex_float*
z
,
lapack_int
ldz
);
lapack_int LAPACKE_zhpgv
(
int
matrix_layout
,
lapack_int
itype
,
char
jobz
,
char
uplo
,
lapack_int
n
,
lapack_complex_double*
ap
,
lapack_complex_double*
bp
,
double*
w
,
lapack_complex_double*
z
,
lapack_int
ldz
);
Include Files
- mkl.h
Description
The routine computes all the 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.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'
- 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).
- ldz
- The leading dimension of the output arrayz;. Ifldz≥1,jobz='V'.ldz≥max(1,n)
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.
- w
- Array, size at least max(1,n).If, contains the eigenvalues in ascending order.info= 0
- z
- Arrayz(size max(1,.ldz*n))If, then ifjobz='V',info= 0zcontains the matrixZof eigenvectors. The eigenvectors are normalized as follows:iforitype= 12,;Z*HB*Z= Iif,itype= 3;Z*inv(HB)*Z= IIf, thenjobz='N'zis 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 chpev
/zhpev
returned an error code:If ,
info
= i
≤
n
chpev
/zhpev
failed to converge, and i
off-diagonal elements of an intermediate tridiagonal did not converge to zero;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.