hemm¶
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
). The operation is defined as:
If (left_right
= side::left
),
If (left_right
= side::right
),
where:
alpha
andbeta
are scalarsA
is eitherm
xm
orn
xn
hermitian matrixB
andC
arem
xn
matrices
hemm
supports the following precisions:
T |
---|
|
|
hemm (Buffer Version)¶
Syntax¶
namespace oneapi::mkl::blas::column_major {
void hemm(sycl::queue &queue,
oneapi::mkl::side left_right,
oneapi::mkl::uplo upper_lower,
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 hemm(sycl::queue &queue,
oneapi::mkl::side left_right,
oneapi::mkl::uplo upper_lower,
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 matrix
A
is on the left side or right side of the multiplication. See Data Types for more details.- upper_lower
Specifies whether matrix
A
is upper or lower triangular. See Data Types for more details.- m
Number of rows of matrix
B
and matrixC
. Must be at least zero.- n
Number of columns of matrix
B
and matrixC
. Must be at least zero.- alpha
Scaling factor for matrix-matrix product.
- a
Buffer holding input matrix
A
. Size of the buffer must be at leastlda
*m
ifleft_right
=side::left
orlda
*n
ifleft_right
=side::right
. See Matrix Storage for more details.- lda
Leading dimension of matrix
A
. Must be at leastm
ifleft_right
=side::left
or at leastn
ifleft_right
=side::right
. Must be positive.- b
Buffer holding input matrix
B
. Size of the buffer must be at leastldb
*n
if column major layout or at leastldb
*m
if row major layout is used. See Matrix Storage for more details.- ldb
Leading dimension of matrix
B
. Must be at leastm
if column major layout or at leastn
if row major layout is used. Must be positive.- beta
Scaling factor for matrix
C
.- c
Buffer holding input/output matrix
C
. Size of the buffer must be at leastldc
*n
if column major layout or at leastldc
*m
if row major layout is used. See Matrix Storage for more details.- ldc
Leading dimension of matrix
C
. Must be at leastm
if column major layout or at leastn
if row major layout is used. Must be positive.
Output Parameters¶
- c
Output buffer overwritten by
alpha
*A
*B
+beta
*C
ifleft_right
=side::left
oralpha
*B
*A
+beta
*C
ifleft_right
=side::right
.
Note
If beta
= 0, matrix C
does not need to be initialized before calling hemm
.
hemm (Buffer Version)¶
Syntax¶
namespace oneapi::mkl::blas::column_major {
sycl::event hemm(sycl::queue &queue,
oneapi::mkl::side left_right,
oneapi::mkl::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 std::vector<sycl::event> &dependencies = {})
}
namespace oneapi::mkl::blas::row_major {
sycl::event hemm(sycl::queue &queue,
oneapi::mkl::side left_right,
oneapi::mkl::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 std::vector<sycl::event> &dependencies = {})
}
Input Parameters¶
- queue
The queue where the routine should be executed.
- left_right
Specifies whether matrix
A
is on the left side or right side of the multiplication. See Data Types for more details.- upper_lower
Specifies whether matrix
A
is upper or lower triangular. See Data Types for more details.- m
Number of rows of matrix
B
and matrixC
. Must be at least zero.- n
Number of columns of matrix
B
and matrixC
. Must be at least zero.- alpha
Scaling factor for matrix-matrix product.
- a
Pointer to input matrix
A
. Size of the array must be at leastlda
*m
ifleft_right
=side::left
orlda
*n
ifleft_right
=side::right
. See Matrix Storage for more details.- lda
Leading dimension of matrix
A
. Must be at leastm
ifleft_right
=side::left
or at leastn
ifleft_right
=side::right
. Must be positive.- b
Pointer to input matrix
B
. Size of the array must be at leastldb
*n
if column major layout or at leastldb
*m
if row major layout is used. See Matrix Storage for more details.- ldb
Leading dimension of matrix
B
. Must be at leastm
if column major layout or at leastn
if row major layout is used. Must be positive.- beta
Scaling factor for matrix
C
.- c
Pointer to input/output matrix
C
. Size of the array must be at leastldc
*n
if column major layout or at leastldc
*m
if row major layout is used. See Matrix Storage for more details.- ldc
Leading dimension of matrix
C
. Must be at leastm
if column major layout or at leastn
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 output matrix overwritten by
alpha
*A
*B
+beta
*C
ifleft_right
=side::left
oralpha
*B
*A
+beta
*C
ifleft_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.