gemmt#
Computes a matrix-matrix product with general matrices, but updates only the upper or lower triangular part of the result matrix.
Description
The gemmt routines compute a scalar-matrix-matrix product and add the result to the upper or lower part of a scalar-matrix product, with general matrices. 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 n x k, op(B) is k x n, and
C is n x n.
gemmt supports the following precisions.
T
float
double
std::complex<float>
std::complex<double>
gemmt (Buffer Version)#
Syntax
namespace oneapi::mkl::blas::column_major {
void gemmt(sycl::queue &queue,
onemkl::uplo upper_lower,
onemkl::transpose transa,
onemkl::transpose transb,
std::int64_t n,
std::int64_t k,
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 gemmt(sycl::queue &queue,
onemkl::uplo upper_lower,
onemkl::transpose transa,
onemkl::transpose transb,
std::int64_t n,
std::int64_t k,
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.
- upper_lower
Specifies whether
C’s data is stored in its upper or lower triangle. See oneMKL Defined Datatypes for more details.- transa
Specifies op(
A), the transposition operation applied toA. See oneMKL Defined Datatypes for more details.- transb
Specifies op(
B), the transposition operation applied toB. See oneMKL Defined Datatypes for more details.- n
Number of rows of op(
A), columns of op(B), and columns and rows ofC. Must be at least zero.- k
Number of columns of op(
A) and rows of op(B). Must be at least zero.- alpha
Scaling factor for the matrix-matrix product.
- a
Buffer holding the input matrix
A.Anot transposedAtransposedColumn major
Ais ann-by-kmatrix so the arrayamust have size at leastlda*k.Ais ank-by-nmatrix so the arrayamust have size at leastlda*nRow major
Ais ann-by-kmatrix so the arrayamust have size at leastlda*n.Ais ank-by-nmatrix so the arrayamust have size at leastlda*k.See Matrix Storage for more details.
- lda
The leading dimension of
A. It must be positive.Anot transposedAtransposedColumn major
ldamust be at leastn.ldamust be at leastk.Row major
ldamust be at leastk.ldamust be at leastn.- b
Buffer holding the input matrix
B.Bnot transposedBtransposedColumn major
Bis ank-by-nmatrix so the arraybmust have size at leastldb*n.Bis ann-by-kmatrix so the arraybmust have size at leastldb*kRow major
Bis ank-by-nmatrix so the arraybmust have size at leastldb*k.Bis ann-by-kmatrix so the arraybmust have size at leastldb*n.See Matrix Storage for more details.
- ldb
The leading dimension of
B. It must be positive.Bnot transposedBtransposedColumn major
ldbmust be at leastk.ldbmust be at leastn.Row major
ldbmust be at leastn.ldbmust be at leastk.- beta
Scaling factor for matrix
C.- c
Buffer holding the input/output matrix
C. Must have size at leastldc*n. See Matrix Storage for more details.- ldc
Leading dimension of
C. Must be positive and at leastm.
Output Parameters
- c
Output buffer, overwritten by the upper or lower triangular part of
alpha* op(A)*op(B) +beta*C.
Notes
If beta = 0, matrix C does not need to be initialized
before calling gemmt.
gemmt (USM Version)#
Syntax
namespace oneapi::mkl::blas::column_major {
sycl::event gemmt(sycl::queue &queue,
onemkl::uplo upper_lower,
onemkl::transpose transa,
onemkl::transpose transb,
std::int64_t n,
std::int64_t k,
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 gemmt(sycl::queue &queue,
onemkl::uplo upper_lower,
onemkl::transpose transa,
onemkl::transpose transb,
std::int64_t n,
std::int64_t k,
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.
- upper_lower
Specifies whether
C’s data is stored in its upper or lower triangle. See oneMKL Defined Datatypes for more details.- transa
Specifies op(
A), the transposition operation applied toA. See oneMKL Defined Datatypes for more details.- transb
Specifies op(
B), the transposition operation applied toB. See oneMKL Defined Datatypes for more details.- n
Number of columns of op(
A), columns of op(B), and columns ofC. Must be at least zero.- k
Number of columns of op(
A) and rows of op(B). Must be at least zero.- alpha
Scaling factor for the matrix-matrix product.
- a
Pointer to input matrix
A.Anot transposedAtransposedColumn major
Ais ann-by-kmatrix so the arrayamust have size at leastlda*k.Ais ank-by-nmatrix so the arrayamust have size at leastlda*nRow major
Ais ann-by-kmatrix so the arrayamust have size at leastlda*n.Ais ank-by-nmatrix so the arrayamust have size at leastlda*kSee Matrix Storage for more details.
- lda
The leading dimension of
A. It must be positive.Anot transposedAtransposedColumn major
ldamust be at leastn.ldamust be at leastk.Row major
ldamust be at leastk.ldamust be at leastn.- b
Pointer to input matrix
B.Bnot transposedBtransposedColumn major
Bis ank-by-nmatrix so the arraybmust have size at leastldb*n.Bis ann-by-kmatrix so the arraybmust have size at leastldb*kRow major
Bis ank-by-nmatrix so the arraybmust have size at leastldb*k.Bis ann-by-kmatrix so the arraybmust have size at leastldb*nSee Matrix Storage for more details.
- ldb
The leading dimension of
B. It must be positive.Bnot transposedBtransposedColumn major
ldbmust be at leastk.ldbmust be at leastn.Row major
ldbmust be at leastn.ldbmust be at leastk.- beta
Scaling factor for matrix
C.- c
Pointer to input/output matrix
C. Must have size at leastldc*n. See Matrix Storage for more details.- ldc
Leading dimension of
C. Must be positive and at leastm.- 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 the upper or lower triangular part of
alpha* op(A)*op(B) +beta*C.
Notes
If beta = 0, matrix C does not need to be initialized
before calling gemmt.
Return Values
Output event to wait on to ensure computation is complete.
Parent topic: BLAS-like Extensions