potrs_batch#
Computes the LU factorizations of a batch of general matrices.
Description
potrs_batch
supports the following precisions.
T
float
double
std::complex<float>
std::complex<double>
potrs_batch (Buffer Version)#
Description
The buffer version of potrs_batch
supports only the strided API.
Strided API
The routine solves for \(X_i\) the systems of linear equations \(A_iX_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 \in \{1...batch\_size\}\):\(A_i = U_i^TU_i\) for real data, \(A_i = U_i^HU_i\) for complex data ifuplo = mkl::uplo::upper
,\(A_i = L_iL_i^T\) for real data, \(A_i = L_iL_i^H\) for complex data ifuplo = mkl::uplo::lower
,where \(L_i\) is a lower triangular matrix and \(U_i\) is upper triangular.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 call to the Strided API of the potrf_batch (Buffer Version) function.
Syntax
namespace oneapi::mkl::lapack {
void potrs_batch(cl::sycl::queue &queue, mkl::uplo uplo, std::int64_t n, std::int64_t nrhs, cl::sycl::buffer<T> &a, std::int64_t lda, std::int64_t stride_a, cl::sycl::buffer<T> &b, std::int64_t ldb, std::int64_t stride_b, std::int64_t batch_size, cl::sycl::buffer<T> &scratchpad, std::int64_t scratchpad_size)
}
Input Parameters
- queue
Device queue where calculations will be performed.
- uplo
- Indicates how the input matrices have been factored:If
uplo = mkl::uplo::upper
, the upper triangle \(U_i\) of \(A_i\) is stored, where \(A_i = U_i^TU_i\) for real data, \(A_i = U_i^HU_i\) for complex data.Ifuplo = mkl::uplo::lower
, the upper triangle \(L_i\) of \(A_i\) is stored, where \(A_i = L_iL_i^T\) for real data, \(A_i = L_iL_i^H\) for complex data. - n
The order of matrices \(A_i\) (\(0 \le n\)).
- nrhs
The number of right-hand sides (\(0 \le \text{nrhs}\)).
- a
Array containing batch of factorizations of the matrices \(A_i\), as returned by the Strided API of the potrf_batch (Buffer Version) function.
- lda
Leading dimension of \(A_i\).
- stride_a
Stride between the beginnings of matrices inside the batch array
a
.- b
Array containing batch of matrices \(B_i\) whose columns are the right-hand sides for the systems of equations.
- ldb
Leading dimension of \(B_i\).
- stride_b
Stride between the beginnings of matrices \(B_i\) inside the batch array
b
.- batch_size
Number of problems in a batch.
- 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 the Strided API of the potrs_batch_scratchpad_size function.
Output Parameters
- b
Solution matrices \(X_i\).
Throws
This routine shall throw the following exceptions if the associated condition is detected. An implementation may throw additional implementation-specific exception(s) in case of error conditions not covered here.
oneapi::mkl::lapack::batch_error
oneapi::mkl::unsupported_device
oneapi::mkl::lapack::invalid_argument
The
info
code of the problem can be obtained by info() method of exception object:If
info = -n
, the \(n\)-th parameter had an illegal value.If
info
equals to value passed as scratchpad size, and detail() returns non zero, then passed scratchpad is of insufficient size, and required size should be not less then value returned by detail() method of exception object.If
info
is not zero and detail() returns zero, then there were some errors for some of the problems in the supplied batch andinfo
code contains the number of failed calculations in a batch.If
info
is zero, then for some of the matrices diagonal element of the Cholesky factor is zero, and the solve could not be completed. The indices of such matrices in the batch can be obtained with ids() method of the exception object. The indices of first zero diagonal elements in these matrices can be obtained by exceptions() method of exception object.
potrs_batch (USM Version)#
Description
The USM version of potrs_batch
supports the group API and strided API.
Group 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, const T * const *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 = {})
}
Input Parameters
- queue
Device queue where calculations will be performed.
- uplo
- Array of
group_count
\(\text{uplo}_g\) parameters.Each of \(\text{uplo}_g\) indicates whether the upper or lower triangular parts of the input matrices are provided:If \(\text{uplo}_g\) ismkl::uplo::upper
, input matrices from arraya
belonging to group \(g\) store the upper triangular parts,If \(\text{uplo}_g\) ismkl::uplo::lower
, input matrices from arraya
belonging to group \(g\) store the lower triangular parts. - n
- Array of
group_count
\(n_g\) parameters.Each \(n_g\) specifies the order of the input matrices from arraya
belonging to group \(g\). - nrhs
- Array of
group_count
\(\text{nrhs}_g\) parameters.Each \(\text{nrhs}_g\) specifies the number of right-hand sides supplied for group \(g\) in corresponding part of arrayb
. - a
Array of
batch_size
pointers to Cholesky factored matrices \(A_i\) as returned by the Group API of the potrf_batch (USM Version) function.- lda
- Array of
group_count
\(\text{lda}_g\) parameters.Each \(\text{lda}_g\) specifies the leading dimensions of the matrices froma
belonging to group \(g\). - b
Array of
batch_size
pointers to right-hand side matrices \(B_i\), each of size \(\text{ldb}_g \cdot \text{nrhs}_g\), where \(g\) is an index of group corresponding to \(B_i\).- ldb
- Array of
group_count
\(\text{ldb}_g\) parameters.Each \(\text{ldb}_g\) specifies the leading dimensions of the matrices fromb
belonging to group \(g\). - group_count
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 the Group API of the potrs_batch_scratchpad_size function.- events
List of events to wait for before starting computation. Defaults to empty list.
Output Parameters
- b
Solution matrices \(X_i\).
Return Values
Output event to wait on to ensure computation is complete.
Throws
This routine shall throw the following exceptions if the associated condition is detected. An implementation may throw additional implementation-specific exception(s) in case of error conditions not covered here.
oneapi::mkl::lapack::batch_error
oneapi::mkl::unsupported_device
oneapi::mkl::lapack::invalid_argument
The
info
code of the problem can be obtained by info() method of exception object:If
info = -n
, the n-th parameter had an illegal value.If
info
equals to value passed as scratchpad size, and detail() returns non zero, then passed scratchpad is of insufficient size, and required size should be not less then value returned by detail() method of exception object.If
info
is not zero and detail() returns zero, then there were some errors for some of the problems in the supplied batch andinfo
code contains the number of failed calculations in a batch.If
info
is zero, then for some of the matrices diagonal element of the Cholesky factor is zero, and the solve could not be completed. The indices of such matrices in the batch can be obtained with ids() method of the exception object. The indices of first zero diagonal elements in these matrices can be obtained by exceptions() method of exception object.
Strided API
The routine solves for \(X_i\) the systems of linear equations \(A_iX_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 \in \{1...batch\_size\}\):\(A_i = U_i^TU_i\) for real data, \(A_i = U_i^HU_i\) for complex data ifuplo = mkl::uplo::upper
,\(A_i = L_iL_i^T\) for real data, \(A_i = L_iL_i^H\) for complex data ifuplo = mkl::uplo::lower
,where \(L_i\) is a lower triangular matrix and \(U_i\) is upper triangular.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 call to the Strided API of the potrf_batch (USM Version) function.
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, const T *a, std::int64_t lda, std::int64_t stride_a, T *b, std::int64_t ldb, std::int64_t stride_b, std::int64_t batch_size, T *scratchpad, std::int64_t scratchpad_size, const std::vector<cl::sycl::event> &events = {})
};
Input Parameters
- queue
Device queue where calculations will be performed.
- uplo
- Indicates how the input matrices have been factored:If
uplo = mkl::uplo::upper
, the upper triangle \(U_i\) of \(A_i\) is stored, where \(A_i = U_i^TU_i\) for real data, \(A_i = U_i^HU_i\) for complex data.Ifuplo = mkl::uplo::lower
, the upper triangle \(L_i\) of \(A_i\) is stored, where \(A_i = L_iL_i^T\) for real data, \(A_i = L_iL_i^H\) for complex data. - n
The order of matrices \(A_i\) (\(0 \le n\)).
- nrhs
The number of right-hand sides (\(0 \le nrhs\)).
- a
Array containing batch of factorizations of the matrices \(A_i\), as returned by the Strided API of the potrf_batch (USM Version) function.
- lda
Leading dimension of \(A_i\).
- stride_a
Stride between the beginnings of matrices inside the batch array
a
.- b
Array containing batch of matrices \(B_i\) whose columns are the right-hand sides for the systems of equations.
- ldb
Leading dimension of \(B_i\).
- stride_b
Stride between the beginnings of matrices \(B_i\) inside the batch array
b
.- batch_size
Number of problems in a batch.
- 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 the Strided API of the potrs_batch_scratchpad_size function.- events
List of events to wait for before starting computation. Defaults to empty list.
Output Parameters
- b
Solution matrices \(X_i\).
Return Values
Output event to wait on to ensure computation is complete.
Throws
This routine shall throw the following exceptions if the associated condition is detected. An implementation may throw additional implementation-specific exception(s) in case of error conditions not covered here.
oneapi::mkl::lapack::batch_error
oneapi::mkl::unsupported_device
oneapi::mkl::lapack::invalid_argument
Exception is thrown in case of problems during calculations. The info code of the problem can be obtained by info() method of exception object:
If
info = -n
, the \(n\)-th parameter had an illegal value.If
info
equals to value passed as scratchpad size, and detail() returns non zero, then passed scratchpad is of insufficient size, and required size should be not less then value returned by detail() method of exception object.If
info
is not zero and detail() returns zero, then there were some errors for some of the problems in the supplied batch andinfo
code contains the number of failed calculations in a batch.If
info
is zero, then for some of the matrices diagonal element of the Cholesky factor is zero, and the solve could not be completed. The indices of such matrices in the batch can be obtained with ids() method of the exception object. The indices of first zero diagonal elements in these matrices can be obtained by exceptions() method of exception object.
Parent topic: LAPACK-like Extensions Routines