gemm_batch#
Computes a group of gemm operations.
Description
The gemm_batch routines are batched versions of gemm, performing
multiple gemm operations in a single call. Each gemm
operation perform a matrix-matrix product with general matrices.
gemm_batch supports the following precisions.
T
half
float
double
std::complex<float>
std::complex<double>
gemm_batch (Buffer Version)#
Description
The buffer version of gemm_batch supports only the strided API.
The strided API operation is defined as:
for i = 0 … batch_size – 1
A, B and C are matrices at offset i * stridea, i * strideb, i * stridec in a, b and c.
C := alpha * op(A) * op(B) + beta * C
end for
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 k, op(B) is
k x n, and C is m x n.
The a, b and c buffers contain all the input matrices. The stride
between matrices is given by the stride parameter. The total number
of matrices in a, b and c buffers is given by the batch_size parameter.
Strided API
Syntax
namespace oneapi::mkl::blas::column_major {
void gemm_batch(sycl::queue &queue,
onemkl::transpose transa,
onemkl::transpose transb,
std::int64_t m,
std::int64_t n,
std::int64_t k,
T alpha,
sycl::buffer<T,1> &a,
std::int64_t lda,
std::int64_t stridea,
sycl::buffer<T,1> &b,
std::int64_t ldb,
std::int64_t strideb,
T beta,
sycl::buffer<T,1> &c,
std::int64_t ldc,
std::int64_t stridec,
std::int64_t batch_size)
}
namespace oneapi::mkl::blas::row_major {
void gemm_batch(sycl::queue &queue,
onemkl::transpose transa,
onemkl::transpose transb,
std::int64_t m,
std::int64_t n,
std::int64_t k,
T alpha,
sycl::buffer<T,1> &a,
std::int64_t lda,
std::int64_t stridea,
sycl::buffer<T,1> &b,
std::int64_t ldb,
std::int64_t strideb,
T beta,
sycl::buffer<T,1> &c,
std::int64_t ldc,
std::int64_t stridec,
std::int64_t batch_size)
}
Input Parameters
- queue
The queue where the routine should be executed.
- transa
Specifies op(
A) the transposition operation applied to the matricesA. See oneMKL Defined Datatypes for more details.- transb
Specifies op(
B) the transposition operation applied to the matricesB. See oneMKL Defined Datatypes for more details.- m
Number of rows of op(
A) andC. Must be at least zero.- n
Number of columns of op(
B) andC. 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 products.
- a
Buffer holding the input matrices
Awith sizestridea*batch_size.- lda
The leading dimension of the matrices
A. It must be positive.Anot transposedAtransposedColumn major
ldamust be at leastm.ldamust be at leastk.Row major
ldamust be at leastk.ldamust be at leastm.- stridea
Stride between different
Amatrices.- b
Buffer holding the input matrices
Bwith sizestrideb*batch_size.- ldb
The leading dimension of the matrices``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.- strideb
Stride between different
Bmatrices.- beta
Scaling factor for the matrices
C.- c
Buffer holding input/output matrices
Cwith sizestridec*batch_size.- ldc
The leading dimension of the matrices
C. It must be positive and at leastmif column major layout is used to store matrices or at leastnif row major layout is used to store matrices.- stridec
Stride between different
Cmatrices. Must be at leastldc*n.- batch_size
Specifies the number of matrix multiply operations to perform.
Output Parameters
- c
Output buffer, overwritten by
batch_sizematrix multiply operations of the formalpha* op(A)*op(B) +beta*C.
Notes
If beta = 0, matrix C does not need to be initialized before
calling gemm_batch.
gemm_batch (USM Version)#
Description
The USM version of gemm_batch supports the group API and strided API.
The group API operation is defined as:
idx = 0
for i = 0 … group_count – 1
for j = 0 … group_size – 1
A, B, and C are matrices in a[idx], b[idx] and c[idx]
C := alpha[i] * op(A) * op(B) + beta[i] * C
idx = idx + 1
end for
end for
The strided API operation is defined as
for i = 0 … batch_size – 1
A, B and C are matrices at offset i * stridea, i * strideb, i * stridec in a, b and c.
C := alpha * op(A) * op(B) + beta * C
end for
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 k, op(B) is k x n, and C is m x n.
For group API, a, b and c arrays contain the pointers for all the input matrices.
The total number of matrices in a, b and c are given by:
For strided API, a, b, c arrays contain all the input matrices. The total number of matrices
in a, b and c are given by the batch_size parameter.
Group API
Syntax
namespace oneapi::mkl::blas::column_major {
sycl::event gemm_batch(sycl::queue &queue,
onemkl::transpose *transa,
onemkl::transpose *transb,
std::int64_t *m,
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,
std::int64_t group_count,
std::int64_t *group_size,
const std::vector<sycl::event> &dependencies = {})
}
namespace oneapi::mkl::blas::row_major {
sycl::event gemm_batch(sycl::queue &queue,
onemkl::transpose *transa,
onemkl::transpose *transb,
std::int64_t *m,
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,
std::int64_t group_count,
std::int64_t *group_size,
const std::vector<sycl::event> &dependencies = {})
}
Input Parameters
- queue
The queue where the routine should be executed.
- transa
Array of
group_countonemkl::transposevalues.transa[i]specifies the form of op(A) used in the matrix multiplication in groupi. See oneMKL Defined Datatypes for more details.- transb
Array of
group_countonemkl::transposevalues.transb[i]specifies the form of op(B) used in the matrix multiplication in groupi. See oneMKL Defined Datatypes for more details.- m
Array of
group_countintegers.m[i]specifies the number of rows of op(A) andCfor every matrix in groupi. All entries must be at least zero.- n
Array of
group_countintegers.n[i]specifies the number of columns of op(B) andCfor every matrix in groupi. All entries must be at least zero.- k
Array of
group_countintegers.k[i]specifies the number of columns of op(A) and rows of op(B) for every matrix in groupi. All entries must be at least zero.- alpha
Array of
group_countscalar elements.alpha[i]specifies the scaling factor for every matrix-matrix product in groupi.- a
Array of pointers to input matrices
Awith sizetotal_batch_count.See Matrix Storage for more details.
- lda
Array of
group_countintegers.lda[i]specifies the leading dimension ofAfor every matrix in groupi. All entries must be positive.Anot transposedAtransposedColumn major
lda[i]must be at leastm[i].lda[i]must be at leastk[i].Row major
lda[i]must be at leastk[i].lda[i]must be at leastm[i].- b
Array of pointers to input matrices
Bwith sizetotal_batch_count.See Matrix Storage for more details.
- ldb
Array of
group_countintegers.ldb[i]specifies the leading dimension ofBfor every matrix in groupi. All entries must be positive.Bnot transposedBtransposedColumn major
ldb[i]must be at leastk[i].ldb[i]must be at leastn[i].Row major
ldb[i]must be at leastn[i].ldb[i]must be at leastk[i].- beta
Array of
group_countscalar elements.beta[i]specifies the scaling factor for matrixCfor every matrix in groupi.- c
Array of pointers to input/output matrices
Cwith sizetotal_batch_count.See Matrix Storage for more details.
- ldc
Array of
group_countintegers.ldc[i]specifies the leading dimension ofCfor every matrix in groupi. All entries must be positive andldc[i]must be at leastm[i]if column major layout is used to store matrices or at leastn[i]if row major layout is used to store matrices.- group_count
Specifies the number of groups. Must be at least 0.
- group_size
Array of
group_countintegers.group_size[i]specifies the number of matrix multiply products in groupi. All entries must be at least 0.- dependencies
List of events to wait for before starting computation, if any. If omitted, defaults to no dependencies.
Output Parameters
- c
Overwritten by the
m[i]-by-n[i]matrix calculated by (alpha[i]* op(A)*op(B) +beta[i]*C) for groupi.
Notes
If beta = 0, matrix C does not need to be initialized
before calling gemm_batch.
Return Values
Output event to wait on to ensure computation is complete.
Strided API
Syntax
namespace oneapi::mkl::blas::column_major {
sycl::event gemm_batch(sycl::queue &queue,
onemkl::transpose transa,
onemkl::transpose transb,
std::int64_t m,
std::int64_t n,
std::int64_t k,
T alpha,
const T *a,
std::int64_t lda,
std::int64_t stridea,
const T *b,
std::int64_t ldb,
std::int64_t strideb,
T beta,
T *c,
std::int64_t ldc,
std::int64_t stridec,
std::int64_t batch_size,
const std::vector<sycl::event> &dependencies = {})
}
namespace oneapi::mkl::blas::row_major {
sycl::event gemm_batch(sycl::queue &queue,
onemkl::transpose transa,
onemkl::transpose transb,
std::int64_t m,
std::int64_t n,
std::int64_t k,
T alpha,
const T *a,
std::int64_t lda,
std::int64_t stridea,
const T *b,
std::int64_t ldb,
std::int64_t strideb,
T beta,
T *c,
std::int64_t ldc,
std::int64_t stridec,
std::int64_t batch_size,
const std::vector<sycl::event> &dependencies = {})
}
Input Parameters
- queue
The queue where the routine should be executed.
- transa
Specifies op(
A) the transposition operation applied to the matricesA. See oneMKL Defined Datatypes for more details.- transb
Specifies op(
B) the transposition operation applied to the matricesB. See oneMKL Defined Datatypes for more details.- m
Number of rows of op(
A) andC. Must be at least zero.- n
Number of columns of op(
B) andC. 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 products.
- a
Pointer to input matrices
Awith sizestridea*batch_size.- lda
The leading dimension of the matrices
A. It must be positive.Anot transposedAtransposedColumn major
ldamust be at leastm.ldamust be at leastk.Row major
ldamust be at leastk.ldamust be at leastm.- stridea
Stride between different
Amatrices.- b
Pointer to input matrices
Bwith sizestrideb*batch_size.- ldb
The leading dimension of the matrices``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.- strideb
Stride between different
Bmatrices.- beta
Scaling factor for the matrices
C.- c
Pointer to input/output matrices
Cwith sizestridec*batch_size.- ldc
The leading dimension of the matrices
C. It must be positive and at leastmif column major layout is used to store matrices or at leastnif row major layout is used to store matrices.- stridec
Stride between different
Cmatrices.- batch_size
Specifies the number of matrix multiply operations to perform.
- dependencies
List of events to wait for before starting computation, if any. If omitted, defaults to no dependencies.
Output Parameters
- c
Output matrices, overwritten by
batch_sizematrix multiply operations of the formalpha* op(A)*op(B) +beta*C.
Notes
If beta = 0, matrix C does not need to be initialized before
calling gemm_batch.
Return Values
Output event to wait on to ensure computation is complete.
Parent topic: BLAS-like Extensions