getrs_batch (Group Version)¶
Solves a batch of system sof 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.
Syntax
-
cl::sycl::event
getrs_batch
(cl::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 cl::sycl::vector_class<cl::sycl::event> &events = {})¶
Function supports the following precisions and devices.
T |
Devices supported |
---|---|
|
Host, CPU, and GPU |
|
Host, CPU, and GPU |
|
Host, CPU, and GPU |
|
Host, CPU, and GPU |
Description
The routine solves for X
i (iϵ{1...batch_size}
) the
following systems of linear equations:
A
i * X
i = B
i If
trans = mkl::transpose::notrans
A
iT * X
i = B
i If
trans = mkl::transpose::trans
A
iH * X
i = B
i If
trans = mkl::transpose::conjtrans
Before calling this routine you must call getrf_batch (Group
Version) to
compute the LU factorization of A
1.
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.
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 groupg
:If trans = mkl::transpose::nontrans, then
A
i*X
i =B
i is solved forX
i.If trans = mkl::transpose::trans, then
A
iT*X
i =B
i is solved forX
i.If trans = mkl::transpose::conjtrans, then
A
iH*X
i =B
i is solved forX
i.- n
Array of
group_count
parametersn
g specifying the order of the matricesA
i and the number of rows in matricesB
i (0 ≤n
g) belonging to groupg
.- nrhs
Array of
group_count
parameters nrhsg specifying the number of right hand sides(0≤nrhs)
for groupg
.- a
Array of
batch_size
pointers to factorization of the matricesA
i, as returned by getrf_batch (Group Version).- lda
Array of
group_count
parameters ldag specifying the leading dimension ofA
i from groupg
.- ipiv
The ipiv array, as returned by getrf_batch (Group Version).
- b
The array containing @
batch_size
pointers to the matricesB
i whose columns are the right-hand sides for the systems of equations.- ldb
Array of
group_count
parameters ldbg specifying the leading dimensions ofB
i in the groupg
.- 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 parametersg
. 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.
- events
List of events to wait for before starting computation. Defaults to empty list.
Output Parameters
- b
Overwritten by the solution matrices
X
i.
Exceptions
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 |
Return Values
Output event to wait on to ensure computation is complete.