Visible to Intel only — GUID: GUID-58B539AE-836B-4C66-8232-EAB117174908
Visible to Intel only — GUID: GUID-58B539AE-836B-4C66-8232-EAB117174908
getrs_batch (Group Version)
Solves a batch of systems of linear equations with a batch of LU-factored square coefficient matrices, with multiple right-hand sides. This routine belongs to the oneapi::mkl::lapack namespace.
Description
The routine solves for Xi (iϵ{1...batch_size}) the following systems of linear equations:
Ai * Xi = Bi If trans = mkl::transpose::notrans
AiT * Xi = Bi If trans = mkl::transpose::trans
AiH * Xi = Bi If trans = mkl::transpose::conjtrans
Before calling this routine you must call getrf_batch (Group Version) to compute the LU factorization of Ai.
The total number of problems to solve, batch_size, is a sum of sizes of all of the groups of parameters as provided by``group_sizes`` array.
API
Syntax
namespace oneapi::mkl::lapack { sycl::event getrs_batch(sycl::queue &queue, mkl::transpose *trans, std::int64_t *n, std::int64_t *nrhs, T **a, std::int64_t *lda, std::int64_t **ipiv, 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.
- trans
-
Array of group_count parameters transg indicating the form of the equations for the group g:
If trans = mkl::transpose::nontrans, then Ai*Xi = Bi is solved for Xi.
If trans = mkl::transpose::trans, then AiT*Xi = Bi is solved for Xi.
If trans = mkl::transpose::conjtrans, then AiH*Xi = Bi is solved for Xi.
- n
-
Array of group_count parameters ng specifying the order of the matrices Ai and the number of rows in matrices Bi (0 ≤ ng) belonging to group g.
- nrhs
-
Array of group_count parameters nrhsg specifying the number of right hand sides(0 ≤ nrhs)for group g.
- a
-
Array of batch_size pointers to factorization of the matricesAi, as returned by getrf_batch (Group Version).
- lda
-
Array of group_count parameters ldag specifying the leading dimension of Ai from group g.
- ipiv
-
The ipiv array, as returned by getrf_batch (Group Version).
- b
-
The array containing @ batch_size pointers to the matrices Bi whose columns are the right-hand sides for the systems of equations.
- ldb
-
Array of group_count parameters ldbg specifying the leading dimensions of Bi in the 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 getrs_batch_scratchpad_size (Group Version).
- events
-
List of events to wait for before starting computation. Defaults to empty list.
Output Parameters
- b
-
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.