Visible to Intel only — GUID: GUID-3210E054-C515-4EF4-A7AC-15AA35DA0D30
Visible to Intel only — GUID: GUID-3210E054-C515-4EF4-A7AC-15AA35DA0D30
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.
Description
The routine forms the Cholesky factorizations of a symmetric positive-definite or, for complex data, Hermitian positive-definite matrices Ai, iϵ{1...batch_size}:
Ai = UiT * Ui for real data, Ai = UiH * Ui for complex data. if uplog = mkl::uplo::upper,
Ai = LiT * Li for real data, Ai = LiH * Li for complex data if uplog = mkl::uplo::lower
Where Li is a lower triangular matrix and Ui is an upper triangular matrix, g is an index of group of parameters corresponding to Ai, 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
namespace oneapi::mkl::lapack { sycl::event potrf_batch(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<sycl::event> &events = {}) }
This function supports the following precisions and devices:
T |
Devices supported |
---|---|
float |
CPU and GPU |
double |
CPU and GPU |
std::complex<float> |
CPU and GPU |
std::complex<double> |
CPU and GPU |
Input Parameters
- queue
-
Device queue where calculations will be performed.
- uplo
-
Array of group_countuplog parameters.
Each of uplog indicates whether the upper or lower triangular parts of the input matrices are provided:
If uplog=mkl::uplo::upper, input matrices from array a belonging to group g store the upper triangular parts.
If uplog=mkl::uplo::lower, input matrices from array a belonging to group g store the lower triangular parts.
- n
-
Array of group_count parameters ng.
Each ng specifies the order of the input matrices from array a belonging to group g.
- a
-
Array of batch_size pointers to input matrices Ai, each being of size ldag*ng (g is an index of the group to which Ai belongs) and holding either upper or lower triangular part as specified by uplog.
- lda
-
Array of group_count parameters ldag.
Each ldag 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 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 Ui or Li, as specified by uplog from the corresponding group of parameters.
Exceptions
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 Ui 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 Ui matrices using the infos() method of the exception object. |
Return Values
Output event to wait on to ensure computation is complete.