GEMM - SBS

General matrix multiplication functions for computing \( C \leftarrow \alpha A B + \beta C \) with stripe-block-stipe distribution.

  ------                     ------         ------
  |    |                     |    |         |    |
  |    |                     |    |         |    |
  ------                     ------         ------
  |    |       -------       |    |         |    |
  |    |       |  |  |       |    |         |    |
  ------   *   -------   +   ------   -->   ------
  |    |       |  |  |       |    |         |    |
  |    |       -------       |    |         |    |
  ------          B          ------         ------
  |    |                     |    |         |    |
  |    |                     |    |         |    |
  ------                     ------         ------
    A                          C              C

Functions

void pgemm_sbs(int mLocal, int n, int k, float alpha, const float *A, int lda, const float *B, int ldb, int bRowOffset, int bColOffset, MatrixDistribution &distB, float beta, float *C, int ldc, Context &ctx)

Computes a distributed general matrix multiplication of the form \( C \leftarrow \alpha A B + \beta C \) in single precision.

\(A\) and \(C\) are only split along the row dimension (stripes), while \(B\) can be distributed as any supported MatrixDistribution type.

Parameters
  • mLocal[in] Number rows of \(A\) and \(C\) stored at calling MPI rank. This number may differ for each rank.

  • n[in] Number of columns of \(B\).

  • k[in] Number of columns of \(C\) and rows of \(B\).

  • alpha[in] Scaling of multiplication of \(A^H\) and \(B\)

  • A[in] Pointer to matrix \(A\).

  • lda[in] Leading dimension of \(A\) with lda \(\geq\) kLocal.

  • B[in] Pointer to matrix \(B\).

  • ldb[in] Leading dimension of \(B\) with ldb \(\geq\) loc(k), where loc(k) is the number of locally stored rows of \(B\).

  • bRowOffset[in] Row offset in the global matrix \(B\), identifying the first row of the submatrix \(B\).

  • bColOffset[in] Column offset in the global matrix \(B\), identifying the first coloumn of the submatrix \(B\).

  • distB[in] Matrix distribution of global matrix \(B\).

  • beta[in] Scaling of \(C\) before summation.

  • C[out] Pointer to matrix \(C\).

  • ldc[in] Leading dimension of \(C\) with ldC \(\geq\) mLocal.

  • ctx[in] Context, which provides configuration settings and reusable resources.

void pgemm_sbs(int mLocal, int n, int k, double alpha, const double *A, int lda, const double *B, int ldb, int bRowOffset, int bColOffset, MatrixDistribution &distB, double beta, double *C, int ldc, Context &ctx)

Computes a distributed general matrix multiplication of the form \( C \leftarrow \alpha A B + \beta C \) in double precision.

See documentation above.

void pgemm_sbs(int mLocal, int n, int k, std::complex<float> alpha, const std::complex<float> *A, int lda, const std::complex<float> *B, int ldb, int bRowOffset, int bColOffset, MatrixDistribution &distB, std::complex<float> beta, std::complex<float> *C, int ldc, Context &ctx)

Computes a distributed general matrix multiplication of the form \( C \leftarrow \alpha A B + \beta C \) in double precision.

See documentation above.

void pgemm_sbs(int mLocal, int n, int k, std::complex<double> alpha, const std::complex<double> *A, int lda, const std::complex<double> *B, int ldb, int bRowOffset, int bColOffset, MatrixDistribution &distB, std::complex<double> beta, std::complex<double> *C, int ldc, Context &ctx)

Computes a distributed general matrix multiplication of the form \( C \leftarrow \alpha A B + \beta C \) in double precision.

See documentation above.