trtrs#

Solves a system of linear equations with a triangular coefficient matrix, with multiple right-hand sides.

Description

trtrs supports the following precisions.

T

float

double

std::complex<float>

std::complex<double>

The routine solves for \(X\) the following systems of linear equations with a triangular matrix \(A\), with multiple right-hand sides stored in \(B\):

\(AX = B\)

if transa =transpose::nontrans,

\(A^TX = B\)

if transa =transpose::trans,

\(A^HX = B\)

if transa =transpose::conjtrans (for complex matrices only).

trtrs (Buffer Version)#

Syntax

namespace oneapi::mkl::lapack {
  void trtrs(cl::sycl::queue &queue, oneapi::mkl::uplo upper_lower, oneapi::mkl::transpose transa, oneapi::mkl::diag unit_diag, std::int64_t n, std::int64_t nrhs, cl::sycl::buffer<T,1> &a, std::int64_t lda, cl::sycl::buffer<T,1> &b, std::int64_t ldb, cl::sycl::buffer<T,1> &scratchpad, std::int64_t scratchpad_size)
}

Input Parameters

queue

The queue where the routine should be executed.

upper_lower

Indicates whether \(A\) is upper or lower triangular:

If upper_lower = uplo::upper, then \(A\) is upper triangular.

If upper_lower = uplo::lower, then \(A\) is lower triangular.

transa

If transa = transpose::nontrans, then \(AX = B\) is solved for \(X\).

If transa = transpose::trans, then \(A^{T}X = B\) is solved for \(X\).

If transa = transpose::conjtrans, then \(A^{H}X = B\) is solved for \(X\).

unit_diag

If unit_diag = diag::nonunit, then \(A\) is not a unit triangular matrix.

If unit_diag = diag::unit, then \(A\) is unit triangular: diagonal elements of \(A\) are assumed to be 1 and not referenced in the array a.

n

The order of \(A\); the number of rows in \(B\); \(n \ge 0\).

nrhs

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

a

Buffer containing the matrix \(A\). The second dimension of a must be at least \(\max(1,n)\).

lda

The leading dimension of a; \(\text{lda} \ge \max(1, n)\).

b

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

ldb

The leading dimension of b; \(\text{ldb} \ge \max(1, n)\).

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 trtrs_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.

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::host_bad_alloc

oneapi::mkl::device_bad_alloc

oneapi::mkl::unimplemented

oneapi::mkl::unsupported_device

oneapi::mkl::lapack::invalid_argument

oneapi::mkl::lapack::computation_error

Exception is thrown in case of problems during calculations. The info code of the problem can be obtained by info() method of exception object:

If \(\text{info}=-i\), the \(i\)-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 not be less than value return by detail() method of exception object.

trtrs (USM Version)#

Syntax

namespace oneapi::mkl::lapack {
  cl::sycl::event trtrs(cl::sycl::queue &queue, oneapi::mkl::uplo upper_lower, oneapi::mkl::transpose transa, oneapi::mkl::diag unit_diag, std::int64_t n, std::int64_t nrhs, const T *a, std::int64_t lda, T *b, std::int64_t ldb, T *scratchpad, std::int64_t scratchpad_size, const std::vector<cl::sycl::event> &events = {})
}

Input Parameters

queue

The queue where the routine should be executed.

upper_lower

Indicates whether \(A\) is upper or lower triangular:

If upper_lower = uplo::upper, then \(A\) is upper triangular.

If upper_lower = uplo::lower, then \(A\) is lower triangular.

transa

If transa = transpose::nontrans, then \(AX = B\) is solved for \(X\).

If transa = transpose::trans, then \(A^{T}X = B\) is solved for \(X\).

If transa = transpose::conjtrans, then \(A^{H}X = B\) is solved for \(X\).

unit_diag

If unit_diag = diag::nonunit, then \(A\) is not a unit triangular matrix.

If unit_diag = diag::unit, then \(A\) is unit triangular: diagonal elements of \(A\) are assumed to be 1 and not referenced in the array a.

n

The order of \(A\); the number of rows in \(B\); \(n \ge 0\).

nrhs

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

a

Array containing the matrix \(A\). The second dimension of a must be at least \(\max(1,n)\).

lda

The leading dimension of a; \(\text{lda} \ge \max(1, n)\).

b

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

ldb

The leading dimension of b; \(\text{ldb} \ge \max(1, n)\).

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 trtrs_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.

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::host_bad_alloc

oneapi::mkl::device_bad_alloc

oneapi::mkl::unimplemented

oneapi::mkl::unsupported_device

oneapi::mkl::lapack::invalid_argument

oneapi::mkl::lapack::computation_error

Exception is thrown in case of problems during calculations. The info code of the problem can be obtained by info() method of exception object:

If \(\text{info}=-i\), the \(i\)-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 not be less than value return by detail() method of exception object.

Return Values

Output event to wait on to ensure computation is complete.

Parent topic: LAPACK Linear Equation Routines