You can not select more than 25 topics Topics must start with a chinese character,a letter or number, can include dashes ('-') and can be up to 35 characters long.

gemm_vec.c 24 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686
  1. /*
  2. * Copyright (c) IBM Corporation 2020.
  3. * All rights reserved.
  4. *
  5. * Redistribution and use in source and binary forms, with or without
  6. * modification, are permitted provided that the following conditions are
  7. * met:
  8. *
  9. * 1. Redistributions of source code must retain the above copyright
  10. * notice, this list of conditions and the following disclaimer.
  11. *
  12. * 2. Redistributions in binary form must reproduce the above copyright
  13. * notice, this list of conditions and the following disclaimer in
  14. * the documentation and/or other materials provided with the
  15. * distribution.
  16. * 3. Neither the name of the OpenBLAS project nor the names of
  17. * its contributors may be used to endorse or promote products
  18. * derived from this software without specific prior written
  19. * permission.
  20. *
  21. * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
  22. * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
  23. * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
  24. * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
  25. * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL
  26. * DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR
  27. * SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
  28. * CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY,
  29. * OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE
  30. * USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
  31. */
  32. #include "common.h"
  33. #include <vecintrin.h>
  34. #include <stdbool.h>
  35. #include <stdio.h>
  36. #include <stdlib.h>
  37. #ifdef COMPLEX
  38. #error "Handling for complex numbers is not supported in this kernel"
  39. #endif
  40. #ifdef DOUBLE
  41. #define UNROLL_M DGEMM_DEFAULT_UNROLL_M
  42. #define UNROLL_N DGEMM_DEFAULT_UNROLL_N
  43. #else
  44. #define UNROLL_M SGEMM_DEFAULT_UNROLL_M
  45. #define UNROLL_N SGEMM_DEFAULT_UNROLL_N
  46. #endif
  47. static const size_t unroll_m = UNROLL_M;
  48. static const size_t unroll_n = UNROLL_N;
  49. /* Handling of triangular matrices */
  50. #ifdef TRMMKERNEL
  51. static const bool trmm = true;
  52. static const bool left =
  53. #ifdef LEFT
  54. true;
  55. #else
  56. false;
  57. #endif
  58. static const bool backwards =
  59. #if defined(LEFT) != defined(TRANSA)
  60. true;
  61. #else
  62. false;
  63. #endif
  64. #else
  65. static const bool trmm = false;
  66. static const bool left = false;
  67. static const bool backwards = false;
  68. #endif /* TRMMKERNEL */
  69. /*
  70. * Background:
  71. *
  72. * The algorithm of GotoBLAS / OpenBLAS breaks down the matrix multiplication
  73. * problem by splitting all matrices into partitions multiple times, so that the
  74. * submatrices fit into the L1 or L2 caches. As a result, each multiplication of
  75. * submatrices can stream data fast from L1 and L2 caches. Inbetween, it copies
  76. * and rearranges the submatrices to enable contiguous memory accesses to
  77. * improve locality in both caches and TLBs.
  78. *
  79. * At the heart of the algorithm is this kernel, which multiplies, a "Block
  80. * matrix" A (small dimensions) with a "Panel matrix" B (number of rows is
  81. * small) and adds the result into a "Panel matrix" C; GotoBLAS calls this
  82. * operation GEBP. This kernel further partitions GEBP twice, such that (1)
  83. * submatrices of C and B fit into the L1 caches (GEBP_column_block) and (2) a
  84. * block of C fits into the registers, while multiplying panels from A and B
  85. * streamed from the L2 and L1 cache, respectively (GEBP_block).
  86. *
  87. *
  88. * Algorithm GEBP(A, B, C, m, n, k, alpha):
  89. *
  90. * The problem is calculating C += alpha * (A * B)
  91. * C is an m x n matrix, A is an m x k matrix, B is an k x n matrix.
  92. *
  93. * - C is in column-major-order, with an offset of ldc to the element in the
  94. * next column (same row).
  95. * - A is in row-major-order yet stores SGEMM_UNROLL_M elements of each column
  96. * contiguously while walking along rows.
  97. * - B is in column-major-order but packs SGEMM_UNROLL_N elements of a row
  98. * contiguously.
  99. * If the numbers of rows and columns are not multiples of SGEMM_UNROLL_M or
  100. * SGEMM_UNROLL_N, the remaining elements are arranged in blocks with power-of-2
  101. * dimensions (e.g., 5 remaining columns would be in a block-of-4 and a
  102. * block-of-1).
  103. *
  104. * Note that packing A and B into that form is taken care of by the caller in
  105. * driver/level3/level3.c (actually done by "copy kernels").
  106. *
  107. * Steps:
  108. * - Partition C and B into blocks of n_r (SGEMM_UNROLL_N) columns, C_j and B_j.
  109. * Now, B_j should fit into the L1 cache.
  110. * - For each partition, calculate C_j += alpha * (A * B_j) by
  111. * (1) Calculate C_aux := A * B_j (see below)
  112. * (2) unpack C_j = C_j + alpha * C_aux
  113. *
  114. *
  115. * Algorithm for Calculating C_aux:
  116. *
  117. * - Further partition C_aux and A into groups of m_r (SGEMM_UNROLL_M) rows,
  118. * such that the m_r x n_r-submatrix of C_aux can be held in registers. Each
  119. * submatrix of C_aux can be calculated independently, and the registers are
  120. * added back into C_j.
  121. *
  122. * - For each row-block of C_aux:
  123. * (uses a row block of A and full B_j)
  124. * - stream over all columns of A, multiply with elements from B and
  125. * accumulate in registers. (use different inner-kernels to exploit
  126. * vectorization for varying block sizes)
  127. * - add alpha * row block of C_aux back into C_j.
  128. *
  129. * Note that there are additional mechanics for handling triangular matrices,
  130. * calculating B := alpha (A * B) where either of the matrices A or B can be
  131. * triangular. In case of A, the macro "LEFT" is defined. In addition, A can
  132. * optionally be transposed.
  133. * The code effectively skips an "offset" number of columns in A and rows of B
  134. * in each block, to save unnecessary work by exploiting the triangular nature.
  135. * To handle all cases, the code discerns (1) a "left" mode when A is triangular
  136. * and (2) "forward" / "backwards" modes where only the first "offset"
  137. * columns/rows of A/B are used or where the first "offset" columns/rows are
  138. * skipped, respectively.
  139. *
  140. * Reference:
  141. *
  142. * The summary above is based on staring at various kernel implementations and:
  143. * K. Goto and R. A. Van de Geijn, Anatomy of High-Performance Matrix
  144. * Multiplication, in ACM Transactions of Mathematical Software, Vol. 34, No.
  145. * 3, May 2008.
  146. */
  147. #define VLEN_BYTES 16
  148. #define VLEN_FLOATS (VLEN_BYTES / sizeof(FLOAT))
  149. typedef FLOAT vector_float __attribute__ ((vector_size (16)));
  150. /**
  151. * Load a vector into register, and hint on 8-byte alignment to improve
  152. * performance. gcc-9 and newer will create these hints by itself. For older
  153. * compiler versions, use inline assembly to explicitly express the hint.
  154. * Provide explicit hex encoding to cater for binutils versions that do not know
  155. * about vector-load with alignment hints yet.
  156. *
  157. * Note that, for block sizes where we apply vectorization, vectors in A will
  158. * always be 8-byte aligned.
  159. */
  160. static inline vector_float vec_load_hinted(FLOAT const *restrict a) {
  161. vector_float const *restrict addr = (vector_float const *restrict)a;
  162. vector_float y;
  163. #if __GNUC__ < 9
  164. // hex-encode vl %[out],%[addr],3
  165. asm(".insn vrx,0xe70000003006,%[out],%[addr],3"
  166. : [ out ] "=v"(y)
  167. : [ addr ] "R"(*addr));
  168. #else
  169. y = *addr;
  170. #endif
  171. return y;
  172. }
  173. /**
  174. * Calculate for a row-block in C_i of size ROWSxCOLS using vector intrinsics.
  175. *
  176. * @param[in] A Pointer current block of input matrix A.
  177. * @param[in] k Number of columns in A.
  178. * @param[in] B Pointer current block of input matrix B.
  179. * @param[inout] C Pointer current block of output matrix C.
  180. * @param[in] ldc Offset between elements in adjacent columns in C.
  181. * @param[in] alpha Scalar factor.
  182. */
  183. #define VECTOR_BLOCK(ROWS, COLS) \
  184. static inline void GEBP_block_##ROWS##_##COLS( \
  185. FLOAT const *restrict A, BLASLONG bk, FLOAT const *restrict B, \
  186. FLOAT *restrict C, BLASLONG ldc, FLOAT alpha) { \
  187. _Static_assert( \
  188. ROWS % VLEN_FLOATS == 0, \
  189. "rows in block must be multiples of vector length"); \
  190. vector_float Caux[ROWS / VLEN_FLOATS][COLS]; \
  191. \
  192. for (BLASLONG i = 0; i < ROWS / VLEN_FLOATS; i++) { \
  193. vector_float A0 = \
  194. vec_load_hinted(A + i * VLEN_FLOATS); \
  195. for (BLASLONG j = 0; j < COLS; j++) \
  196. Caux[i][j] = A0 * B[j]; \
  197. } \
  198. \
  199. /* \
  200. * Stream over the row-block of A, which is packed \
  201. * column-by-column, multiply by coefficients in B and add up \
  202. * into temporaries Caux (which the compiler will hold in \
  203. * registers). Vectorization: Multiply column vectors from A \
  204. * with scalars from B and add up in column vectors of Caux. \
  205. * That equates to unrolling the loop over rows (in i) and \
  206. * executing each unrolled iteration as a vector element. \
  207. */ \
  208. for (BLASLONG k = 1; k < bk; k++) { \
  209. for (BLASLONG i = 0; i < ROWS / VLEN_FLOATS; i++) { \
  210. vector_float Ak = vec_load_hinted( \
  211. A + i * VLEN_FLOATS + k * ROWS); \
  212. \
  213. for (BLASLONG j = 0; j < COLS; j++) \
  214. Caux[i][j] += Ak * B[j + k * COLS]; \
  215. } \
  216. } \
  217. \
  218. /* \
  219. * Unpack row-block of C_aux into outer C_i, multiply by \
  220. * alpha and add up. \
  221. */ \
  222. for (BLASLONG j = 0; j < COLS; j++) { \
  223. for (BLASLONG i = 0; i < ROWS / VLEN_FLOATS; i++) { \
  224. vector_float *C_ij = \
  225. (vector_float *)(C + i * VLEN_FLOATS + \
  226. j * ldc); \
  227. if (trmm) { \
  228. *C_ij = alpha * Caux[i][j]; \
  229. } else { \
  230. *C_ij += alpha * Caux[i][j]; \
  231. } \
  232. } \
  233. } \
  234. }
  235. #if UNROLL_M == 16
  236. VECTOR_BLOCK(16, 2)
  237. VECTOR_BLOCK(16, 1)
  238. #endif
  239. #if UNROLL_N == 8
  240. VECTOR_BLOCK(8, 8)
  241. VECTOR_BLOCK(4, 8)
  242. #endif
  243. #ifndef DOUBLE
  244. VECTOR_BLOCK(8, 4)
  245. #endif
  246. VECTOR_BLOCK(8, 2)
  247. VECTOR_BLOCK(8, 1)
  248. VECTOR_BLOCK(4, 4)
  249. VECTOR_BLOCK(4, 2)
  250. VECTOR_BLOCK(4, 1)
  251. #ifdef DOUBLE
  252. VECTOR_BLOCK(2, 4)
  253. VECTOR_BLOCK(2, 2)
  254. VECTOR_BLOCK(2, 1)
  255. #endif
  256. /**
  257. * Calculate a row-block that fits 4x4 vector registers using a loop
  258. * unrolled-by-2 with explicit interleaving to better overlap loads and
  259. * computation.
  260. * This function fits 16x4 blocks for SGEMM and 8x4 blocks for DGEMM.
  261. */
  262. #ifdef DOUBLE
  263. static inline void GEBP_block_8_4(
  264. #else // float
  265. static inline void GEBP_block_16_4(
  266. #endif
  267. FLOAT const *restrict A, BLASLONG bk, FLOAT const *restrict B,
  268. FLOAT *restrict C, BLASLONG ldc, FLOAT alpha) {
  269. #define VEC_ROWS 4
  270. #define VEC_COLS 4
  271. #define ROWS VEC_ROWS * VLEN_FLOATS
  272. #define COLS (VEC_COLS)
  273. /*
  274. * Hold intermediate results in vector registers.
  275. * Since we need to force the compiler's hand in places, we need to use
  276. * individual variables in contrast to the generic implementation's
  277. * arrays.
  278. */
  279. #define INIT_ROW_OF_C(ROW) \
  280. vector_float A##ROW = vec_load_hinted(A + ROW * VLEN_FLOATS); \
  281. vector_float C_##ROW##_0 = A##ROW * B[0]; \
  282. vector_float C_##ROW##_1 = A##ROW * B[1]; \
  283. vector_float C_##ROW##_2 = A##ROW * B[2]; \
  284. vector_float C_##ROW##_3 = A##ROW * B[3];
  285. INIT_ROW_OF_C(0)
  286. INIT_ROW_OF_C(1)
  287. INIT_ROW_OF_C(2)
  288. INIT_ROW_OF_C(3)
  289. #undef INIT_ROW_OF_C
  290. if (bk > 1) {
  291. BLASLONG k = 1;
  292. vector_float Ak[VEC_ROWS], Aknext[VEC_ROWS];
  293. vector_float Bk[VEC_COLS], Bknext[VEC_COLS];
  294. /*
  295. * Note that in several places, we enforce an instruction
  296. * sequence that we identified empirically by utilizing dummy
  297. * asm statements.
  298. */
  299. for (BLASLONG j = 0; j < VEC_COLS; j++)
  300. Bk[j] = vec_splats(B[j + k * COLS]);
  301. asm("");
  302. for (BLASLONG i = 0; i < VEC_ROWS; i++)
  303. Ak[i] = vec_load_hinted(A + i * VLEN_FLOATS + k * ROWS);
  304. for (; k < (bk - 2); k += 2) {
  305. /*
  306. * Load inputs for (k+1) into registers.
  307. * Loading from B first is advantageous.
  308. */
  309. for (BLASLONG j = 0; j < VEC_COLS; j++)
  310. Bknext[j] = vec_splats(B[j + (k + 1) * COLS]);
  311. asm("");
  312. for (BLASLONG i = 0; i < VEC_ROWS; i++)
  313. Aknext[i] = vec_load_hinted(A + i * VLEN_FLOATS +
  314. (k + 1) * ROWS);
  315. /*
  316. * To achieve better instruction-level parallelism,
  317. * make sure to first load input data for (k+1) before
  318. * initiating compute for k. We enforce that ordering
  319. * with a pseudo asm statement.
  320. * Note that we need to massage this particular "barrier"
  321. * depending on the gcc version.
  322. */
  323. #if __GNUC__ > 7
  324. #define BARRIER_READ_BEFORE_COMPUTE(SUFFIX) \
  325. do { \
  326. asm("" \
  327. : "+v"(C_0_0), "+v"(C_0_1), "+v"(C_0_2), "+v"(C_0_3), "+v"(C_1_0), \
  328. "+v"(C_1_1), "+v"(C_1_2), "+v"(C_1_3) \
  329. : "v"(B##SUFFIX[0]), "v"(B##SUFFIX[1]), "v"(B##SUFFIX[2]), \
  330. "v"(B##SUFFIX[3]), "v"(A##SUFFIX[0]), "v"(A##SUFFIX[1]), \
  331. "v"(A##SUFFIX[2]), "v"(A##SUFFIX[3])); \
  332. asm("" \
  333. : "+v"(C_2_0), "+v"(C_2_1), "+v"(C_2_2), "+v"(C_2_3), "+v"(C_3_0), \
  334. "+v"(C_3_1), "+v"(C_3_2), "+v"(C_3_3) \
  335. : "v"(B##SUFFIX[0]), "v"(B##SUFFIX[1]), "v"(B##SUFFIX[2]), \
  336. "v"(B##SUFFIX[3]), "v"(A##SUFFIX[0]), "v"(A##SUFFIX[1]), \
  337. "v"(A##SUFFIX[2]), "v"(A##SUFFIX[3])); \
  338. } while (0)
  339. #else // __GNUC__ <= 7
  340. #define BARRIER_READ_BEFORE_COMPUTE(SUFFIX) \
  341. do { \
  342. asm(""); \
  343. } while (0)
  344. #endif
  345. BARRIER_READ_BEFORE_COMPUTE(knext);
  346. /* Compute for (k) */
  347. C_0_0 += Ak[0] * Bk[0];
  348. C_1_0 += Ak[1] * Bk[0];
  349. C_2_0 += Ak[2] * Bk[0];
  350. C_3_0 += Ak[3] * Bk[0];
  351. C_0_1 += Ak[0] * Bk[1];
  352. C_1_1 += Ak[1] * Bk[1];
  353. C_2_1 += Ak[2] * Bk[1];
  354. C_3_1 += Ak[3] * Bk[1];
  355. C_0_2 += Ak[0] * Bk[2];
  356. C_1_2 += Ak[1] * Bk[2];
  357. C_2_2 += Ak[2] * Bk[2];
  358. C_3_2 += Ak[3] * Bk[2];
  359. C_0_3 += Ak[0] * Bk[3];
  360. C_1_3 += Ak[1] * Bk[3];
  361. C_2_3 += Ak[2] * Bk[3];
  362. C_3_3 += Ak[3] * Bk[3];
  363. asm("");
  364. /*
  365. * Load inputs for (k+2) into registers.
  366. * First load from B.
  367. */
  368. for (BLASLONG j = 0; j < VEC_COLS; j++)
  369. Bk[j] = vec_splats(B[j + (k + 2) * COLS]);
  370. asm("");
  371. for (BLASLONG i = 0; i < VEC_ROWS; i++)
  372. Ak[i] = vec_load_hinted(A + i * VLEN_FLOATS + (k + 2) * ROWS);
  373. /*
  374. * As above, make sure to first schedule the loads for (k+2)
  375. * before compute for (k+1).
  376. */
  377. BARRIER_READ_BEFORE_COMPUTE(k);
  378. /* Compute on (k+1) */
  379. C_0_0 += Aknext[0] * Bknext[0];
  380. C_1_0 += Aknext[1] * Bknext[0];
  381. C_2_0 += Aknext[2] * Bknext[0];
  382. C_3_0 += Aknext[3] * Bknext[0];
  383. C_0_1 += Aknext[0] * Bknext[1];
  384. C_1_1 += Aknext[1] * Bknext[1];
  385. C_2_1 += Aknext[2] * Bknext[1];
  386. C_3_1 += Aknext[3] * Bknext[1];
  387. C_0_2 += Aknext[0] * Bknext[2];
  388. C_1_2 += Aknext[1] * Bknext[2];
  389. C_2_2 += Aknext[2] * Bknext[2];
  390. C_3_2 += Aknext[3] * Bknext[2];
  391. C_0_3 += Aknext[0] * Bknext[3];
  392. C_1_3 += Aknext[1] * Bknext[3];
  393. C_2_3 += Aknext[2] * Bknext[3];
  394. C_3_3 += Aknext[3] * Bknext[3];
  395. }
  396. /* Wrapup remaining k's */
  397. for (; k < bk; k++) {
  398. vector_float Ak;
  399. #define COMPUTE_WRAPUP_ROW(ROW) \
  400. Ak = vec_load_hinted(A + ROW * VLEN_FLOATS + k * ROWS); \
  401. C_##ROW##_0 += Ak * B[0 + k * COLS]; \
  402. C_##ROW##_1 += Ak * B[1 + k * COLS]; \
  403. C_##ROW##_2 += Ak * B[2 + k * COLS]; \
  404. C_##ROW##_3 += Ak * B[3 + k * COLS];
  405. COMPUTE_WRAPUP_ROW(0)
  406. COMPUTE_WRAPUP_ROW(1)
  407. COMPUTE_WRAPUP_ROW(2)
  408. COMPUTE_WRAPUP_ROW(3)
  409. #undef COMPUTE_WRAPUP_ROW
  410. }
  411. }
  412. /*
  413. * Unpack row-block of C_aux into outer C_i, multiply by
  414. * alpha and add up (or assign for TRMM).
  415. */
  416. #define WRITE_BACK_C(ROW, COL) \
  417. do { \
  418. vector_float *Cij = \
  419. (vector_float *)(C + ROW * VLEN_FLOATS + COL * ldc); \
  420. if (trmm) { \
  421. *Cij = alpha * C_##ROW##_##COL; \
  422. } else { \
  423. *Cij += alpha * C_##ROW##_##COL; \
  424. } \
  425. } while (0)
  426. WRITE_BACK_C(0, 0); WRITE_BACK_C(0, 1); WRITE_BACK_C(0, 2); WRITE_BACK_C(0, 3);
  427. WRITE_BACK_C(1, 0); WRITE_BACK_C(1, 1); WRITE_BACK_C(1, 2); WRITE_BACK_C(1, 3);
  428. WRITE_BACK_C(2, 0); WRITE_BACK_C(2, 1); WRITE_BACK_C(2, 2); WRITE_BACK_C(2, 3);
  429. WRITE_BACK_C(3, 0); WRITE_BACK_C(3, 1); WRITE_BACK_C(3, 2); WRITE_BACK_C(3, 3);
  430. #undef WRITE_BACK_C
  431. #undef ROWS
  432. #undef VEC_ROWS
  433. #undef COLS
  434. #undef VEC_COLS
  435. #undef BARRIER_READ_BEFORE_COMPUTE
  436. }
  437. /**
  438. * Handle calculation for row blocks in C_i of any size by dispatching into
  439. * macro-defined (inline) functions or by deferring to a simple generic
  440. * implementation. Note that the compiler can remove this awkward-looking
  441. * dispatching code while inlineing.
  442. *
  443. * @param[in] m Number of rows in block C_i.
  444. * @param[in] n Number of columns in block C_i.
  445. * @param[in] first_row Index of first row of the block C_i (relative to C).
  446. * @param[in] A Pointer to input matrix A (note: all of it).
  447. * @param[in] k Number of columns in A and rows in B.
  448. * @param[in] B Pointer to current column block (panel) of input matrix B.
  449. * @param[inout] C Pointer to current column block (panel) of output matrix C.
  450. * @param[in] ldc Offset between elements in adjacent columns in C.
  451. * @param[in] alpha Scalar factor.
  452. * @param[in] offset Number of columns of A and rows of B to skip (for triangular matrices).
  453. * @param[in] off Running offset for handling triangular matrices.
  454. */
  455. static inline void GEBP_block(BLASLONG m, BLASLONG n,
  456. BLASLONG first_row,
  457. const FLOAT * restrict A, BLASLONG k,
  458. const FLOAT * restrict B,
  459. FLOAT *restrict C, BLASLONG ldc,
  460. FLOAT alpha,
  461. BLASLONG offset, BLASLONG off)
  462. {
  463. if (trmm && left)
  464. off = offset + first_row;
  465. A += first_row * k;
  466. C += first_row;
  467. if (trmm) {
  468. if (backwards) {
  469. A += off * m;
  470. B += off * n;
  471. k -= off;
  472. } else {
  473. if (left) {
  474. k = off + m;
  475. } else {
  476. k = off + n;
  477. }
  478. }
  479. }
  480. #define BLOCK(bm, bn) \
  481. if (m == bm && n == bn) { \
  482. GEBP_block_##bm##_##bn(A, k, B, C, ldc, alpha); \
  483. return; \
  484. }
  485. #if UNROLL_M == 16
  486. BLOCK(16, 4); BLOCK(16, 2); BLOCK(16, 1);
  487. #endif
  488. #if UNROLL_N == 8
  489. BLOCK(8, 8); BLOCK(4, 8);
  490. #endif
  491. BLOCK(8, 4); BLOCK(8, 2); BLOCK(8, 1);
  492. BLOCK(4, 4); BLOCK(4, 2); BLOCK(4, 1);
  493. #ifdef DOUBLE
  494. BLOCK(2, 4);
  495. BLOCK(2, 2);
  496. #endif
  497. #undef BLOCK
  498. /* simple implementation for smaller block sizes: */
  499. FLOAT Caux[m][n] __attribute__ ((aligned (16)));
  500. /*
  501. * Peel off first iteration (i.e., column of A) for initializing Caux
  502. */
  503. for (BLASLONG i = 0; i < m; i++)
  504. for (BLASLONG j = 0; j < n; j++)
  505. Caux[i][j] = A[i] * B[j];
  506. for (BLASLONG kk = 1; kk < k; kk++)
  507. for (BLASLONG i = 0; i < m; i++)
  508. for (BLASLONG j = 0; j < n; j++)
  509. Caux[i][j] += A[i + kk * m] * B[j + kk * n];
  510. for (BLASLONG i = 0; i < m; i++)
  511. for (BLASLONG j = 0; j < n; j++)
  512. if (trmm) {
  513. C[i + j * ldc] = alpha * Caux[i][j];
  514. } else {
  515. C[i + j * ldc] += alpha * Caux[i][j];
  516. }
  517. }
  518. /**
  519. * Handle a column block (panel) of C and B while calculating C += alpha(A * B).
  520. *
  521. * @param[in] num_cols Number of columns in the block (in C and B).
  522. * @param[in] first_col First column of the current block (in C and B).
  523. * @param[in] A Pointer to input matrix A.
  524. * @param[in] bk Number of columns in A and rows in B.
  525. * @param[in] B Pointer to input matrix B (note: all of it).
  526. * @param[in] bm Number of rows in C and A.
  527. * @param[inout] C Pointer to output matrix C (note: all of it).
  528. * @param[in] ldc Offset between elements in adjacent columns in C.
  529. * @param[in] alpha Scalar factor.
  530. * @param[in] offset Number of columns of A and rows of B to skip (for triangular matrices).
  531. */
  532. static inline void GEBP_column_block(BLASLONG num_cols, BLASLONG first_col,
  533. const FLOAT *restrict A, BLASLONG bk,
  534. const FLOAT *restrict B, BLASLONG bm,
  535. FLOAT *restrict C, BLASLONG ldc,
  536. FLOAT alpha,
  537. BLASLONG const offset) {
  538. FLOAT *restrict C_i = C + first_col * ldc;
  539. /*
  540. * B is in column-order with n_r packed row elements, which does
  541. * not matter -- we always move in full such blocks of
  542. * column*pack
  543. */
  544. const FLOAT *restrict B_i = B + first_col * bk;
  545. BLASLONG off = 0;
  546. if (trmm) {
  547. if (left) {
  548. off = offset;
  549. } else {
  550. off = -offset + first_col;
  551. }
  552. }
  553. /*
  554. * Calculate C_aux := A * B_j
  555. * then unpack C_i += alpha * C_aux.
  556. *
  557. * For that purpose, further partition C_aux and A into blocks
  558. * of m_r (unroll_m) rows, or powers-of-2 if smaller.
  559. */
  560. BLASLONG row = 0;
  561. for (BLASLONG block_size = unroll_m; block_size > 0; block_size /= 2)
  562. for (; bm - row >= block_size; row += block_size)
  563. GEBP_block(block_size, num_cols, row, A, bk, B_i, C_i,
  564. ldc, alpha, offset, off);
  565. }
  566. /**
  567. * Inner kernel for matrix-matrix multiplication. C += alpha (A * B)
  568. * where C is an m-by-n matrix, A is m-by-k and B is k-by-n. Note that A, B, and
  569. * C are pointers to submatrices of the actual matrices.
  570. *
  571. * For triangular matrix multiplication, calculate B := alpha (A * B) where A
  572. * or B can be triangular (in case of A, the macro LEFT will be defined).
  573. *
  574. * @param[in] bm Number of rows in C and A.
  575. * @param[in] bn Number of columns in C and B.
  576. * @param[in] bk Number of columns in A and rows in B.
  577. * @param[in] alpha Scalar factor.
  578. * @param[in] ba Pointer to input matrix A.
  579. * @param[in] bb Pointer to input matrix B.
  580. * @param[inout] C Pointer to output matrix C.
  581. * @param[in] ldc Offset between elements in adjacent columns in C.
  582. * @param[in] offset Number of columns of A and rows of B to skip (for triangular matrices).
  583. * @returns 0 on success.
  584. */
  585. int CNAME(BLASLONG bm, BLASLONG bn, BLASLONG bk, FLOAT alpha,
  586. FLOAT *restrict ba, FLOAT *restrict bb,
  587. FLOAT *restrict C, BLASLONG ldc
  588. #ifdef TRMMKERNEL
  589. , BLASLONG offset
  590. #endif
  591. )
  592. {
  593. if ( (bm == 0) || (bn == 0) || (bk == 0) || (alpha == ZERO))
  594. return 0;
  595. /*
  596. * interface code allocates buffers for ba and bb at page
  597. * granularity (i.e., using mmap(MAP_ANONYMOUS), so enable the compiler
  598. * to make use of the fact in vector load operations.
  599. */
  600. ba = __builtin_assume_aligned(ba, 16);
  601. bb = __builtin_assume_aligned(bb, 16);
  602. /*
  603. * Use offset and off even when compiled as SGEMMKERNEL to simplify
  604. * function signatures and function calls.
  605. */
  606. #ifndef TRMMKERNEL
  607. BLASLONG const offset = 0;
  608. #endif
  609. /*
  610. * Partition B and C into blocks of n_r (unroll_n) columns, called B_i
  611. * and C_i. For each partition, calculate C_i += alpha * (A * B_j).
  612. *
  613. * For remaining columns that do not fill up a block of n_r, iteratively
  614. * use smaller block sizes of powers of 2.
  615. */
  616. BLASLONG col = 0;
  617. for (BLASLONG block_size = unroll_n; block_size > 0; block_size /= 2)
  618. for (; bn - col >= block_size; col += block_size)
  619. GEBP_column_block(block_size, col, ba, bk, bb, bm, C, ldc, alpha, offset);
  620. return 0;
  621. }