potrs#

Solves a system of linear equations with a Cholesky-factored symmetric (Hermitian) positive-definite coefficient matrix.

Description

potrs supports the following precisions.

T

float

double

std::complex<float>

std::complex<double>

The routine solves for \(X\) the system of linear equations \(AX = B\) with a symmetric positive-definite or, for complex data, Hermitian positive-definite matrix \(A\), given the Cholesky factorization of \(A\):

\(A = U^TU\) for real data, \(A = U^HU\) for complex data

if upper_lower=oneapi::mkl::uplo::upper

\(A = LL^T\) for real data, \(A = LL^H\) for complex data

if upper_lower=oneapi::mkl::uplo::lower

where \(L\) is a lower triangular matrix and \(U\) is upper triangular. The system is solved with multiple right-hand sides stored in the columns of the matrix \(B\).

Before calling this routine, you must call potrf to compute the Cholesky factorization of \(A\).

potrs (Buffer Version)#

Syntax

namespace oneapi::mkl::lapack {
  void potrs(sycl::queue &queue, oneapi::mkl::uplo upper_lower, std::int64_t n, std::int64_t nrhs, sycl::buffer<T,1> &a, std::int64_t lda, sycl::buffer<T,1> &b, std::int64_t ldb, 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 has been factored:

If upper_lower = oneapi::mkl::uplo::upper, the upper triangle \(U\) of \(A\) is stored, where \(A\) = \(U^{T}`U\) for real data, \(A\) = \(U^{H}U\) for complex data.

If upper_lower = oneapi::mkl::uplo::lower, the lower triangle \(L\) of \(A\) is stored, where \(A\) = \(LL^{T}\) for real data, \(A\) = \(LL^{H}\) for complex data.

n

The order of matrix \(A\) (\(0 \le n\)).

nrhs

The number of right-hand sides (\(0 \le \text{nrhs}\)).

a

Buffer containing 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.

b

The array b contains the matrix \(B\) whose columns are the right-hand sides for the systems of equations. The second dimension of b must be at least \(\max(1,\text{nrhs})\).

ldb

The leading dimension of b.

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 potrs_scratchpad_size function.

Output Parameters

b

Overwritten by the solution matrix \(X\).

scratchpad

Buffer holding scratchpad memory to be used by routine for storing intermediate results.

potrs (USM Version)#

Syntax

namespace oneapi::mkl::lapack {
  sycl::event potrs(sycl::queue &queue, oneapi::mkl::uplo upper_lower, std::int64_t n, std::int64_t nrhs, T *a, std::int64_t lda, T *b, std::int64_t ldb, 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 has been factored:

If upper_lower = oneapi::mkl::uplo::upper, the upper triangle \(U\) of \(A\) is stored, where \(A\) = \(U^{T}U\) for real data, \(A\) = \(U^{H}U\) for complex data.

If upper_lower = oneapi::mkl::uplo::lower, the lower triangle \(L\) of \(A\) is stored, where \(A\) = \(LL^{T}\) for real data, \(A\) = \(LL^{H}\) for complex data.

n

The order of matrix \(A\) (\(0 \le n\)).

nrhs

The number of right-hand sides (\(0 \le \text{nrhs}\)).

a

Pointer to array containing 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.

b

The array b contains the matrix \(B\) whose columns are the right-hand sides for the systems of equations. The second dimension of b must be at least \(\max(1,\text{nrhs})\).

ldb

The leading dimension of b.

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 potrs_scratchpad_size function.

events

List of events to wait for before starting computation. Defaults to empty list.

Output Parameters

b

Overwritten by the solution matrix \(X\).

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