hegvd (USM Version)¶
Computes all the eigenvalues, and optionally, the eigenvectors of a
complex generalized Hermitian positive-definite eigenproblem using a
divide and conquer method. This routine belongs to the
oneapi::mkl::lapack
namespace.
Syntax
-
cl::sycl::event
hegvd
(cl::sycl::queue &queue, std::int64_t itype, mkl::job jobz, mkl::uplo uplo, std::int64_t n, T *a, std::int64_t lda, T *b, std::int64_t ldb, T *w, T *scratchpad, std::int64_t scratchpad_size, const cl::sycl::vector_class<cl::sycl::event> &events = {})¶
hegvd
(USM version) supports the following precision and
devices.
T |
Devices Supported |
---|---|
|
Host, CPU, GPU |
|
Host, CPU, GPU |
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
, or B*A*x = λ*x
.
Here A
and B
are assumed to be Hermitian and B
is also
positive definite.
It uses a divide and conquer algorithm.
Input Parameters
- queue
Device queue where calculations will be performed.
- itype
Must be 1 or 2 or 3. Specifies the problem type to be solved:
if itype
= 1
, the problem type isA*x = lambda*B*x;
if itype
= 2
, the problem type isA*B*x = lambda*x;
if itype
= 3
, the problem type isB*A*x = lambda*x
.- jobz
Must be
job::novec
orjob::vec
.If
jobz = job::novec
, then only eigenvalues are computed.If
jobz = job::vec
, then eigenvalues and eigenvectors are computed.- uplo
Must be
uplo::upper
oruplo::lower
.If
uplo = uplo::upper
, a and b store the upper triangular part ofA
andB
.If
uplo = uplo::lower
, a and b stores the lower triangular part ofA
andB
.- n
The order of the matrices
A
andB
(0≤n
).- a
Pointer to the array of size
a(lda,*)
containing the upper or lower triangle of the Hermitian matrixA
, as specified by uplo.The second dimension of a must be at least
max(1, n)
.- lda
The leading dimension of a; at least
max(1,n)
.- b
Pointer to the array of size
b(ldb,*)
containing the upper or lower triangle of the Hermitian matrixB
, as specified by uplo.The second dimension of b must be at least
max(1, n)
.- ldb
The leading dimension of b; at least
max(1,n)
.- scratchpad
Pointer to scratchpad memory to be used by the routine for storing intermediate results.
- scratchpad_size
Size of scratchpad memory as a number of floating point elements of type
T
. Size should not be less than the value returned by the hegvd_scratchpad_size function.- events
List of events to wait for before starting computation. Defaults to empty list.
Output Parameters
- a
On exit, if
jobz = job::vec
, then ifinfo = 0
, a contains the matrixZ
of eigenvectors. The eigenvectors are normalized as follows:if itype
= 1
or2
,Z
H*B*Z = I
;if itype
= 3
,Z
H*inv(B)*Z = I
;If
jobz = job::novec
, then on exit the upper triangle (ifuplo = uplo::upper
) or the lower triangle (ifuplo = uplo::lower
) ofA
, including the diagonal, is destroyed.- b
On exit, if
info≤n
, the part of b containing the matrix is overwritten by the triangular factorU
orL
from the Cholesky factorizationB = U
H*U
orB
=L
*L
H.- w
Pointer to the array of size at least n. If
info = 0
, contains the eigenvalues of the matrixA
in ascending order. See also info.
Exceptions
mkl::lapack::exception |
This exception is thrown when problems occur during calculations. You can obtain the info code of the problem using the info() method of the exception object: If |
Return Values
Output event to wait on to ensure computation is complete.