trsm¶
Solves a triangular matrix equation (forward or backward solve).
Description
The trsm
routines solve one of the following matrix equations:
or
where:
op(A
) is one of op(A
) = A
, or op(A
) =
A
T, or op(A
) = A
H,
alpha
is a scalar,
A
is a triangular matrix, and
B
and X
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. On return, the matrix B
is overwritten by the solution matrix X
.
trsm
supports the following precisions.
T
float
double
std::complex<float>
std::complex<double>
trsm (Buffer Version)¶
Syntax
namespace oneapi::mkl::blas::column_major {
void trsm(sycl::queue &queue,
onemkl::side left_right,
onemkl::uplo upper_lower,
onemkl::transpose transa,
onemkl::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,
onemkl::side left_right,
onemkl::uplo upper_lower,
onemkl::transpose transa,
onemkl::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
multipliesX
on the left (side::left
) or on the right (side::right
). See oneMKL defined datatypes for more details.- uplo
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 toA
. 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 ofm
must be at least zero.- n
Specifies the number of columns of
B
. The value ofn
must be at least zero.- alpha
Scaling factor for the solution.
- a
Buffer holding input matrix
A
. Must have size at leastlda
*m
ifleft_right
=side::left
, orlda
*n
ifleft_right
=side::right
. See Matrix Storage for more details.- lda
Leading dimension of
A
. Must be at leastm
ifleft_right
=side::left
, and at leastn
ifleft_right
=side::right
. Must be positive.- b
Buffer holding input/output matrix
B
. Must have size at leastldb
*n
if column major layout is used to store matrices or at leastldb
*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 leastm
if column major layout is used to store matrices or at leastn
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.
trsm (USM Version)¶
Syntax
namespace oneapi::mkl::blas::column_major {
sycl::event trsm(sycl::queue &queue,
onemkl::side left_right,
onemkl::uplo upper_lower,
onemkl::transpose transa,
onemkl::diag unit_diag,
std::int64_t m,
std::int64_t n,
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,
onemkl::side left_right,
onemkl::uplo upper_lower,
onemkl::transpose transa,
onemkl::diag unit_diag,
std::int64_t m,
std::int64_t n,
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
multipliesX
on the left (side::left
) or on the right (side::right
). See oneMKL defined datatypes for more details.- uplo
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 toA
. 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 ofm
must be at least zero.- n
Specifies the number of columns of
B
. The value ofn
must be at least zero.- alpha
Scaling factor for the solution.
- a
Pointer to input matrix
A
. Must have size at leastlda
*m
ifleft_right
=side::left
, orlda
*n
ifleft_right
=side::right
. See Matrix Storage for more details.- lda
Leading dimension of
A
. Must be at leastm
ifleft_right
=side::left
, and at leastn
ifleft_right
=side::right
. Must be positive.- b
Pointer to input/output matrix
B
. Must have size at leastldb
*n
if column major layout is used to store matrices or at leastldb
*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 leastm
if column major layout is used to store matrices or at leastn
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::unsupported_device
Parent topic: BLAS Level 3 Routines