hemm (USM Version)¶
Computes a matrix-matrix product where one input matrix is Hermitian and one is general.
Description¶
The hemm
routines compute a scalar-matrix-matrix product and add the
result to a scalar-matrix product, where one of the matrices in the
multiplication is Hermitian. The argument left_right
determines
if the Hermitian matrix, A
, is on the left of the multiplication
(left_right
= side::left
) or on the right (left_right
=
side::right
). Depending on left_right
, the operation is
defined as
or
where:
alpha
andbeta
are scalars,A
is a Hermitian matrix, eitherm
-by-m
orn
-by-n
matrices,B
andC
arem
-by-n
matrices.
API¶
Syntax¶
namespace oneapi::mkl::blas::column_major {
sycl::event hemm(sycl::queue &queue,
onemkl::side left_right,
onemkl::uplo upper_lower,
std::int64_t m,
std::int64_t n,
T alpha,
const T* a,
std::int64_t lda,
const T* b,
std::int64_t ldb,
T beta,
T* c,
std::int64_t ldc,
const sycl::vector_class<sycl::event> &dependencies = {})
}
namespace oneapi::mkl::blas::row_major {
sycl::event hemm(sycl::queue &queue,
onemkl::side left_right,
onemkl::uplo upper_lower,
std::int64_t m,
std::int64_t n,
T alpha,
const T* a,
std::int64_t lda,
const T* b,
std::int64_t ldb,
T beta,
T* c,
std::int64_t ldc,
const sycl::vector_class<sycl::event> &dependencies = {})
}
The USM version of hemm
supports the following precisions and devices:
T |
Devices Supported |
---|---|
|
Host, CPU, and GPU |
|
Host, CPU, and GPU |
Input Parameters¶
- exec_queue
The queue where the routine should be executed.
- left_right
Specifies whether
A
is on the left side of the multiplication (side::left
) or on the right side (side::right
). See Data Types for more details.- uplo
Specifies whether
A
’s data is stored in its upper or lower triangle. See Data Types for more details.- m
Specifies the number of rows of the matrix
B
andC
.The value of
m
must be at least zero.- n
Specifies the number of columns of the matrix
B
andC
.The value of
n
must be at least zero.- alpha
Scaling factor for the matrix-matrix product.
- a
Pointer to input matrix
A
. Must have size at leastlda
*m
ifA
is on the left of the multiplication, orlda
*n
ifA
is on the right. See Matrix Storage for more details.- lda
Leading dimension of
A
. Must be at leastm
ifA
is on the left of the multiplication, or at leastn
ifA
is on the right. Must be positive.- b
Pointer to input matrix
B
. It must have a size of 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 column major layout is used to store matrices.- beta
Scaling factor for matrix
C
.- c
Pointer to input/output matrix
C
. It must have a size of at leastldc*n
if column major layout is used to store matrices or at leastldc*m
if row major layout is used to store matrices. See Matrix Storage for more details.- ldc
Leading dimension of
C
. It must be positive and at leastm
if column major layout is used to store matrices or at leastn
if column 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¶
- c
Pointer to the output matrix, overwritten by
alpha
*A
*B
+beta
*C
(left_right
=side::left
) oralpha
*B
*A
+beta
*C
(left_right
=side::right
).
Note
If beta
= 0, matrix C
does not need to be initialized before calling hemm
.
Return Values¶
Output event to wait on to ensure computation is complete.