Visible to Intel only — GUID: GUID-B92B2BB7-9440-4157-98B8-85621B149A0F
Visible to Intel only — GUID: GUID-B92B2BB7-9440-4157-98B8-85621B149A0F
potrs_batch (Group Version)
Solves a batch of systems of linear equations with a Cholesky-factored symmetric (Hermitian) positive-definite coefficient matrices. This routine belongs to the oneapi::mkl::lapack namespace.
Description
The routine solves for Xi the systems of linear equations Ai*Xi = Bi with a symmetric positive-definite or, for complex data, Hermitian positive-definite matrices Ai, given the Cholesky factorization of Ai, i ϵ{1...batch_size} :
Ai = UiT*Ui for real data, Ai = UiH*Ui for complex data if uplog=mkl::uplo::upper
Ai = Li*LiT for real data, Ai = Li*LiH for complex data if uplog=mkl::uplo::lower
where Li is a lower triangular matrix and Ui is upper triangular, 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.
The systems are solved with multiple right-hand sides stored in the columns of the matrices Bi.
Before calling this routine, matrices Ai should be factorized by a call to potrf_batch (Group Version).
API
Syntax
namespace oneapi::mkl::lapack { sycl::event potrs_batch(sycl::queue &queue, mkl::uplo *uplo, std::int64_t *n, std::int64_t *nrhs, T **a, std::int64_t *lda, T **b, std::int64_t *ldb, 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 (ng ≥ 0).
- nrhs
-
Array of group_count parameters nrhsg parameters.
Each nrhsg specifies the number of right-hand sides supplied for group g in corresponding part of array b (nrhsg ≥ 0).
- a
-
Array of batch_size pointers to Cholesky factored matrices Ai as returned by potrf_batch (Group Version).
- lda
-
Array of group_count parameters ldag.
Each ldag specifies the leading dimension of matrices from a belonging to group g (ldag ≥ max(1, ng)).
- b
-
Array of batch_size pointers to right-hand side matrices Bi, each of size ldbg*nrhsg, where g is an index of group corresponding to Bi.
- ldb
-
Array of group_count parameters ldbg.
Each ldbg specifies the leading dimension of matrices from b belonging to group g (ldbg ≥ max(1, ng)).
- 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 potrs_batch_scratchpad_size (Group Version).
- events
-
List of events to wait for before starting computation. Defaults to empty list.
Output Parameters
- b
-
The matrices pointed to by array b are overwritten by the solution matrices Xi.
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.