Visible to Intel only — GUID: GUID-2137C746-4A4B-4AB3-B7B8-0E84C4490306
Visible to Intel only — GUID: GUID-2137C746-4A4B-4AB3-B7B8-0E84C4490306
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 and beta are scalars
A is either m x m or n x n hermitian matrix
B and C are m x n matrices
hemm supports the following precisions:
T |
---|
std::complex<float> |
std::complex<double> |
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, compute_mode mode = compute_mode::unset) }
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, compute_mode mode = compute_mode::unset) }
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 matrix C. Must be at least zero.
- n
-
Number of columns of matrix B and matrix C. 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 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 matrix 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. Size of the buffer must be 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. Must be at least m if column major layout or at least n 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 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.
- mode
-
Optional. Compute mode settings. See Compute Modes for more details.
Output Parameters
- c
-
Output buffer overwritten by alpha * A * B + beta * C if left_right = side::left or alpha * B * A + beta * C if left_right = side::right.
hemm (USM 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, compute_mode mode = compute_mode::unset, 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, compute_mode mode = compute_mode::unset, 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 matrix C. Must be at least zero.
- n
-
Number of columns of matrix B and matrix C. 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 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 matrix 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. Size of the array must be 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. Must be at least m if column major layout or at least n 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 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.
- mode
-
Optional. Compute mode settings. See Compute Modes for more details.
- dependencies
-
Optional. List of events to wait for before starting computation, if any. If omitted, defaults to no dependencies.
mode and dependencies may be omitted independently; it is not necessary to specify mode in order to provide dependencies.
Output Parameters
- c
-
Pointer to output matrix overwritten by alpha * A * B + beta * C if left_right = side::left or alpha * B * A + beta * C if left_right = side::right.
Return Values
Output event to wait on to ensure computation is complete.