.. _hegvd-usm-version: 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. .. contents:: :local: :depth: 1 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. API *** Syntax ------ .. code-block:: cpp namespace oneapi::mkl::lapack { 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 std::vector &events = {}) } ``hegvd`` (USM version) supports the following precision and devices. .. list-table:: :header-rows: 1 * - T - Devices Supported * - ``std::complex`` - Host, CPU, GPU * - ``std::complex`` - Host, CPU, GPU 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 is ``A*x = lambda*B*x;`` if itype\ ``= 2``, the problem type is ``A*B*x = lambda*x;`` if itype\ ``= 3``, the problem type is ``B*A*x = lambda*x``. jobz Must be ``job::novec`` or ``job::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`` or ``uplo::lower``. If ``uplo = uplo::upper``, a and b store the upper triangular part of ``A`` and ``B``. If ``uplo = uplo::lower``, a and b stores the lower triangular part of ``A`` and ``B``. n The order of the matrices ``A`` and ``B`` (``0≤n``). a Pointer to the array of size ``a(lda,*)`` containing the upper or lower triangle of the Hermitian matrix ``A``, 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 matrix ``B``, 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 :ref:`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 if ``info = 0``, a contains the matrix ``Z`` of eigenvectors. The eigenvectors are normalized as follows: if itype\ ``= 1`` or ``2``, ``Z``\ :sup:`H`\ ``*B*Z = I``; if itype\ ``= 3``, ``Z``\ :sup:`H`\ \*\ ``inv(B)*Z = I``; If ``jobz = job::novec``, then on exit the upper triangle (if ``uplo = uplo::upper``) or the lower triangle (if ``uplo = uplo::lower``) of ``A``, including the diagonal, is destroyed. b On exit, if ``info≤n``, the part of b containing the matrix is overwritten by the triangular factor ``U`` or ``L`` from the Cholesky factorization ``B = U``\ :sup:`H`\ ``*U``\ or ``B`` = ``L``\ \*\ ``L``\ :sup:`H`. w Pointer to the array of size at least n. If ``info = 0``, contains the eigenvalues of the matrix ``A`` in ascending order. See also info. Exceptions ---------- .. tabularcolumns:: |\Y{0.3}|\Y{0.7}| .. list-table:: :header-rows: 1 * - Exception - Description * - ``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 ``info = -i``, the ``i``-th parameter had an illegal value. For ``info ≤ n``: If ``info = i``, and ``jobz = job::novec``, then the algorithm failed to converge; ``i`` indicates the number of off-diagonal elements of an intermediate tridiagonal form which did not converge to zero. If ``info = i``, and ``jobz = job:vec``, then the algorithm failed to compute an eigenvalue while working on the submatrix lying in rows and columns ``info/(n+1)`` through ``mod(info,n+1)``. For ``info > n``: If ``info = n + i``, for ``1 ≤ i ≤ n``, then the leading minor of order ``i`` of ``B`` is not positive-definite. The factorization of ``B`` could not be completed and no eigenvalues or eigenvectors were computed. If ``info`` is equal to the value passed as scratchpad size, and detail() returns non zero, then the passed scratchpad has an insufficient size, and the required size should not be less than the value returned by the detail() method of the exception object. Return Values ------------- Output event to wait on to ensure computation is complete.