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
X
i
the systems of linear
equations A
i
*X
i
= B
i
with
a symmetric positive-definite or, for complex data, Hermitian
positive-definite matrices A
i
, given the Cholesky
factorization of A
i
, i ϵ{1...batch_size}
:- Ai=UiT*Uifor real data,Ai=UiH*Uifor complex data if uplog=mkl::uplo::upper
- Ai=Li*LiTfor real data,Ai=Li*LiHfor complex data if uplog=mkl::uplo::lower
where
L
i
is a lower triangular matrix and
U
i
is upper triangular, g
is an index of group of
parameters corresponding to A
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.The systems are solved with multiple right-hand sides stored in the
columns of the matrices
B
i
.Before calling this routine, matrices
A
i
should be
factorized by a call to potrf_batch (Group Version).API
Syntax
namespace oneapi::mkl::lapack {
cl::sycl::event potrs_batch(cl::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<cl::sycl::event> &events = {})
}
Function supports the following precisions and devices.
T | Devices supported |
---|---|
float | Host, CPU, and GPU |
double | Host, CPU, and GPU |
std::complex<float> | Host, CPU, and GPU |
std::complex<double> | Host, CPU, and GPU |
Input Parameters
- queue
- Device queue where calculations will be performed.
- uplo
- Array ofgroup_countuplogparameters.Each of uplogindicates whether the upper or lower triangular parts of the input matrices are provided:If uplog=mkl::uplo::upper, input matrices from arrayabelonging to groupgstore the upper triangular parts.If uplog=mkl::uplo::lower, input matrices from arrayabelonging to groupgstore the lower triangular parts.
- n
- Array ofgroup_countparametersng.Eachngspecifies the order of the input matrices from arrayabelonging to groupg.
- nrhs
- Array ofgroup_countparameters nrhsgparameters.Each nrhsgspecifies the number of right-hand sides supplied for groupgin corresponding part of arrayb.
- a
- Array ofbatch_sizepointers to Cholesky factored matricesAias returned by potrf_batch (Group Version).
- lda
- Array ofgroup_countparameters ldag.Each ldagspecifies the leading dimension of matrices fromabelonging to groupg.
- b
- Array ofbatch_sizepointers to right-hand side matricesBi, each of size ldbg*nrhsg, wheregis an index of group corresponding toBi.
- ldb
- Array ofgroup_countparameters ldbg.Each ldbgspecifies the leading dimension of matrices frombbelonging to 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 indexgspecifies 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 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 matricesXi.
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 U 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 i matrices using the infos() method of the exception object. |
Return Values
Output event to wait on to ensure computation is complete.