trsm#

Solves a triangular matrix equation (forward or backward solve).

Description

The trsm routines solves a triangular matrix equations. There are two operations available, an in-place operation and an out-of-place operation. The in-place operation solves for X in:

\[op(A)*X = alpha*B\]

or

\[X*op(A) = alpha*B\]

The out-of-place operation solves for X and then adds that solution to a scaled matrix C:

\[op(A)*X = alpha*B, C \leftarrow X + beta*C\]

or

\[X*op(A) = alpha*B, C \leftarrow X + beta*C\]

where:

op(A) is one of op(A) = A, or op(A) = AT, or op(A) = AH,

alpha and beta are scalars,

A is a triangular matrix, and

B, X, and C are m x n general matrices.

A is either m x m or n x n, depending on whether it multiplies X on the left or right.

For the in-place operation, the matrix B is overwritten by the solution matrix X on return. For the out-of-place operation, B remains untouched and the solution is added to a scaled C matrix.

trsm supports the following precisions.

T

float

double

std::complex<float>

std::complex<double>

trsm (Buffer Version)#

In-place API#

Syntax

namespace oneapi::mkl::blas::column_major {
    void trsm(sycl::queue &queue,
              oneapi::mkl::side left_right,
              oneapi::mkl::uplo upper_lower,
              oneapi::mkl::transpose transa,
              oneapi::mkl::diag unit_diag,
              std::int64_t m,
              std::int64_t n,
              T alpha,
              sycl::buffer<T,1> &a,
              std::int64_t lda,
              sycl::buffer<T,1> &b,
              std::int64_t ldb)
}
namespace oneapi::mkl::blas::row_major {
    void trsm(sycl::queue &queue,
              oneapi::mkl::side left_right,
              oneapi::mkl::uplo upper_lower,
              oneapi::mkl::transpose transa,
              oneapi::mkl::diag unit_diag,
              std::int64_t m,
              std::int64_t n,
              T alpha,
              sycl::buffer<T,1> &a,
              std::int64_t lda,
              sycl::buffer<T,1> &b,
              std::int64_t ldb)
}

Input Parameters

queue

The queue where the routine should be executed.

left_right

Specifies whether A multiplies X on the left (side::left) or on the right (side::right). See oneMKL defined datatypes for more details.

upper_lower

Specifies whether the matrix A is upper or lower triangular. See oneMKL defined datatypes for more details.

trans

Specifies op(A), the transposition operation applied to A. See oneMKL defined datatypes for more details.

unit_diag

Specifies whether A is assumed to be unit triangular (all diagonal elements are 1). See oneMKL defined datatypes for more details.

m

Specifies the number of rows of B. The value of m must be at least zero.

n

Specifies the number of columns of B. The value of n must be at least zero.

alpha

Scaling factor for the solution.

a

Buffer holding input matrix A. Must have size at least lda*m if left_right = side::left, or lda*n if left_right = side::right. See Matrix Storage for more details.

lda

Leading dimension of A. Must be at least m if left_right = side::left, and at least n if left_right = side::right. Must be positive.

b

Buffer holding input/output matrix B. Must have size at least ldb*n if column major layout is used to store matrices or at least ldb*m if row major layout is used to store matrices. See Matrix Storage for more details.

ldb

Leading dimension of B. It must be positive and at least m if column major layout is used to store matrices or at least n if row major layout is used to store matrices.

Output Parameters

b

Output buffer. Overwritten by the solution matrix X.

Notes

If alpha = 0, matrix B is set to zero, and A and B do not need to be initialized at entry.

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

oneapi::mkl::unsupported_device

oneapi::mkl::host_bad_alloc

oneapi::mkl::device_bad_alloc

oneapi::mkl::unimplemented

Out-of-place API#

Syntax

namespace oneapi::mkl::blas::column_major {
    void trsm(sycl::queue &queue,
              oneapi::mkl::side left_right,
              oneapi::mkl::uplo upper_lower,
              oneapi::mkl::transpose trans,
              oneapi::mkl::diag unit_diag,
              std::int64_t m,
              std::int64_t n,
              T alpha,
              sycl::buffer<T,1> &a,
              std::int64_t lda,
              sycl::buffer<T,1> &b,
              std::int64_t ldb,
              T beta,
              sycl::buffer<T,1> &c,
              std::int64_t ldc)
}
namespace oneapi::mkl::blas::row_major {
    void trsm(sycl::queue &queue,
              oneapi::mkl::side left_right,
              oneapi::mkl::uplo upper_lower,
              oneapi::mkl::transpose trans,
              oneapi::mkl::diag unit_diag,
              std::int64_t m,
              std::int64_t n,
              T alpha,
              sycl::buffer<T,1> &a,
              std::int64_t lda,
              sycl::buffer<T,1> &b,
              std::int64_t ldb,
              T beta,
              sycl::buffer<T,1> &c,
              std::int64_t ldc)
}

Input Parameters

queue

The queue where the routine should be executed.

left_right

Specifies whether A is on the left side of the matrix solve (side::left) or on the right side (side::right). See oneMKL defined datatypes for more details.

upper_lower

Specifies whether the matrix A is upper or lower triangular. See oneMKL defined datatypes for more details.

trans

Specifies op(A), the transposition operation applied to matrix A. See oneMKL defined datatypes for more details.

unit_diag

Specifies whether A is assumed to be unit triangular (all diagonal elements are 1). See oneMKL defined datatypes for more details.

m

Specifices the number of rows of B. The value of m must be at least zero.

n

Specifies the Number of columns of B. The value of n must be at least zero.

alpha

Scaling factor for the solution.

a

Buffer holding input matrix A. Must have size at least lda*m if left_right = side::left or lda*n if left_right = side::right. See Matrix Storage for more details.

lda

Leading dimension of A. Must be at least m if left_right = side::left or at least n if left_right = side::right. Must be positive.

b

Buffer holding input matrix B. Must have size at least ldb*n if column major layout or at least ldb*m if row major layout is used. See Matrix Storage for more details.

ldb

Leading dimension of matrix B. It must be positive and at least m if column major layout or at least n if row major layout is used.

beta

Scaling factor for matrix C.

c

Buffer holding input/output matrix C. Size of the buffer must be at least ldc*n if column major layout or at least ldc*m if row major layout is used. See Matrix Storage for more details.

ldc

Leading dimension of matrix C. Must be at least m if column major layout or at least n if row major layout is used. Must be positive.

Output Parameters

c

Output buffer overwritten by solution matrix X + beta*C.

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

oneapi::mkl::unsupported_device

oneapi::mkl::host_bad_alloc

oneapi::mkl::device_bad_alloc

oneapi::mkl::unimplemented

trsm (USM Version)#

In-place API#

Syntax

namespace oneapi::mkl::blas::column_major {
    sycl::event trsm(sycl::queue &queue,
                     oneapi::mkl::side left_right,
                     oneapi::mkl::uplo upper_lower,
                     oneapi::mkl::transpose transa,
                     oneapi::mkl::diag unit_diag,
                     std::int64_t m,
                     std::int64_t n,
                     value_or_pointer<T> alpha,
                     const T *a,
                     std::int64_t lda,
                     T *b,
                     std::int64_t ldb,
                     const std::vector<sycl::event> &dependencies = {})
}
namespace oneapi::mkl::blas::row_major {
    sycl::event trsm(sycl::queue &queue,
                     oneapi::mkl::side left_right,
                     oneapi::mkl::uplo upper_lower,
                     oneapi::mkl::transpose transa,
                     oneapi::mkl::diag unit_diag,
                     std::int64_t m,
                     std::int64_t n,
                     value_or_pointer<T> alpha,
                     const T *a,
                     std::int64_t lda,
                     T *b,
                     std::int64_t ldb,
                     const std::vector<sycl::event> &dependencies = {})
}

Input Parameters

queue

The queue where the routine should be executed.

left_right

Specifies whether A multiplies X on the left (side::left) or on the right (side::right). See oneMKL defined datatypes for more details.

upper_lower

Specifies whether the matrix A is upper or lower triangular. See oneMKL defined datatypes for more details.

transa

Specifies op(A), the transposition operation applied to A. See oneMKL defined datatypes for more details.

unit_diag

Specifies whether A is assumed to be unit triangular (all diagonal elements are 1). See oneMKL defined datatypes for more details.

m

Specifies the number of rows of B. The value of m must be at least zero.

n

Specifies the number of columns of B. The value of n must be at least zero.

alpha

Scaling factor for the solution. See Scalar Arguments in BLAS for more details.

a

Pointer to input matrix A. Must have size at least lda*m if left_right = side::left, or lda*n if left_right = side::right. See Matrix Storage for more details.

lda

Leading dimension of A. Must be at least m if left_right = side::left, and at least n if left_right = side::right. Must be positive.

b

Pointer to input/output matrix B. Must have size at least ldb*n if column major layout is used to store matrices or at least ldb*m if row major layout is used to store matrices. See Matrix Storage for more details.

ldb

Leading dimension of B. It must be positive and at least m if column major layout is used to store matrices or at least n if row major layout is used to store matrices.

dependencies

List of events to wait for before starting computation, if any. If omitted, defaults to no dependencies.

Output Parameters

b

Pointer to the output matrix. Overwritten by the solution matrix X.

Notes

If alpha = 0, matrix B is set to zero, and A and B do not need to be initialized at entry.

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::mkl::invalid_argument

oneapi::mkl::unsupported_device

oneapi::mkl::host_bad_alloc

oneapi::mkl::device_bad_alloc

oneapi::mkl::unimplemented

Out-of-place API#

Syntax

namespace oneapi::mkl::blas::column_major {
    void trsm(sycl::queue &queue,
              oneapi::mkl::side left_right,
              oneapi::mkl::uplo upper_lower,
              oneapi::mkl::transpose trans,
              oneapi::mkl::diag unit_diag,
              std::int64_t m,
              std::int64_t n,
              value_or_pointer<T> alpha,
              const T *a,
              std::int64_t lda,
              const T *b,
              std::int64_t ldb,
              value_or_pointer<T> beta,
              T *c,
              std::int64_t ldc,
              const std::vector<sycl::event> &dependencies = {})
}
namespace oneapi::mkl::blas::row_major {
    void trsm(sycl::queue &queue,
              oneapi::mkl::side left_right,
              oneapi::mkl::uplo upper_lower,
              oneapi::mkl::transpose trans,
              oneapi::mkl::diag unit_diag,
              std::int64_t m,
              std::int64_t n,
              value_or_pointer<T> alpha,
              const T *a,
              std::int64_t lda,
              const T *b,
              std::int64_t ldb,
              value_or_pointer<T> beta,
              T *c,
              std::int64_t ldc,
              const std::vector<sycl::event> &dependencies = {})
}

Input Parameters

queue

The queue where the routine should be executed.

left_right

Specifies whether A is on the left side of the matrix solve (side::left) or on the right side (side::right). See oneMKL defined datatypes for more details.

upper_lower

Specifies whether the matrix A is upper or lower triangular. See oneMKL defined datatypes for more details.

trans

Specifies op(A), the transposition operation applied to matrix A. See oneMKL defined datatypes for more details.

unit_diag

Specifies whether A is assumed to be unit triangular (all diagonal elements are 1). See oneMKL defined datatypes for more details.

m

Specifices the number of rows of B. The value of m must be at least zero.

n

Specifies the Number of columns of B. The value of n must be at least zero.

alpha

Scaling factor for the solution. See Scalar Arguments in BLAS for more details.

a

Pointer to input matrix A. Must have size at least lda*m if left_right = side::left or lda*n if left_right = side::right. See Matrix Storage for more details.

lda

Leading dimension of A. Must be at least m if left_right = side::left or at least n if left_right = side::right. Must be positive.

b

Pointer to input matrix B. Must have size at least ldb*n if column major layout or at least ldb*m if row major layout is used. See Matrix Storage for more details.

ldb

Leading dimension of matrix B. It must be positive and at least m if column major layout or at least n if row major layout is used.

beta

Scaling factor for matrix C. See Scalar Arguments in BLAS for more details.

c

Pointer to input/output matrix C. Must have size at least ldc*n if column major layout or at least ldc*m if row major layout is used. See Matrix Storage for more details.

ldc

Leading dimension of matrix C. Must be at least m if column major layout or at least n if row major layout is used. Must be positive.

dependencies

List of events to wait for before starting computation, if any. If omitted, defaults to no dependencies.

Output Parameters

c

Pointer to the output matrix, overwritten by the solution matrix X + beta*C.

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

oneapi::mkl::unsupported_device

oneapi::mkl::host_bad_alloc

oneapi::mkl::device_bad_alloc

oneapi::mkl::unimplemented

Parent topic: BLAS Level 3 Routines