potri#
Computes the inverse of a symmetric (Hermitian) positive-definite matrix using the Cholesky factorization.
Description
potri
supports the following precisions.
T
float
double
std::complex<float>
std::complex<double>
The routine computes the inverse \(A^{-1}\) of a symmetric positive definite or, for complex flavors, Hermitian positive-definite matrix \(A\). Before calling this routine, call potrf to factorize \(A\).
potri (Buffer Version)#
Syntax
namespace oneapi::mkl::lapack {
void potri(sycl::queue &queue, oneapi::mkl::uplo upper_lower, std::int64_t n, sycl::buffer<T,1> &a, std::int64_t lda, sycl::buffer<T,1> &scratchpad, std::int64_t scratchpad_size)
}
Input Parameters
- queue
The queue where the routine should be executed.
- upper_lower
Indicates how the input matrix \(A\) has been factored:
If
upper_lower = oneapi::mkl::uplo::upper
, the upper triangle of \(A\) is stored.If
upper_lower = oneapi::mkl::uplo::lower
, the lower triangle of \(A\) is stored.- n
Specifies the order of the matrix \(A\) (\(0 \le n\)).
- a
Contains the factorization of the matrix \(A\), as returned by potrf. The second dimension of
a
must be at least \(\max(1, n)\).- lda
The leading dimension of
a
.- 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 potri_scratchpad_size function.
Output Parameters
- a
Overwritten by the upper or lower triangle of the inverse of \(A\). Specified by
upper_lower
.- scratchpad
Buffer holding scratchpad memory to be used by routine for storing intermediate results.
potri (USM Version)#
Syntax
namespace oneapi::mkl::lapack {
sycl::event potri(sycl::queue &queue, oneapi::mkl::uplo upper_lower, std::int64_t n, T *a, std::int64_t lda, T *scratchpad, std::int64_t scratchpad_size, const std::vector<sycl::event> &events = {})
}
Input Parameters
- queue
The queue where the routine should be executed.
- upper_lower
Indicates how the input matrix \(A\) has been factored:
If
upper_lower = oneapi::mkl::uplo::upper
, the upper triangle of \(A\) is stored.If
upper_lower = oneapi::mkl::uplo::lower
, the lower triangle of \(A\) is stored.- n
Specifies the order of the matrix \(A\) (\(0 \le n\)).
- a
Contains the factorization of the matrix \(A\), as returned by potrf. The second dimension of
a
must be at least \(\max(1, n)\).- lda
The leading dimension of
a
.- 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 potri_scratchpad_size function.- events
List of events to wait for before starting computation. Defaults to empty list.
Output Parameters
- a
Overwritten by the upper or lower triangle of the inverse of \(A\). Specified by
upper_lower
.- scratchpad
Pointer to scratchpad memory to be used by routine for storing intermediate results.
Return Values
Output event to wait on to ensure computation is complete.
Parent topic: LAPACK Linear Equation Routines