Visible to Intel only — GUID: GUID-691EA08D-CBCC-4D30-9D0E-2DDBDBD66352
Visible to Intel only — GUID: GUID-691EA08D-CBCC-4D30-9D0E-2DDBDBD66352
omatadd
Computes a sum of two general matrices, with optional transposes.
Description
The omatadd routine performs an out-of-place scaled matrix addition with optional transposes in the arguments. The operation is defined as:
where:
op(X) is one of op(X) = X, or op(X) = XT, or op(X) = XH
alpha and beta are scalars
A, B and C are matrices
op(A) is m x n matrix
op(B) is m x n matrix
C is m x n matrix
In general, A, B, and C must not overlap in memory, with the exception of the following in-place operations:
A and C may point to the same memory if op(A) is non-transpose and lda = ldc;
B and C may point to the same memory if op(B) is non-transpose and ldb = ldc.
omatadd supports the following precisions:
T |
---|
float |
double |
std::complex<float> |
std::complex<double> |
omatadd (Buffer Version)
Syntax
namespace oneapi::mkl::blas::column_major {
void omatadd(sycl::queue &queue,
oneapi::mkl::transpose transa,
oneapi::mkl::transpose transb,
std::int64_t m,
std::int64_t n,
T alpha,
sycl::buffer<T, 1> &a,
std::int64_t lda,
T beta,
sycl::buffer<T, 1> &b,
std::int64_t ldb,
sycl::buffer<T, 1> &c,
std::int64_t ldc)
}
namespace oneapi::mkl::blas::row_major {
void omatadd(sycl::queue &queue,
oneapi::mkl::transpose transa,
oneapi::mkl::transpose transb,
std::int64_t m,
std::int64_t n,
T alpha,
sycl::buffer<T, 1> &a,
std::int64_t lda,
T beta,
sycl::buffer<T, 1> &b,
std::int64_t ldb,
sycl::buffer<T, 1> &c,
std::int64_t ldc)
}
Input Parameters
- queue
-
The queue where the routine will be executed.
- transa
-
Specifies op(A), the transposition operation applied to the matrix A.
- transb
-
Specifies op(B), the transposition operation applied to the matrix B.
- m
-
Number of rows for the result matrix C. Must be at least zero.
- n
-
Number of columns for the result matrix C. Must be at least zero.
- alpha
-
Scaling factor for the matrix A.
- a
-
Buffer holding the input matrix A.
transa = transpose::nontrans
transa = transpose::trans or transa = transpose::conjtrans
Column major
A is m x n matrix. Size of array a must be at least lda * n
A is n x m matrix. Size of array a must be at least lda * m
Row major
A is m x n matrix. Size of array a must be at least lda * m
A is n x m matrix. Size of array a must be at least lda * n
- lda
-
Leading dimension of matrix A. Must be positive
transa = transpose::nontrans
transa = transpose::trans or transa = transpose::conjtrans
Column major
Must be at least m
Must be at least n
Row major
Must be at least n
Must be at least m
- beta
-
Scaling factor for the matrix B.
- b
-
Buffer holding the input matrix B.
transb = transpose::nontrans
transb = transpose::trans or transb = transpose::conjtrans
Column major
B is m x n matrix. Size of array b must be at least ldb * n
B is n x m matrix. Size of array b must be at least ldb * m
Row major
B is m x n matrix. Size of array b must be at least ldb * m
B is n x m matrix. Size of array b must be at least ldb * n
- ldb
-
Leading dimension of marix B. Must be positive.
transb = transpose::nontrans
transb = transpose::trans or transb = transpose::conjtrans
Column major
Must be at least m
Must be at least n
Row major
Must be at least n
Must be at least m
- ldc
-
Leading dimension of the matrix C. If matrices are stored using column major layout, ldc must be at least m. If matrices are stored using row major layout, ldc must be at least n. Must be positive.
Output Parameters
- c
-
Output buffer overwritten by alpha * op(A) + beta * op(B).
Column major
C is m x n matrix. Size of array c must be at least ldc * n
Row major
C is m x n matrix. Size of array c must be at least ldc * m
omatadd (USM Version)
Syntax
namespace oneapi::mkl::blas::column_major {
sycl::event omatadd(sycl::queue &queue,
oneapi::mkl::transpose transa,
oneapi::mkl::transpose transb,
std::int64_t m,
std::int64_t n,
oneapi::mkl::value_or_pointer<T> alpha,
const T *a,
std::int64_t lda,
oneapi::mkl::value_or_pointer<T> beta,
const T *b,
std::int64_t ldb,
T *c,
std::int64_t ldc,
const std::vector<sycl::event> &dependencies = {});
}
namespace oneapi::mkl::blas::row_major {
sycl::event omatadd(sycl::queue &queue,
oneapi::mkl::transpose transa,
oneapi::mkl::transpose transb,
std::int64_t m,
std::int64_t n,
oneapi::mkl::value_or_pointer<T> alpha,
const T *a,
std::int64_t lda,
oneapi::mkl::value_or_pointer<T> beta,
const T *b,
std::int64_t ldb,
T *c,
std::int64_t ldc,
const std::vector<sycl::event> &dependencies = {});
}
Input Parameters
- queue
-
The queue where the routine will be executed.
- transa
-
Specifies op(A), the transposition operation applied to matrix A. See Data Types for more details.
- transb
-
Specifies op(B), the transposition operation applied to matrix B. See Data Types for more details.
- m
-
Number of rows for the result matrix C. Must be at least zero.
- n
-
Number of columns for the result matrix C. Must be at least zero.
- alpha
-
Scaling factor for the matrix A. See Scalar Arguments for more information on the value_or_pointer data type.
- a
-
Pointer to input matrix A. See Matrix Storage for more details.
transa = transpose::nontrans
transa = transpose::trans or transa = transpose::conjtrans
Column major
A is m x n matrix. Size of array a must be at least lda * n
A is n x m matrix. Size of array a must be at least lda * m
Row major
A is m x n matrix. Size of array a must be at least lda * m
A is n x m matrix. Size of array a must be at least lda * n
- lda
-
Leading dimension of matrix A. Must be positive
transa = transpose::nontrans
transa = transpose::trans or transa = transpose::conjtrans
Column major
Must be at least m
Must be at least n
Row major
Must be at least n
Must be at least m
- beta
-
Scaling factor for the matrix B. See Scalar Arguments for more information on the value_or_pointer data type.
- b
-
Pointer to input matrix B. See Matrix Storage for more details.
transb = transpose::nontrans
transb = transpose::trans or transb = transpose::conjtrans
Column major
B is m x n matrix. Size of array b must be at least ldb * n
B is n x m matrix. Size of array b must be at least ldb * m
Row major
B is m x n matrix. Size of array b must be at least ldb * m
B is n x m matrix. Size of array b must be at least ldb * n
- ldb
-
Leading dimension of marix B. Must be positive.
transb = transpose::nontrans
transb = transpose::trans or transb = transpose::conjtrans
Column major
Must be at least m
Must be at least n
Row major
Must be at least n
Must be at least m
- ldc
-
Leading dimension of the matrix C. If matrices are stored using column major layout, ldc must be at least m. If matrices are stored using row major layout, ldc must be at least n. 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 * op(A) + beta * op(B).
Column major
C is m x n matrix. Size of array c must be at least ldc * n
Row major
C is m x n matrix. Size of array c must be at least ldc * m
Return Values
Output event to wait for to ensure computation is complete.