(git:ed6f26b)
Loading...
Searching...
No Matches
dbm_multiply.c
Go to the documentation of this file.
1/*----------------------------------------------------------------------------*/
2/* CP2K: A general program to perform molecular dynamics simulations */
3/* Copyright 2000-2025 CP2K developers group <https://cp2k.org> */
4/* */
5/* SPDX-License-Identifier: BSD-3-Clause */
6/*----------------------------------------------------------------------------*/
7
8#include <assert.h>
9#include <limits.h>
10#include <omp.h>
11#include <stdlib.h>
12#include <string.h>
13
14#include "../offload/offload_runtime.h"
15#include "dbm_hyperparams.h"
16#include "dbm_internal.h"
17#include "dbm_library.h"
18#include "dbm_mempool.h"
19#include "dbm_multiply.h"
20#include "dbm_multiply_comm.h"
21#include "dbm_multiply_cpu.h"
22#include "dbm_multiply_gpu.h"
23
24#if defined(__LIBXSMM)
25#include <libxsmm.h>
26#endif
27
28#if !defined(DBM_VALIDATE_AGAINST_LIBXSMM) && 0
29#define DBM_VALIDATE_AGAINST_LIBXSMM
30#endif
31
32/*******************************************************************************
33 * \brief Private routine for computing the max filter threshold for each row.
34 * \author Ole Schuett
35 ******************************************************************************/
36static float *compute_rows_max_eps(const bool trans, const dbm_matrix_t *matrix,
37 const double filter_eps) {
38 const int nrows = (trans) ? matrix->ncols : matrix->nrows;
39 int *nblocks_per_row = calloc(nrows, sizeof(int));
40 float *row_max_eps = malloc(nrows * sizeof(float));
41 assert(row_max_eps != NULL);
42
43#pragma omp parallel
44 {
45#pragma omp for
46 for (int ishard = 0; ishard < dbm_get_num_shards(matrix); ishard++) {
47 dbm_shard_t *shard = &matrix->shards[ishard];
48 for (int iblock = 0; iblock < shard->nblocks; iblock++) {
49 const dbm_block_t *blk = &shard->blocks[iblock];
50 const int row = (trans) ? blk->col : blk->row;
51#pragma omp atomic
52 nblocks_per_row[row]++;
53 }
54 }
55#pragma omp single
56 dbm_mpi_sum_int(nblocks_per_row, nrows, matrix->dist->comm);
57#pragma omp barrier
58#pragma omp for
59 for (int i = 0; i < nrows; i++) {
60 const float f =
61 ((float)filter_eps) / ((float)imax(1, nblocks_per_row[i]));
62 row_max_eps[i] = f * f;
63 }
64 } // end of omp parallel region
65
66 free(nblocks_per_row);
67 return row_max_eps; // Ownership of row_max_eps transfers to caller.
68}
69
70/*******************************************************************************
71 * \brief Private struct for storing the context of the multiplication backend.
72 * \author Ole Schuett
73 ******************************************************************************/
74typedef struct {
75#if defined(__OFFLOAD) && !defined(__NO_OFFLOAD_DBM)
76 dbm_multiply_gpu_context_t gpu;
77#endif
79
80/*******************************************************************************
81 * \brief Private routine for initializing the multiplication backend.
82 * \author Ole Schuett
83 ******************************************************************************/
85 backend_context_t *ctx = calloc(1, sizeof(backend_context_t));
86
87#if defined(__OFFLOAD) && !defined(__NO_OFFLOAD_DBM)
88 dbm_multiply_gpu_start(DBM_MAX_BATCH_SIZE, dbm_get_num_shards(matrix_c),
89 matrix_c->shards, &ctx->gpu);
90#else
91 (void)matrix_c; // mark as used
92#endif
93
94 return ctx;
95}
96
97/*******************************************************************************
98 * \brief Private routine for handing newly arrived packs to the backend.
99 * \author Ole Schuett
100 ******************************************************************************/
101static void backend_upload_packs(const dbm_pack_t *pack_a,
102 const dbm_pack_t *pack_b,
103 backend_context_t *ctx) {
104
105#if defined(__OFFLOAD) && !defined(__NO_OFFLOAD_DBM)
106 dbm_multiply_gpu_upload_packs(pack_a, pack_b, &ctx->gpu);
107#else
108 (void)pack_a; // mark as used
109 (void)pack_b;
110 (void)ctx;
111#endif
112}
113
114/*******************************************************************************
115 * \brief Private routine for sending a batch to the multiplication backend.
116 * \author Ole Schuett
117 ******************************************************************************/
118static void backend_process_batch(const int ntasks, dbm_task_t batch[ntasks],
119 const double alpha, const dbm_pack_t *pack_a,
120 const dbm_pack_t *pack_b, const int kshard,
121 dbm_shard_t *shard_c,
122 backend_context_t *ctx) {
123#if defined(__OFFLOAD) && !defined(__NO_OFFLOAD_DBM)
124 dbm_multiply_gpu_process_batch(ntasks, batch, alpha, kshard, &ctx->gpu);
125#if defined(DBM_VALIDATE_AGAINST_LIBXSMM) && defined(__LIBXSMM)
126 dbm_shard_gpu_t *const shard_g = &ctx->gpu.shards_c_dev[kshard];
127 dbm_shard_t shard_r;
129 /* start transferring GPU result to host */
130 assert(shard_c->data_size == shard_g->data_size);
131 dbm_shard_init(&shard_r);
132 dbm_shard_copy(&shard_r, shard_c);
133 offloadMemcpyAsyncDtoH(shard_c->data, shard_g->data,
134 shard_c->data_size * sizeof(double), shard_g->stream);
135 dbm_multiply_cpu_process_batch(ntasks, batch, alpha, pack_a, pack_b,
136 &shard_r);
137 /* finish transferring GPU result to host */
138 offloadStreamSynchronize(shard_g->stream);
139 libxsmm_matdiff_info diff;
140 libxsmm_matdiff_clear(&diff);
141 for (int itask = 0; itask < ntasks; ++itask) {
142 const dbm_task_t task = batch[itask];
143 const double *const tst = &shard_c->data[task.offset_c];
144 const double *const ref = &shard_r.data[task.offset_c];
145 libxsmm_matdiff_info d;
146 if (EXIT_SUCCESS == libxsmm_matdiff(&d, LIBXSMM_DATATYPE(double), task.m,
147 task.n, ref, tst, NULL /*ldref*/,
148 NULL /*ldtst*/)) {
149 libxsmm_matdiff_reduce(&diff, &d);
150 }
151 }
152 const char *const maxeps_env = getenv("DBM_MULTIPLY_MAXEPS");
153 const double maxeps = (NULL == maxeps_env ? 1E-13 : fabs(atof(maxeps_env)));
154 const double epsilon = libxsmm_matdiff_epsilon(&diff);
155 if (maxeps < epsilon) {
156 if (LIBXSMM_NOTNAN(diff.v_tst)) {
157 fprintf(stderr, "INFO ACC/LIBDBM: diff=%g (|%g-%g|=%g)\n", epsilon,
158 diff.v_ref, diff.v_tst, diff.linf_abs);
159 } else {
160 fprintf(stderr, "INFO ACC/LIBDBM: diff=%g\n", epsilon);
161 }
162 }
163 dbm_shard_release(&shard_r);
164#else
165 (void)pack_a;
166 (void)pack_b;
167 (void)shard_c; // mark as used
168#endif
169#else
170 (void)kshard;
171 (void)ctx; // mark as used
172 dbm_multiply_cpu_process_batch(ntasks, batch, alpha, pack_a, pack_b, shard_c);
173#endif
174}
175
176/*******************************************************************************
177 * \brief Private routine for downloading results of the multiplication backend.
178 * \author Ole Schuett
179 ******************************************************************************/
181#if defined(__OFFLOAD) && !defined(__NO_OFFLOAD_DBM)
182 dbm_multiply_gpu_download_results(&ctx->gpu);
183#else
184 (void)ctx; // mark as used
185#endif
186}
187
188/*******************************************************************************
189 * \brief Private routine for shutting down the multiplication backend.
190 * \author Ole Schuett
191 ******************************************************************************/
193#if defined(__OFFLOAD) && !defined(__NO_OFFLOAD_DBM)
194 dbm_multiply_gpu_stop(&ctx->gpu);
195#endif
196 free(ctx);
197}
198
199/*******************************************************************************
200 * \brief Private routine for multipling two packs.
201 * \author Ole Schuett
202 ******************************************************************************/
203static void multiply_packs(const bool transa, const bool transb,
204 const double alpha, const dbm_pack_t *pack_a,
205 const dbm_pack_t *pack_b,
206 const dbm_matrix_t *matrix_a,
207 const dbm_matrix_t *matrix_b, dbm_matrix_t *matrix_c,
208 const bool retain_sparsity,
209 const float *rows_max_eps, int64_t *flop,
210 backend_context_t *ctx) {
211 const float alpha2 = alpha * alpha;
212 int64_t flop_sum = 0;
213
214 const int nshard_rows = matrix_c->dist->rows.nshards;
215 const int nshard_cols = matrix_c->dist->cols.nshards;
216 int *shard_row_start = calloc(nshard_rows, sizeof(int));
217 int *shard_col_start = calloc(nshard_cols, sizeof(int));
218 assert(NULL != shard_row_start && NULL != shard_col_start);
219
220 const int *sum_index_sizes_a =
221 (transa) ? matrix_a->row_sizes : matrix_a->col_sizes;
222 const int *sum_index_sizes_b =
223 (transb) ? matrix_b->col_sizes : matrix_b->row_sizes;
224 const int *free_index_sizes_a =
225 (transa) ? matrix_a->col_sizes : matrix_a->row_sizes;
226 const int *free_index_sizes_b =
227 (transb) ? matrix_b->row_sizes : matrix_b->col_sizes;
228
229#pragma omp parallel reduction(+ : flop_sum)
230 {
231 // Thread-private array covering given work in piece-wise fashion.
232 dbm_task_t *batch =
234
235 // Blocks are ordered first by shard. Creating lookup tables of boundaries.
236#pragma omp for nowait
237 for (int iblock = 1; iblock < pack_a->nblocks; iblock++) {
238 const int shard_row = pack_a->blocks[iblock].free_index % nshard_rows;
239 const int prev_shard_row =
240 pack_a->blocks[iblock - 1].free_index % nshard_rows;
241 if (prev_shard_row != shard_row) {
242 shard_row_start[shard_row] = iblock;
243 }
244 }
245#pragma omp for
246 for (int jblock = 1; jblock < pack_b->nblocks; jblock++) {
247 const int shard_col = pack_b->blocks[jblock].free_index % nshard_cols;
248 const int prev_shard_col =
249 pack_b->blocks[jblock - 1].free_index % nshard_cols;
250 if (prev_shard_col != shard_col) {
251 shard_col_start[shard_col] = jblock;
252 }
253 }
254
255#pragma omp for collapse(2) DBM_OMP_SCHEDULE
256 for (int shard_row = 0; shard_row < nshard_rows; shard_row++) {
257 for (int shard_col = 0; shard_col < nshard_cols; shard_col++) {
258 const int ishard = shard_row * nshard_cols + shard_col;
259 dbm_shard_t *shard_c = &matrix_c->shards[ishard];
260 int ntasks = 0;
261
262 // Use a merge-join to find pairs of blocks with matching sum indices.
263 // This utilizes that blocks within a shard are ordered by sum_index.
264 const int iblock_start = shard_row_start[shard_row];
265 int jblock_start = shard_col_start[shard_col];
266 for (int iblock = iblock_start; iblock < pack_a->nblocks; iblock++) {
267 const dbm_pack_block_t *blk_a = &pack_a->blocks[iblock];
268 if (blk_a->free_index % nshard_rows != shard_row) {
269 break;
270 }
271 for (int jblock = jblock_start; jblock < pack_b->nblocks; jblock++) {
272 const dbm_pack_block_t *blk_b = &pack_b->blocks[jblock];
273 if (blk_b->free_index % nshard_cols != shard_col) {
274 break;
275 }
276 if (blk_a->sum_index < blk_b->sum_index) {
277 break;
278 }
279 if (blk_a->sum_index > blk_b->sum_index) {
280 jblock_start++;
281 continue;
282 }
283 // Found block pair with blk_a->sum_index == blk_b->sum_index.
284
285 // Check norms.
286 const float result_norm = alpha2 * blk_a->norm * blk_b->norm;
287 if (result_norm < rows_max_eps[blk_a->free_index]) {
288 continue;
289 }
290
291 // Check block sizes.
292 const int m = free_index_sizes_a[blk_a->free_index];
293 const int n = free_index_sizes_b[blk_b->free_index];
294 const int k = sum_index_sizes_a[blk_a->sum_index];
295 assert(m == matrix_c->row_sizes[blk_a->free_index]);
296 assert(n == matrix_c->col_sizes[blk_b->free_index]);
297 assert(k == sum_index_sizes_b[blk_b->sum_index]);
298
299 // Get C block.
300 const int row = blk_a->free_index, col = blk_b->free_index;
301 dbm_block_t *blk_c = dbm_shard_lookup(shard_c, row, col);
302 if (blk_c == NULL && retain_sparsity) {
303 continue;
304 } else if (blk_c == NULL) {
305 assert(dbm_get_shard_index(matrix_c, row, col) == ishard);
306 assert(dbm_get_stored_coordinates(matrix_c, row, col) ==
307 matrix_c->dist->my_rank);
308 blk_c = dbm_shard_promise_new_block(shard_c, row, col, m * n);
309 }
310
311 // Count flops.
312 const int64_t task_flops = 2LL * m * n * k;
313 if (task_flops == 0) {
314 continue;
315 }
316 flop_sum += task_flops;
318
319 // Add block multiplication to batch.
320 batch[ntasks].m = m;
321 batch[ntasks].n = n;
322 batch[ntasks].k = k;
323 batch[ntasks].offset_a = blk_a->offset;
324 batch[ntasks].offset_b = blk_b->offset;
325 batch[ntasks].offset_c = blk_c->offset;
326 ++ntasks;
327
328 if (ntasks == DBM_MAX_BATCH_SIZE) {
329 backend_process_batch(ntasks, batch, alpha, pack_a, pack_b,
330 ishard, shard_c, ctx);
331 ntasks = 0;
332 }
333 }
334 }
335 backend_process_batch(ntasks, batch, alpha, pack_a, pack_b, ishard,
336 shard_c, ctx);
337 }
338 }
339
341 }
342
343 free(shard_row_start);
344 free(shard_col_start);
345
346 *flop += flop_sum;
347}
348
349/*******************************************************************************
350 * \brief Performs a multiplication of two dbm_matrix_t matrices.
351 * See dbm_matrix.h for details.
352 * \author Ole Schuett
353 ******************************************************************************/
354void dbm_multiply(const bool transa, const bool transb, const double alpha,
355 const dbm_matrix_t *matrix_a, const dbm_matrix_t *matrix_b,
356 const double beta, dbm_matrix_t *matrix_c,
357 const bool retain_sparsity, const double filter_eps,
358 int64_t *flop) {
359
360 assert(omp_get_num_threads() == 1);
361
362 // Throughout the matrix multiplication code the "sum_index" and "free_index"
363 // denote the summation (aka dummy) and free index from the Einstein notation.
364 const int num_sum_index_a = (transa) ? matrix_a->nrows : matrix_a->ncols;
365 const int num_sum_index_b = (transb) ? matrix_b->ncols : matrix_b->nrows;
366 const int num_free_index_a = (transa) ? matrix_a->ncols : matrix_a->nrows;
367 const int num_free_index_b = (transb) ? matrix_b->nrows : matrix_b->ncols;
368
369 // Sanity check matrix dimensions.
370 assert(num_sum_index_a == num_sum_index_b);
371 assert(num_free_index_a == matrix_c->nrows);
372 assert(num_free_index_b == matrix_c->ncols);
373
374 // Prepare matrix_c.
375 dbm_scale(matrix_c, beta);
376
377 // Start uploading matrix_c to the GPU.
378 backend_context_t *ctx = backend_start(matrix_c);
379
380 // Compute filter thresholds for each row.
381 float *rows_max_eps = compute_rows_max_eps(transa, matrix_a, filter_eps);
382
383 // Redistribute matrix_a and matrix_b across MPI ranks.
384 dbm_comm_iterator_t *iter =
385 dbm_comm_iterator_start(transa, transb, matrix_a, matrix_b, matrix_c);
386
387 // Main loop.
388 *flop = 0;
389 dbm_pack_t *pack_a, *pack_b;
390 while (dbm_comm_iterator_next(iter, &pack_a, &pack_b)) {
391 backend_upload_packs(pack_a, pack_b, ctx);
392 multiply_packs(transa, transb, alpha, pack_a, pack_b, matrix_a, matrix_b,
393 matrix_c, retain_sparsity, rows_max_eps, flop, ctx);
394 }
395
396 // Start downloading matrix_c from the GPU.
398
399 // Wait for all other MPI ranks to complete, then release ressources.
401 free(rows_max_eps);
402 backend_stop(ctx);
403
404 // Final filter pass.
405 dbm_filter(matrix_c, filter_eps);
406}
407
408// EOF
#define DBM_MAX_BATCH_SIZE
static int imax(int x, int y)
Returns the larger of two given integers (missing from the C standard)
void dbm_library_counter_increment(const int m, const int n, const int k)
Add given block multiplication to stats. This routine is thread-safe.
static int dbm_get_shard_index(const dbm_matrix_t *matrix, const int row, const int col)
Internal routine for getting a block's shard index.
Definition dbm_matrix.h:240
static int dbm_get_num_shards(const dbm_matrix_t *matrix)
Internal routine that returns the number of shards for given matrix.
Definition dbm_matrix.h:232
void dbm_mempool_host_free(const void *memory)
Internal routine for releasing memory back to the pool.
void * dbm_mempool_host_malloc(size_t size)
Internal routine for allocating host memory from the pool.
void dbm_mpi_sum_int(int *values, const int count, const dbm_mpi_comm_t comm)
Wrapper around MPI_Allreduce for op MPI_SUM and datatype MPI_INT.
Definition dbm_mpi.c:305
static void backend_upload_packs(const dbm_pack_t *pack_a, const dbm_pack_t *pack_b, backend_context_t *ctx)
Private routine for handing newly arrived packs to the backend.
static float * compute_rows_max_eps(const bool trans, const dbm_matrix_t *matrix, const double filter_eps)
Private routine for computing the max filter threshold for each row.
static backend_context_t * backend_start(const dbm_matrix_t *matrix_c)
Private routine for initializing the multiplication backend.
static void backend_process_batch(const int ntasks, dbm_task_t batch[ntasks], const double alpha, const dbm_pack_t *pack_a, const dbm_pack_t *pack_b, const int kshard, dbm_shard_t *shard_c, backend_context_t *ctx)
Private routine for sending a batch to the multiplication backend.
static void backend_stop(backend_context_t *ctx)
Private routine for shutting down the multiplication backend.
static void backend_download_results(backend_context_t *ctx)
Private routine for downloading results of the multiplication backend.
static void multiply_packs(const bool transa, const bool transb, const double alpha, const dbm_pack_t *pack_a, const dbm_pack_t *pack_b, const dbm_matrix_t *matrix_a, const dbm_matrix_t *matrix_b, dbm_matrix_t *matrix_c, const bool retain_sparsity, const float *rows_max_eps, int64_t *flop, backend_context_t *ctx)
Private routine for multipling two packs.
dbm_comm_iterator_t * dbm_comm_iterator_start(const bool transa, const bool transb, const dbm_matrix_t *matrix_a, const dbm_matrix_t *matrix_b, const dbm_matrix_t *matrix_c)
Internal routine for creating a communication iterator.
void dbm_comm_iterator_stop(dbm_comm_iterator_t *iter)
Internal routine for releasing the given communication iterator.
bool dbm_comm_iterator_next(dbm_comm_iterator_t *iter, dbm_pack_t **pack_a, dbm_pack_t **pack_b)
Internal routine for retriving next pair of packs from given iterator.
void dbm_multiply_cpu_process_batch(const int ntasks, dbm_task_t batch[ntasks], const double alpha, const dbm_pack_t *pack_a, const dbm_pack_t *pack_b, dbm_shard_t *shard_c)
Private hash function based on Szudzik's elegant pairing. Using unsigned int to return a positive num...
void dbm_shard_release(dbm_shard_t *shard)
Internal routine for releasing a shard.
Definition dbm_shard.c:125
dbm_block_t * dbm_shard_promise_new_block(dbm_shard_t *shard, const int row, const int col, const int block_size)
Internal routine for allocating the metadata of a new block.
Definition dbm_shard.c:197
dbm_block_t * dbm_shard_lookup(const dbm_shard_t *shard, const int row, const int col)
Internal routine for looking up a block from a shard.
Definition dbm_shard.c:175
void dbm_shard_allocate_promised_blocks(dbm_shard_t *shard)
Internal routine for allocating and zeroing any promised block's data.
Definition dbm_shard.c:231
void dbm_shard_init(dbm_shard_t *shard)
Internal routine for initializing a shard.
Definition dbm_shard.c:64
void dbm_shard_copy(dbm_shard_t *shard_a, const dbm_shard_t *shard_b)
Internal routine for copying content of shard_b into shard_a.
Definition dbm_shard.c:80
static void const int const int i
Private struct for storing the context of the multiplication backend.
Internal struct for storing a block's metadata.
Definition dbm_shard.h:20
Internal struct for storing a communication iterator.
dbm_mpi_comm_t comm
Internal struct for storing a matrix.
Definition dbm_matrix.h:20
int * row_sizes
Definition dbm_matrix.h:25
int * col_sizes
Definition dbm_matrix.h:26
dbm_shard_t * shards
Definition dbm_matrix.h:28
dbm_distribution_t * dist
Definition dbm_matrix.h:21
Internal struct for storing a dbm_block_t plus its norm.
Internal struct for storing a pack - essentially a shard for MPI.
dbm_pack_block_t * blocks
Internal struct for storing a matrix shard.
Definition dbm_shard.h:30
double * data
Definition dbm_shard.h:43
dbm_block_t * blocks
Definition dbm_shard.h:33
int data_size
Definition dbm_shard.h:42
Internal struct for storing a task, ie. a single block multiplication.