.. _potrf_batch-group-version: potrf_batch (Group Version) =========================== Computes the Cholesky factorizations of a batch of symmetric (or Hermitian, for complex data) positive-definite matrices. This routine belongs to the ``oneapi::mkl::lapack`` namespace. .. contents:: :local: :depth: 1 Description *********** The routine forms the Cholesky factorizations of a symmetric positive-definite or, for complex data, Hermitian positive-definite matrices ``A``\ :sub:`i`, ``i``\ ``ϵ{1...batch_size}``: - ``A``\ :sub:`i` = ``U``\ :sub:`i`\ :sup:`T` \* ``U``\ :sub:`i` for real data, ``A``\ :sub:`i` = ``U``\ :sub:`i`\ :sup:`H` \* ``U``\ :sub:`i` for complex data. if ``uplog = mkl::uplo::upper``, - ``A``\ :sub:`i` = ``L``\ :sub:`i`\ :sup:`T` \* ``L``\ :sub:`i` for real data, ``A``\ :sub:`i` = ``L``\ :sub:`i`\ :sup:`H` \* ``L``\ :sub:`i` for complex data if ``uplog = mkl::uplo::lower`` Where ``L``\ :sub:`i` is a lower triangular matrix and ``U``\ :sub:`i` is an upper triangular matrix, ``g`` is an index of group of parameters corresponding to ``A``\ :sub:`i`, and the total number of problems to solve, ``batch_size``, is a sum of sizes for all of the groups of parameters as provided by the ``group_sizes`` array. API *** Syntax ------ .. code-block:: cpp namespace oneapi::mkl::lapack { cl::sycl::event potrf_batch(cl::sycl::queue &queue, mkl::uplo *uplo, std::int64_t *n, T **a, std::int64_t *lda, std::int64_t group_count, std::int64_t *group_sizes, T *scratchpad, std::int64_t scratchpad_size, const std::vector &events = {}) } Function supports the following precisions and devices. .. list-table:: :header-rows: 1 * - T - Devices supported * - ``float`` - Host, CPU, and GPU * - ``double`` - Host, CPU, and GPU * - ``std::complex`` - Host, CPU, and GPU * - ``std::complex`` - Host, CPU, and GPU Input Parameters ---------------- queue Device queue where calculations will be performed. uplo Array of ``group_count``\ uplo\ :sub:`g` parameters. Each of uplo\ :sub:`g` indicates whether the upper or lower triangular parts of the input matrices are provided: If uplo\ :sub:`g`\ =\ ``mkl::uplo::upper``, input matrices from array ``a`` belonging to group ``g`` store the upper triangular parts. If uplo\ :sub:`g`\ =\ ``mkl::uplo::lower``, input matrices from array ``a`` belonging to group ``g`` store the lower triangular parts. n Array of ``group_count`` parameters ``n``\ :sub:`g`. Each ``n``\ :sub:`g` specifies the order of the input matrices from array ``a`` belonging to group ``g``. a Array of ``batch_size`` pointers to input matrices ``A``\ :sub:`i`, each being of size ``lda``\ :sub:`g`\ \*\ ``n``\ :sub:`g` (``g`` is an index of the group to which ``A``\ :sub:`i` belongs) and holding either upper or lower triangular part as specified by uplo\ :sub:`g`. lda Array of ``group_count`` parameters lda\ :sub:`g`. Each lda\ :sub:`g` specifies the leading dimension of matrices from ``a`` belonging to group ``g``. group_count Specifies the number of groups of parameters. Must be at least 0. group_sizes Array of group_count integers. Array element with index ``g`` specifies the number of problems to solve for each of the groups of parameters ``g``. So the total number of problems to solve, ``batch_size``, is a sum of all parameter group sizes. scratchpad Scratchpad memory to be used by 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 then the value returned by :ref:`potrf_batch_scratchpad_size-group-version`. events List of events to wait for before starting computation. Defaults to empty list. Output Parameters ----------------- a The matrices pointed to by array a are overwritten by the Cholesky factors ``U``\ :sub:`i` or ``L``\ :sub:`i`, as specified by uplo\ :sub:`g` from the corresponding group of parameters. Exceptions ---------- .. tabularcolumns:: |\Y{0.3}|\Y{0.7}| .. list-table:: :header-rows: 1 * - Exception - Description * - ``mkl::lapack::batch_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 = -n``, the ``n``-th parameter had an illegal value. If ``info`` equals the value passed as scratchpad size, and detail() returns non-zero, then the passed scratchpad is of insufficient size, and the required size should be not less then value returned by the detail() method of the exception object. If ``info`` is zero, then the diagonal element of some of ``U``\ :sub:`i` is zero, and the solve could not be completed. The indexes of such matrices in the batch can be obtained with the ids() method of the exception object. You can obtain the indexes of the first zero diagonal elements in these ``U``\ :sub:`i` matrices using the infos() method of the exception object. Return Values ------------- Output event to wait on to ensure computation is complete.