gemm_batch#

一般行列の行列-行列積のグループを計算します。

説明

gemm_batch ルーチンは gemm のバッチバージョンであり、1 回の呼び出しで複数の gemm 操作を実行します。各 gemm 演算は、一般行列の行列 - 行列積を計算します。

gemm_batch は次の精度をサポートします。

Ta
(A 行列)
Tb
(B 行列)
Tc
(C 行列)
Ts
(alpha/beta)

sycl::half

sycl::half

sycl::half

sycl::half

sycl::half

sycl::half

float

float

oneapi::mkl::bfloat16

oneapi::mkl::bfloat16

oneapi::mkl::bfloat16

float

oneapi::mkl::bfloat16

oneapi::mkl::bfloat16

float

float

std::int8_t

std::int8_t

std::int32_t

float

std::int8_t

std::int8_t

float

float

float

float

float

float

double

double

double

double

std::complex<float>

std::complex<float>

std::complex<float>

std::complex<float>

std::complex<double>

std::complex<double>

std::complex<double>

std::complex<double>

gemm_batch (バッファーバージョン)#

gemm_batch のバッファーバージョンは、ストライド API のみをサポートします。

ストライド API#

ストライド API 操作は次のように定義されます。

for i = 0 … batch_size – 1 
    A、B、C は、a、b、c 内のオフセット i * stridea、i * strideb、i * stridec にある行列です。
    C = alpha * op(A) * op(B) + beta * C 
end for

説明:

  • op(X) は、op(X) = X、op(X) = XT、または op(X) = XH のいずれかです

  • alphabeta はスカラーです

  • ABC は行列です

  • op(A) は m x k、op(B)は k x nCm x n です。

ストライド API の場合、ab および c バッファーにすべての入力行列が含まれます。行列間のストライドは、ストライド・パラメーターで指定されます。abc バッファー内の行列の合計数は、batch_size パラメーターで指定されます。

構文#

namespace oneapi::mkl::blas::column_major { 
    void gemm_batch(sycl::queue &queue, 
                    oneapi::mkl::transpose transa, 
                    oneapi::mkl::transpose transb, 
                    std::int64_t m, 
                    std::int64_t n, 
                    std::int64_t k, 
                    Ts alpha, 
                    sycl::buffer<Ta,1> &a, 
                    std::int64_t lda, 
                    std::int64_t stridea, 
                    sycl::buffer<Tb,1> &b, 
                    std::int64_t ldb, 
                    std::int64_t strideb, Ts beta, 
                    sycl::buffer<Tc,1> &c, 
                    std::int64_t ldc, 
                    std::int64_t stridec, 
                    std::int64_t batch_size, 
                    compute_mode mode = compute_mode::unset) 
}
namespace oneapi::mkl::blas::row_major { 
    void gemm_batch(sycl::queue &queue, 
                    oneapi::mkl::transpose transa, 
                    oneapi::mkl::transpose transb, 
                    std::int64_t m, 
                    std::int64_t n, 
                    std::int64_t k, 
                    Ts alpha, 
                    sycl::buffer<Ta,1> &a, 
                    std::int64_t lda, 
                    std::int64_t stridea, 
                    sycl::buffer<Tb,1> &b, 
                    std::int64_t ldb, 
                    std::int64_t strideb, 
                    Ts beta, 
                    sycl::buffer<Tc,1> &c, 
                    std::int64_t ldc, 
                    std::int64_t stridec, 
                    std::int64_t batch_size, 
                    compute_mode mode = compute_mode::unset) 
}

入力パラメーター#

queue

ルーチンを実行するキュー。

transa

行列 A に適用される転置演算 op(A) を指定します。詳細はデータタイプを参照してください。

transb

行列 B に適用される転置演算 op(B) を指定します。詳細はデータタイプを参照してください。

m

行列 op(A) と行列 C の行数。最小値はゼロです。

n

行列 op(B) と行列 C の列数。最小値はゼロです。

k

行列 op(A) の列数と行列 op(B) の行数。最小値はゼロです。

alpha

行列 - 行列積のスケーリング係数。

a

入力行列 A を保持するバッファー。バッファーのサイズは stridea * batch_size 以上である必要があります。

lda

行列 A の先頭次元。正である必要があります。

transa = transpose::nontrans

transa = transpose::trans または trans = transpose::conjtrans

列優先

m 以上である必要があります

k 以上である必要があります

行優先

k 以上である必要があります

m 以上である必要があります

stridea

連続する 2 つの行列 A 間のストライド。

transa = transpose::nontrans

transa = transpose::trans または trans = transpose::conjtrans

列優先

lda * k 以上である必要があります

lda * m 以上である必要があります

行優先

lda * m 以上である必要があります

lda * k 以上である必要があります

b

入力行列 B を保持するバッファー。バッファーのサイズは strideb * batch_size 以上である必要があります。

ldb

行列 B の先頭次元。正である必要があります。

transa = transpose::nontrans

transa = transpose::trans または trans = transpose::conjtrans

列優先

k 以上である必要があります

n 以上である必要があります

行優先

n 以上である必要があります

k 以上である必要があります

strideb

連続する 2 つの行列 B 間のストライド。

transa = transpose::nontrans

transa = transpose::trans または trans = transpose::conjtrans

列優先

ldb * n 以上である必要があります

ldb * k 以上である必要があります

行優先

ldb * k 以上である必要があります

ldb * n 以上である必要があります

beta

行列 C のスケーリング係数。

c

入力/出力行列 C を保持するバッファー。バッファーのサイズは stridec * batch_size 以上である必要があります。

ldc

行列 C の先頭次元。正である必要があります。

列優先

m 以上である必要があります

行優先

n 以上である必要があります

stridec

連続する 2 つの行列 C 間のストライド。

列優先

ldc * n 以上である必要があります

行優先

ldc * m 以上である必要があります

batch_size

実行する行列乗算演算の数を指定します。最小値はゼロです。

mode

オプション。計算モードの設定。計算モードを参照してください。

出力パラメーター#

c

出力バッファーは、alpha * op(A) * op(B) + beta * C 形式の batch_size gemm 操作によって上書きされます。

beta = 0 の場合、gemm_batch を呼び出す前に行列 C を初期化する必要はありません。

gemm_batch (USM バージョン)#

gemm_batch の USM バージョンは、グループ API とストライド API をサポートします。

グループ API#

グループ API の整数ポインターの Ti タイプは、std::int64_t または std::int32_t のいずれかです。

グループ API 操作は次のように定義されます。

idx = 0 
    for i = 0 … group_count – 1 
        for j = 0 … group_size – 1 
            A、B、C は a[idx]、b[idx]、c[idx] の行列です。
            C = alpha[i] * op(A) * op(B) + beta[i] * C idx = idx + 1 
        end for 
    end for

説明:

  • op(X) は、op(X) = X、op(X) = XT、または op(X) = XH のいずれかです

  • alphabeta はスカラーです

  • AB、および C は行列です

  • op(A) は m x k、op(B) は k x n、そして Cm x n です。

グループ API の場合、abc 配列にはすべての入力行列へのポインターが含まれます。ab および c の行列の合計数は次のように与えられます。

\[total\_batch\_count = \sum_{i=0}^{group\_count-1}group\_size[i]\]

グループ API は、ポインター入力とスパン入力の両方をサポートします。

ポインターの代わりにスパンを使用する利点は、配列のサイズを可変にでき、実行時にスパンのサイズを照会できることです。出力行列を除く各 gemm パラメーターでは、スパンのサイズは 1、グループの数、または合計バッチサイズにすることができます。出力行列の場合、すべての計算が独立していることを保証するには、スパンのサイズが合計バッチサイズである必要があります。

スパンのサイズに応じて、gemm 計算の各パラメーターは次のように使用されます。

span size = 1

パラメーターはすべての gemm 計算で再利用されます

span size = group_count

パラメーターはすべての gemm グループ内で再利用されます。各グループでこのパラメーターの値が異なります。これはポインター付きの gemm_batch グループ API と同じです。

span size = total_batch_count

gemm 計算では、このパラメーターに異なる値が使用されます。

構文#

namespace oneapi::mkl::blas::column_major { 
    sycl::event gemm_batch(sycl::queue &queue, 
                           const oneapi::mkl::transpose *transa, 
                           const oneapi::mkl::transpose *transb, 
                           const Ti *m, 
                           const Ti *n, 
                           const Ti *k, 
                           const Ts *alpha, 
                           const Ta **a, 
                           const Ti *lda, 
                           const Tb **b, 
                           const Ti *ldb, 
                           const Ts *beta, 
                           Tc **c, 
                           const Ti *ldc, 
                           std::int64_t group_count, 
                           const Ti *group_size, 
                           compute_mode mode = compute_mode::unset, 
                           const std::vector<sycl::event> &dependencies = {}) 

    sycl::event gemm_batch(sycl::queue &queue, 
                           const sycl::span<oneapi::mkl::transpose> &transa, 
                           const sycl::span<oneapi::mkl::transpose> &transb, 
                           const sycl::span<std::int64_t> &m, 
                           const sycl::span<std::int64_t> &n, 
                           const sycl::span<std::int64_t> &k, 
                           const sycl::span<Ts> &alpha, 
                           const sycl::span<const Ta*> &a, 
                           const sycl::span<std::int64_t> &lda, 
                           const sycl::span<const Tb*> &b, 
                           const sycl::span<std::int64_t> &ldb, 
                           const sycl::span<Ts> &beta, 
                           sycl::span<Tc*> &c, 
                           const sycl::span<std::int64_t> &ldc, 
                           size_t group_count, 
                           const sycl::span<size_t> &group_sizes, 
                           compute_mode mode = compute_mode::unset, 
                           const std::vector<sycl::event> &dependencies = {}) 
}
namespace oneapi::mkl::blas::row_major { 
    sycl::event gemm_batch(sycl::queue &queue, 
                           const oneapi::mkl::transpose *transa, 
                           const oneapi::mkl::transpose *transb, 
                           const Ti *m, 
                           const Ti *n, 
                           const Ti *k, 
                           const Ts *alpha, 
                           const Ta **a, 
                           const Ti *lda, 
                           const Tb **b, 
                           const Ti *ldb, 
                           const Ts *beta, 
                           Tc **c, 
                           const Ti *ldc, 
                           std::int64_t group_count, 
                           const Ti *group_size, 
                           compute_mode mode = compute_mode::unset, 
                           const std::vector<sycl::event> &dependencies = {}) 

    sycl::event gemm_batch(sycl::queue &queue, 
                           const sycl::span<oneapi::mkl::transpose> &transa, 
                           const sycl::span<oneapi::mkl::transpose> &transb, 
                           const sycl::span<std::int64_t> &m, 
                           const sycl::span<std::int64_t> &n, 
                           const sycl::span<std::int64_t> &k, 
                           const sycl::span<Ts> &alpha, 
                           const sycl::span<const Ta*> &a, 
                           const sycl::span<std::int64_t> &lda, 
                           const sycl::span<const Tb*> &b, 
                           const sycl::span<std::int64_t> &ldb, 
                           const sycl::span<Ts> &beta, 
                           sycl::span<Tc*> &c, 
                           const sycl::span<std::int64_t> &ldc, 
                           size_t group_count, 
                           const sycl::span<size_t> &group_sizes, 
                           compute_mode mode = compute_mode::unset, 
                           const std::vector<sycl::event> &dependencies = {}) 
}

入力パラメーター#

queue

ルーチンを実行するキュー。

transa

group_count oneapi::mkl::transpose 値の配列またはスパン。transa[i] は グループ i 内の行列 A に適用される転置演算 op(A) を指定します。詳細はデータタイプを参照してください。

transb

group_count oneapi::mkl::transpose 値の配列またはスパン。transb[i] は グループ i 内の行列 B に適用される転置演算 op(B) を指定します。詳細はデータタイプを参照してください。

m

group_count 整数の配列またはスパン。m[i] はグループ i 内の行列 op(A) と行列 C の行数を指定します。すべてのエントリーは 0 以上である必要があります。

n

group_count 整数の配列またはスパン。n[i] はグループ i 内の行列 op(B) と行列 C の列数を指定します。すべてのエントリーは 0 以上である必要があります。

k

group_count 整数の配列またはスパン。k[i] はグループ i 内の行列 op(A) の列数と行列 op(B) の行数を指定します。すべてのエントリーは 0 以上である必要があります。

alpha

group_count スカラー要素は配列またはスパン。alpha[i] はグループ i 内の行列 - 行列積のスケーリング係数を指定します。

a

入力行列 Atotal_batch_count ポインターの配列、または total_batch_count 入力行列 A のスパン。行列ストレージを参照してください。

transa = transpose::nontrans

transa = transpose::trans または trans = transpose::conjtrans

列優先

配列 A[i] のサイズは lda[i] * k[i] 以上でなければなりません

配列 A[i] のサイズは lda[i] * m[i] 以上でなければなりません

行優先

配列 A[i] のサイズは lda[i] * m[i] 以上でなければなりません

配列 A[i] のサイズは lda[i] * k[i] 以上でなければなりません

lda

group_count 整数の配列またはスパン。lda[i] はグループ i 内の行列 A の先頭次元を指定します。正である必要があります。

transa = transpose::nontrans

transa = transpose::trans または trans = transpose::conjtrans

列優先

m[i] 以上である必要があります

k[i] 以上である必要があります

行優先

k[i] 以上である必要があります

m[i] 以上である必要があります

b

入力行列 Btotal_batch_count ポインターの配列、または total_batch_count 入力行列 B のスパン。行列ストレージを参照してください。

transa = transpose::nontrans

transa = transpose::trans または trans = transpose::conjtrans

列優先

配列 B[i] のサイズは ldb[i] * n[i] 以上でなければなりません

配列 B[i] のサイズは ldb[i] * k[i] 以上でなければなりません

行優先

配列 B[i] のサイズは ldb[i] * k[i] 以上でなければなりません

配列 B[i] のサイズは ldb[i] * n[i] 以上でなければなりません

ldb

group_count 整数の配列またはスパン。ldb[i] はグループ i 内の行列 B の先頭次元を指定します。正である必要があります。

transa = transpose::nontrans

transa = transpose::trans または trans = transpose::conjtrans

列優先

k[i] 以上である必要があります

n[i] 以上である必要があります

行優先

n[i] 以上である必要があります

k[i] 以上である必要があります

beta

group_count スカラー要素は配列またはスパン。beta[i] はグループ i 内の行列 C のスケーリング係数を指定します。

c

入力/出力行列 Ctotal_batch_count ポインターの配列、または total_batch_count 入力/出力行列 C のスパン。行列ストレージを参照してください。

列優先

配列 C[i] のサイズは ldc[i] * n[i] 以上でなければなりません

行優先

配列 C[i] のサイズは ldc[i] * m[i] 以上でなければなりません

ldc

group_count 整数の配列またはスパン。ldc[i] はグループ i 内の行列 C の先頭次元を指定します。正である必要があります。

列優先

m[i] 以上である必要があります

行優先

n[i] 以上である必要があります

group_count

グループの数最小値はゼロです。

group_size

group_count 整数の配列またはスパン。group_size[i] はグループ i 内の gemm 操作の数を指定します。group_size 内の各要素は 0 以上である必要があります。

mode

オプション。計算モードの設定。詳細は計算モードを参照してください。

dependencies

オプション。計算を開始する前に待機するイベントのリスト (存在する場合)。省略した場合、依存関係はデフォルトでなくなります。

modedependencies はそれぞれ省略できます。dependencies を提供するのに mode を指定する必要はありません。

出力パラメーター#

c

alpha * op(A) * op(B) + beta * C 形式の total_batch_count gemm 演算によって上書きされる出力行列 C へのポインター配列。

beta = 0 の場合、gemm_batch を呼び出す前に行列 C を初期化する必要はありません。

戻り値#

計算が完了したことを確認するために待機する出力イベント。

#

USM バージョンの gemm_batch の使用方法の例は、oneMKL インストール・ディレクトリーの次の場所にあります。

share/doc/mkl/examples/sycl/blas/source/gemm_batch.cpp

ストライド API#

ストライド API 操作は次のように定義されます。

for i = 0 … batch_size – 1 
    A、B、C は、a、b、c 内のオフセット i * stridea、i * strideb、i * stridec にある行列です。
    C = alpha * op(A) * op(B) + beta * C 
end for

説明:

  • op(X) は、op(X) = X、op(X) = XT、または op(X) = XH のいずれかです

  • alphabeta はスカラーです

  • AB、および C は行列です

  • op(A) は m x k、op(B) は k x n、そして Cm x n です。

ストライド API の場合、ab および c 配列にすべての入力行列が含まれます。行列間のストライドは、ストライド・パラメーターで指定されます。ab および c 配列内の行列の合計数は、batch_size パラメーターで指定されます。

構文#

namespace oneapi::mkl::blas::column_major { 
    sycl::event gemm_batch(sycl::queue &queue, 
                           oneapi::mkl::transpose transa, 
                           oneapi::mkl::transpose transb, 
                           std::int64_t m, 
                           std::int64_t n, 
                           std::int64_t k, 
                           oneapi::mkl::value_or_pointer<Ts> alpha, 
                           const Ta *a, 
                           std::int64_t lda, 
                           std::int64_t stridea, 
                           const Tb *b, 
                           std::int64_t ldb, 
                           std::int64_t strideb, 
                           oneapi::mkl::value_or_pointer<Ts> beta, 
                           Tc *c, 
                           std::int64_t ldc, 
                           std::int64_t stridec, 
                           std::int64_t batch_size, 
                           compute_mode mode = compute_mode::unset, 
                           const std::vector<sycl::event> &dependencies = {}) 
}
namespace oneapi::mkl::blas::row_major { 
    sycl::event gemm_batch(sycl::queue &queue, 
                           oneapi::mkl::transpose transa, 
                           oneapi::mkl::transpose transb, 
                           std::int64_t m, 
                           std::int64_t n,
                           std::int64_t k, 
                           oneapi::mkl::value_or_pointer<Ts> alpha, 
                           const Ta *a, 
                           std::int64_t lda, 
                           std::int64_t stridea, 
                           const Tb *b, 
                           std::int64_t ldb, 
                           std::int64_t strideb, 
                           oneapi::mkl::value_or_pointer<Ts> beta, 
                           Tc *c, 
                           std::int64_t ldc, 
                           std::int64_t stridec, 
                           std::int64_t batch_size, 
                           compute_mode mode = compute_mode::unset, 
                           const std::vector<sycl::event> &dependencies = {}) 
}

入力パラメーター#

queue

ルーチンを実行するキュー。

transa

行列 A に適用される転置演算 op(A) を指定します。詳細はデータタイプを参照してください。

transb

行列 B に適用される転置演算 op(B) を指定します。詳細はデータタイプを参照してください。

m

行列 op(A) と行列 C の行数。最小値はゼロです。

n

行列 op(B) と行列 C の列数。最小値はゼロです。

k

行列 op(A) の列数と行列 op(B) の行数。最小値はゼロです。

alpha

行列 - 行列積のスケーリング係数。value_or_pointer データタイプの詳細については、スカラー引数を参照してください。

a

入力行列 A へのポインター。配列のサイズは stridea * batch_size 以上である必要があります。

lda

行列 A の先頭次元。正である必要があります。

transa = transpose::nontrans

transa = transpose::trans または trans = transpose::conjtrans

列優先

m 以上である必要があります

k 以上である必要があります

行優先

k 以上である必要があります

m 以上である必要があります

stridea

連続する 2 つの行列 A 間のストライド。

transa = transpose::nontrans

transa = transpose::trans または trans = transpose::conjtrans

列優先

lda * k 以上である必要があります

lda * m 以上である必要があります

行優先

lda * m 以上である必要があります

lda * k 以上である必要があります

b

入力行列 B へのポインター。配列のサイズは strideb * batch_size 以上である必要があります。

ldb

行列 B の先頭次元。正である必要があります。

transa = transpose::nontrans

transa = transpose::trans または trans = transpose::conjtrans

列優先

k 以上である必要があります

n 以上である必要があります

行優先

n 以上である必要があります

k 以上である必要があります

strideb

連続する 2 つの行列 B 間のストライド。

transa = transpose::nontrans

transa = transpose::trans または trans = transpose::conjtrans

列優先

ldb * n 以上である必要があります

ldb * k 以上である必要があります

行優先

ldb * k 以上である必要があります

ldb * n 以上である必要があります

beta

行列 C のスケーリング係数。value_or_pointer データタイプの詳細については、スカラー引数を参照してください。

c

入力/出力行列 C へのポインター。配列のサイズは stridec * batch_size 以上である必要があります。

ldc

行列 C の先頭次元。正である必要があります。

列優先

m 以上である必要があります

行優先

n 以上である必要があります

stridec

連続する 2 つの行列 C 間のストライド。

列優先

ldc * n 以上である必要があります

行優先

ldc * m 以上である必要があります

batch_size

実行する行列乗算演算の数を指定します。最小値はゼロです。

mode

オプション。計算モードの設定。詳細は計算モードを参照してください。

dependencies

オプション。計算を開始する前に待機するイベントのリスト (存在する場合)。省略した場合、依存関係はデフォルトでなくなります。

modedependencies はそれぞれ省略できます。dependencies を提供するのに mode を指定する必要はありません。

出力パラメーター#

c

alpha * op(A) * op(B) + beta * C 形式の batch_size gemm 演算によって上書きされる出力行列 C へのポインター。

beta = 0 の場合、gemm_batch を呼び出す前に行列 C を初期化する必要はありません。

戻り値#

計算が完了したことを確認するために待機する出力イベント。