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
) =X
T, or op(X
) =X
Halpha
andbeta
are scalarsA
,B
andC
are matricesop(
A
) ism
xn
matrixop(
B
) ism
xn
matrixC
ism
xn
matrix
In general, A
, B
, and C
must not overlap in memory, with the exception of
the following in-place operations:
A
andC
may point to the same memory ifop(A)
is non-transpose andlda
=ldc
;
B
andC
may point to the same memory ifop(B)
is non-transpose andldb
=ldc
.
omatadd
supports the following precisions:
T |
---|
|
|
|
|
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
ortransa
=transpose::conjtrans
Column major
A
ism
xn
matrix. Size of arraya
must be at leastlda
*n
A
isn
xm
matrix. Size of arraya
must be at leastlda
*m
Row major
A
ism
xn
matrix. Size of arraya
must be at leastlda
*m
A
isn
xm
matrix. Size of arraya
must be at leastlda
*n
- lda
Leading dimension of matrix
A
. Must be positivetransa
=transpose::nontrans
transa
=transpose::trans
ortransa
=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
ortransb
=transpose::conjtrans
Column major
B
ism
xn
matrix. Size of arrayb
must be at leastldb
*n
B
isn
xm
matrix. Size of arrayb
must be at leastldb
*m
Row major
B
ism
xn
matrix. Size of arrayb
must be at leastldb
*m
B
isn
xm
matrix. Size of arrayb
must be at leastldb
*n
- ldb
Leading dimension of marix
B
. Must be positive.transb
=transpose::nontrans
transb
=transpose::trans
ortransb
=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 leastm
. If matrices are stored using row major layout,ldc
must be at leastn
. Must be positive.
Output Parameters¶
- c
Output buffer overwritten by
alpha
* op(A
) +beta
* op(B
).Column major
C
ism
xn
matrix. Size of arrayc
must be at leastldc
*n
Row major
C
ism
xn
matrix. Size of arrayc
must be at leastldc
*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,
T alpha,
const T *a,
std::int64_t lda,
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,
T alpha,
const T *a,
std::int64_t lda,
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 matrixA
. See Data Types for more details.- transb
Specifies op(
B
), the transposition operation applied to matrixB
. 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
.- a
Pointer to input matrix
A
. See Matrix Storage for more details.transa
=transpose::nontrans
transa
=transpose::trans
ortransa
=transpose::conjtrans
Column major
A
ism
xn
matrix. Size of arraya
must be at leastlda
*n
A
isn
xm
matrix. Size of arraya
must be at leastlda
*m
Row major
A
ism
xn
matrix. Size of arraya
must be at leastlda
*m
A
isn
xm
matrix. Size of arraya
must be at leastlda
*n
- lda
Leading dimension of matrix
A
. Must be positivetransa
=transpose::nontrans
transa
=transpose::trans
ortransa
=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
Pointer to input matrix
B
. See Matrix Storage for more details.transb
=transpose::nontrans
transb
=transpose::trans
ortransb
=transpose::conjtrans
Column major
B
ism
xn
matrix. Size of arrayb
must be at leastldb
*n
B
isn
xm
matrix. Size of arrayb
must be at leastldb
*m
Row major
B
ism
xn
matrix. Size of arrayb
must be at leastldb
*m
B
isn
xm
matrix. Size of arrayb
must be at leastldb
*n
- ldb
Leading dimension of marix
B
. Must be positive.transb
=transpose::nontrans
transb
=transpose::trans
ortransb
=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 leastm
. If matrices are stored using row major layout,ldc
must be at leastn
. 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
ism
xn
matrix. Size of arrayc
must be at leastldc
*n
Row major
C
ism
xn
matrix. Size of arrayc
must be at leastldc
*m
Return Values¶
Output event to wait for to ensure computation is complete.