GEMM¶
General matrix multiplication functions for locally computing \( C \leftarrow \alpha OP(A) OP(B) + \beta C \).
Functions
-
void gemm(SplaOperation opA, SplaOperation opB, int m, int n, int k, float alpha, const float *A, int lda, const float *B, int ldb, float beta, float *C, int ldc, Context &ctx)¶
Computes a local general matrix multiplication of the form \( C \leftarrow \alpha op(A) op(B) + \beta C \) in single precision.
If context with processing unit set to GPU, pointers to matrices can be any combination of host and device pointers.
- Parameters
opA – [in] Operation applied when reading matrix \(A\).
opB – [in] Operation applied when reading matrix \(B\).
m – [in] Number of rows of \(OP(A)\)
n – [in] Number of columns of \(OP(B)\)
k – [in] Number rows of \(OP(B)\) and number of columns of \(OP(A)\)
alpha – [in] Scaling of multiplication of \(A^H\) and \(B\)
A – [in] Pointer to matrix \(A\).
lda – [in] Leading dimension of \(A\).
B – [in] Pointer to matrix \(B\).
ldb – [in] Leading dimension of \(B\).
beta – [in] Scaling of \(C\) before summation.
C – [out] Pointer to matrix \(C\).
ldc – [in] Leading dimension of \(C\).
ctx – [in] Context, which provides configuration settings and reusable resources.
-
void gemm(SplaOperation opA, SplaOperation opB, int m, int n, int k, double alpha, const double *A, int lda, const double *B, int ldb, double beta, double *C, int ldc, Context &ctx)¶
Computes a local general matrix multiplication of the form \( C \leftarrow \alpha OP(A) OP(B) + \beta C \) in double precision.
See documentation above.
-
void gemm(SplaOperation opA, SplaOperation opB, int m, int n, int k, std::complex<float> alpha, const std::complex<float> *A, int lda, const std::complex<float> *B, int ldb, std::complex<float> beta, std::complex<float> *C, int ldc, Context &ctx)¶
Computes a local general matrix multiplication of the form \( C \leftarrow \alpha OP(A) OP(B) + \beta C \) in single precision complex types.
See documentation above.
-
void gemm(SplaOperation opA, SplaOperation opB, int m, int n, int k, std::complex<double> alpha, const std::complex<double> *A, int lda, const std::complex<double> *B, int ldb, std::complex<double> beta, std::complex<double> *C, int ldc, Context &ctx)¶
Computes a local general matrix multiplication of the form \( C \leftarrow \alpha OP(A) OP(B) + \beta C \) in double precision complex types.
See documentation above.