potrf_batch¶
Computes the LU factorizations of a batch of general matrices.
Description
potrf_batch
supports the following precisions.
T
float
double
std::complex<float>
std::complex<double>
potrf_batch (Buffer Version)¶
Description
The buffer version of potrf_batch
supports only the strided API.
Strided API
The routine forms the Cholesky factorizations of a symmetric positive-definite or, for complex data, Hermitian positive-definite matrices \(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.
Syntax
namespace oneapi::mkl::lapack {
void potrf_batch(cl::sycl::queue &queue, mkl::uplo uplo, std::int64_t n, cl::sycl::buffer<T> &a, std::int64_t lda, std::int64_t stride_a, 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 whether the upper or lower triangular part of \(A_i\) is stored and how \(A_i\) is factored:If
uplo = mkl::uplo::upper
, the arraya
stores the upper triangular parts of the matrices \(A_i\),Ifuplo = mkl::uplo::lower
, the arraya
stores the lower triangular parts of the matrices \(A_i\). - n
Order of the matrices \(A_i\), (\(0 \le n\)).
- a
Array containing batch of input matrices \(A_i\), each of \(A_i\) being of size \(\text{lda} \cdot n\) and holding either upper or lower triangular parts of the matrices \(A_i\) (see
uplo
).- lda
Leading dimension of \(A_i\).
- stride_a
Stride between the beginnings of matrices \(A_i\) inside the batch.
- 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 potrf_batch_scratchpad_size function.
Output Parameters
- a
Cholesky factors \(U_i\) or \(L_i\), as specified by
uplo
.
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 the leading minors of some of matrices (and therefore some matrices \(A_i\) themselves) are not positive-definite, and the factorizations could not be completed for these matrices from the batch. The indices of such matrices in the batch can be obtained with ids() method of the exception object. The orders of corresponding not positive-definite leading minors of these matrices can be obtained by exceptions() method of exception object.
potrf_batch (USM Version)¶
Description
The USM version of potrf_batch
supports the group API and strided API.
Group API
The routine forms the Cholesky factorizations of symmetric positive-definite or, for complex data, Hermitian positive-definite matrices \(A_i\), \(i \in \{1...batch\_size\}\):\(A_i = U_i^TU_i\) for real data (\(A_i = U_i^HU_i\) for complex), if \(\text{uplo}_g\) ismkl::uplo::upper
,\(A_i = L_iL_i^T\) for real data (\(A_i = L_iL_i^H\) for complex), if \(\text{uplo}_g\) ismkl::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 total number of problems to solve,batch_size
, is a sum of sizes of all of the groups of parameters as provided bygroup_sizes
array
Syntax
namespace oneapi::mkl::lapack {
cl::sycl::event potrf_batch(cl::sycl::queue &queue, mkl::uplo *uplo, std::int64_t *n, T **a, std::int64_t *lda, 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 \(\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 array a belonging to group \(g\).- a
Array of
batch_size
pointers to input matrices \(A_i\), each being of size \(\text{lda}_g \cdot n_g\) (\(g\) is an index of group to which \(A_i\) belongs to) and holding either upper or lower triangular part as specified by \(\text{uplo}_g\).- lda
Array of
group_count
\(\text{lda}_g\) parameters. Each \(\text{lda}_g\) specifies the leading dimensions of the matrices from a 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 potrf_batch_scratchpad_size function.- events
List of events to wait for before starting computation. Defaults to empty list.
Output Parameters
- a
Cholesky factors \(U_i\) or \(L_i\), as specified by \(\text{uplo}_g\) from corresponding group of parameters.
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 the leading minors of some of the input matrices (and therefore some matrices themselves) are not positive-definite, and the factorizations could not be completed for these matrices from the batch. The indices of such matrices in the batch can be obtained with ids() method of the exception object. The orders of corresponding not positive-definite leading minors of these matrices can be obtained by exceptions() method of the exception object.
Strided API
The routine forms the Cholesky factorizations of a symmetric positive-definite or, for complex data, Hermitian positive-definite matrices \(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.
Syntax
namespace oneapi::mkl::lapack {
cl::sycl::event potrf_batch(cl::sycl::queue &queue, mkl::uplo uplo, std::int64_t n, T *a, std::int64_t lda, std::int64_t stride_a, 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 whether the upper or lower triangular part of \(A_i\) is stored and how \(A_i\) is factored:If
uplo = mkl::uplo::upper
, the arraya
stores the upper triangular parts of the matrices \(A_i\),Ifuplo = mkl::uplo::lower
, the arraya
stores the lower triangular parts of the matrices \(A_i\). - n
Order of the matrices \(A_i\), (\(0 \le n\)).
- a
Array containing batch of input matrices \(A_i\), each of \(A_i\) being of size \(\text{lda} \cdot n\) and holding either upper or lower triangular parts of the matrices \(A_i\) (see
uplo
).- lda
Leading dimension of \(A_i\).
- stride_a
Stride between the beginnings of matrices \(A_i\) inside the batch.
- 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 potrf_batch_scratchpad_size function.- events
List of events to wait for before starting computation. Defaults to empty list.
Output Parameters
- a
Cholesky factors \(U_i\) or \(L_i\), as specified by
uplo
.
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 the leading minors of some of matrices (and therefore some matrices \(A_i\) themselves) are not positive-definite, and the factorizations could not be completed for these matrices from the batch. The indices of such matrices in the batch can be obtained with ids() method of the exception object. The orders of corresponding not positive-definite leading minors of these matrices can be obtained by exceptions() method of exception object.
Parent topic: LAPACK-like Extensions Routines