getri (USM Version)#
Computes the inverse of an LU-factored general matrix determined by
getrf (USM Version). This routine belongs to the
oneapi::mkl::lapack namespace.
Description#
The routine computes the inverse inv(A) of a general matrix
A. Before calling this routine, call getrf (USM Version) to
factorize A.
API#
Syntax#
namespace oneapi::mkl::lapack {
sycl::event getri(sycl::queue &queue,
int64_t n,
T *a,
int64_t lda,
const int64_t *ipiv,
T *scratchpad,
int64_t scratchpad_size,
const std::vector<sycl::event> &events = {})
}
getri (USM version) supports the following precisions and
devices:
T |
Devices supported |
|---|---|
|
CPU and GPU |
|
CPU and GPU |
|
CPU and GPU |
|
CPU and GPU |
Input Parameters#
- queue
Device queue where calculations will be performed.
- n
The order of the matrix
A(0 ≤ n).- a
The pointer returned by getrf (USM Version). Must be of size at least
lda*max(1,n).- lda
The leading dimension of
a(n ≤ lda).- ipiv
The array as returned by getrf (USM Version). The dimension of ipiv must be at least
max(1, n).- scratchpad
Pointer to scratchpad memory to be used by the 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 getri_scratchpad_size function.- events
List of events to wait for before starting computation. Defaults to empty list.
Output Parameters#
- a
Overwritten by the n-by-n matrix
inv(A).
Exceptions#
Exception |
Description |
|---|---|
|
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 If |
Return Values#
Output event to wait on to ensure computation is complete.