getri_batch#
Computes the inverses of a batch of LU-factored matrices determined by getrf_batch.
Description
getri_batch supports the following precisions.
T
float
double
std::complex<float>
std::complex<double>
getri_batch (Buffer Version)#
Description
The buffer version of getri_batch supports only the strided API.
Strided API
The routine computes the inverses \(A_i^{-1}\) of general matrices \(A_i\). Before calling this routine, call the Strided API of the getrf_batch (Buffer Version) function to factorize \(A_i\).
Syntax
namespace oneapi::math::lapack {
  void getri_batch(cl::sycl::queue &queue, std::int64_t n, cl::sycl::buffer<T> &a, std::int64_t lda, std::int64_t stride_a, cl::sycl::buffer<std::int64_t> &ipiv, std::int64_t stride_ipiv, 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.
- n
 Order of the matrices \(A_i\) (\(0 \le n\)).
- a
 Result of the Strided API of the getrf_batch (Buffer Version) function.
- lda
 Leading dimension of \(A_i\) (\(n\le \text{lda}\)).
- stride_a
 Stride between the beginnings of matrices \(A_i\) inside the batch array
a.- ipiv
 Arrays returned by the Strided API of the getrf_batch (Buffer Version) function.
- stride_ipiv
 Stride between the beginnings of arrays \(\text{ipiv}_i\) inside the array
ipiv.- 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 than the value returned by the Strided API of the getri_batch_scratchpad_size function.
Output Parameters
- a
 Inverse \(n \times n\) matrices \(A_i^{-1}\).
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::math::lapack::batch_error
oneapi::math::unsupported_device
oneapi::math::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
infoequals 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
infois not zero and detail() returns zero, then there were some errors for some of the problems in the supplied batch andinfocode contains the number of failed calculations in a batch.
getri_batch (USM Version)#
Description
The USM version of getri_batch supports the group API and strided API.
Group API
The routine computes the inverses \(A_i^{-1}\) of general matrices \(A_i\), \(i \in \{1...batch\_size\}\). Before calling this routine, call the Group API of the getrf_batch (USM Version) function to factorize \(A_i\).
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.
Syntax
namespace oneapi::math::lapack {
  cl::sycl::event getri_batch(cl::sycl::queue &queue, std::int64_t *n, T **a, std::int64_t *lda, const std::int64_t * const *ipiv, 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.
- n
 Array of
group_count\(n_g\) parameters specifying the order of the matrices \(A_i\) (\(0 \le n_g\)) belonging to group \(g\).- a
 Result of the Group API of the getrf_batch (USM Version) function.
- lda
 Array of
group_count\(\text{lda}_g\) parameters specifying the leading dimensions of the matrices \(A_i\) (\(n_g \le \text{lda}_g\)) belonging to group \(g\).- ipiv
 Arrays returned by the Group API of the getrf_batch (USM Version) function.
- group_count
 Number of groups of parameters. Must be at least 0.
- group_sizes
 Array of
group_countintegers. 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 than the value returned by the Group API of the getri_batch_scratchpad_size function.- events
 List of events to wait for before starting computation. Defaults to empty list.
Output Parameters
- a
 Inverse \(n_g \times n_g\) matrices \(A_i^{-1}\).
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::math::lapack::batch_error
oneapi::math::unsupported_device
oneapi::math::lapack::invalid_argument
The
infocode of the problem can be obtained by info() method of exception object:If
info = -n, the \(n\)-th parameter had an illegal value.If
infoequals 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
infois not zero and detail() returns zero, then there were some errors for some of the problems in the supplied batch andinfocode contains the number of failed calculations in a batch.
Strided API
The routine computes the inverses \(A_i^{-1}\) of general matrices \(A_i\). Before calling this routine, call the Strided API of the getrf_batch (USM Version) function to factorize \(A_i\).
Syntax
namespace oneapi::math::lapack {
  cl::sycl::event getri_batch(cl::sycl::queue &queue, std::int64_t n, T *a, std::int64_t lda, std::int64_t stride_a, const std::int64_t *ipiv, std::int64_t stride_ipiv, 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.
- n
 Order of the matrices \(A_i\) (\(0 \le n\)).
- a
 Result of the Strided API of the getrf_batch (USM Version) function.
- lda
 Leading dimension of \(A_i\) (\(n \le \text{lda}\)).
- stride_a
 Stride between the beginnings of matrices \(A_i\) inside the batch array
a.- ipiv
 Arrays returned by the Strided API of the getrf_batch (USM Version) function.
- stride_ipiv
 Stride between the beginnings of arrays \(\text{ipiv}_i\) inside the array
ipiv.- 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 than the value returned by the Strided API of the getri_batch_scratchpad_size function.- events
 List of events to wait for before starting computation. Defaults to empty list.
Output Parameters
- a
 Inverse \(n \times n\) matrices \(A_i^{-1}\).
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::math::lapack::batch_error
oneapi::math::unsupported_device
oneapi::math::lapack::invalid_argument
The
infocode of the problem can be obtained by info() method of exception object:If
info = -n, the \(n\)-th parameter had an illegal value.If
infoequals 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
infois not zero and detail() returns zero, then there were some errors for some of the problems in the supplied batch andinfocode contains the number of failed calculations in a batch.
Parent topic: LAPACK-like Extensions Routines