(git:e966546)
Loading...
Searching...
No Matches
dbm_miniapp.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#include "dbm_library.h"
8#include "dbm_matrix.h"
9
10#include "../mpiwrap/cp_mpi.h"
11#include "../offload/offload_library.h"
12#include "../offload/offload_mempool.h"
13
14#include <assert.h>
15#include <omp.h>
16#include <stdint.h>
17#include <stdio.h>
18#include <stdlib.h>
19#include <string.h>
20
21#if defined(__LIBXSMM)
22#include <libxsmm.h>
23#endif
24
25/*******************************************************************************
26 * \brief Wrapper for printf, passed to dbm_library_print_stats.
27 * \author Ole Schuett
28 ******************************************************************************/
29static void print_func(const char *msg, int /*msglen*/, int output_unit) {
30 if (output_unit == 0) { // i.e. my_rank == 0
31 printf("%s", msg);
32 }
33}
34
35/*******************************************************************************
36 * \brief Returns the smaller of the two integers (missing from the C standard).
37 * \author Ole Schuett
38 ******************************************************************************/
39static inline int imin(int x, int y) { return (x < y ? x : y); }
40
41/*******************************************************************************
42 * \brief Private routine for creating a distribution.
43 * \author Ole Schuett
44 ******************************************************************************/
45static dbm_distribution_t *create_dist(const int nrows, const int ncols,
46 const cp_mpi_comm_t comm) {
47 int cart_dims[2], cart_periods[2], cart_coords[2];
48 cp_mpi_cart_get(comm, 2, cart_dims, cart_periods, cart_coords);
49
50 // Create distribution.
51 assert(0 < nrows && 0 < ncols);
52 int *row_dist = malloc(nrows * sizeof(int));
53 int *col_dist = malloc(ncols * sizeof(int));
54 assert(row_dist != NULL && col_dist != NULL);
55 for (int i = 0; i < nrows; i++) {
56 row_dist[i] = i % cart_dims[0];
57 }
58 for (int i = 0; i < ncols; i++) {
59 col_dist[i] = i % cart_dims[1];
60 }
61 const int fortran_comm = cp_mpi_comm_c2f(comm);
62 dbm_distribution_t *dist = NULL;
63 dbm_distribution_new(&dist, fortran_comm, nrows, ncols, row_dist, col_dist);
64 free(row_dist);
65 free(col_dist);
66 return dist;
67}
68
69/*******************************************************************************
70 * \brief Private routine for creating a distribution and an empty matrix.
71 * \author Ole Schuett
72 ******************************************************************************/
73static dbm_matrix_t *
74create_some_matrix(const int nrows, const int ncols, const int nrows_min,
75 const int nrows_max, const int ncols_min,
76 const int ncols_max, const cp_mpi_comm_t comm) {
77 // Create distribution.
78 dbm_distribution_t *dist = create_dist(nrows, ncols, comm);
79
80 // Create matrix.
81 assert(0 < nrows && 0 < ncols);
82 int *row_sizes = malloc(nrows * sizeof(int));
83 int *col_sizes = malloc(ncols * sizeof(int));
84 assert(row_sizes != NULL && col_sizes != NULL);
85 assert(0 < nrows_min && nrows_min <= nrows_max);
86 assert(0 < ncols_min && ncols_min <= ncols_max);
87 if (nrows_min != nrows_max) {
88 const int row_size = nrows_max - nrows_min + 1;
89 for (int i = 0; i < nrows; i++) {
90 row_sizes[i] = rand() % row_size + 1;
91 }
92 } else {
93 for (int i = 0; i < nrows; i++) {
94 row_sizes[i] = nrows_max;
95 }
96 }
97 if (ncols_min != ncols_max) {
98 const int col_size = ncols_max - ncols_min + 1;
99 for (int i = 0; i < ncols; i++) {
100 col_sizes[i] = rand() % col_size + 1;
101 }
102 } else {
103 for (int i = 0; i < ncols; i++) {
104 col_sizes[i] = ncols_max;
105 }
106 }
107 dbm_matrix_t *matrix = NULL;
108 dbm_create(&matrix, dist, "some name", nrows, ncols, row_sizes, col_sizes);
109 dbm_distribution_release(dist);
110 free(row_sizes);
111 free(col_sizes);
112 return matrix;
113}
114
115/*******************************************************************************
116 * \brief Private routine for reserving all blocks of the given matrix.
117 * \author Ole Schuett
118 ******************************************************************************/
119static void reserve_all_blocks(dbm_matrix_t *matrix) {
120 int nrows, ncols;
121 const int *row_sizes, *col_sizes;
122 dbm_get_row_sizes(matrix, &nrows, &row_sizes);
123 dbm_get_col_sizes(matrix, &ncols, &col_sizes);
124
125#pragma omp parallel
126 {
127 int nblocks = 0;
128#pragma omp for collapse(2)
129 for (int row = 0; row < nrows; row++) {
130 for (int col = 0; col < ncols; col++) {
131 if (dbm_get_stored_coordinates(matrix, row, col) ==
132 matrix->dist->my_rank) {
133 ++nblocks;
134 }
135 }
136 }
137 assert(0 < nblocks);
138 int *reserve_row = malloc(nblocks * sizeof(int));
139 int *reserve_col = malloc(nblocks * sizeof(int));
140 assert(reserve_row != NULL && reserve_col != NULL);
141 int iblock = 0;
142#pragma omp for collapse(2)
143 for (int row = 0; row < nrows; row++) {
144 for (int col = 0; col < ncols; col++) {
145 if (dbm_get_stored_coordinates(matrix, row, col) ==
146 matrix->dist->my_rank) {
147 reserve_row[iblock] = row;
148 reserve_col[iblock] = col;
149 iblock++;
150 }
151 }
152 }
153 assert(iblock == nblocks);
154 dbm_reserve_blocks(matrix, nblocks, reserve_row, reserve_col);
155 free(reserve_row);
156 free(reserve_col);
157 }
158}
159
160/*******************************************************************************
161 * \brief Private routine for setting all blocks.
162 * \author Ole Schuett
163 ******************************************************************************/
164static void set_all_blocks(dbm_matrix_t *matrix) {
165#pragma omp parallel
166 {
167 dbm_iterator_t *iter = NULL;
168 dbm_iterator_start(&iter, matrix);
169 while (dbm_iterator_blocks_left(iter)) {
170 int row, col, row_size, col_size;
171 double *block;
172 dbm_iterator_next_block(iter, &row, &col, &block, &row_size, &col_size);
173 const int block_size = row_size * col_size;
174 for (int i = 0; i < block_size; i++) {
175 block[i] = 1.0 / (i + 1);
176 }
177 }
178 dbm_iterator_stop(iter);
179 }
180}
181/*******************************************************************************
182 * \brief Run a benchmark of dbm_multiply with given block sizes.
183 * \author Ole Schuett
184 ******************************************************************************/
185void benchmark_multiply(const int M, const int N, const int K, const int m,
186 const int n, const int k, const cp_mpi_comm_t comm) {
187 dbm_matrix_t *matrix_a = create_some_matrix(M, K, 1, m, k, k, comm);
188 dbm_matrix_t *matrix_b = create_some_matrix(K, N, k, k, 1, n, comm);
189 dbm_distribution_t *dist_c = create_dist(M, N, comm);
190 dbm_matrix_t *matrix_c = NULL, *matrix_d = NULL;
191 dbm_create(&matrix_c, dist_c, "result", M, N, matrix_a->row_sizes,
192 matrix_b->col_sizes);
193 dbm_distribution_release(dist_c);
194
195 reserve_all_blocks(matrix_a);
196 reserve_all_blocks(matrix_b);
197 set_all_blocks(matrix_a);
198 set_all_blocks(matrix_b);
199
200 const char *const verify_env = getenv("DBM_MULTIPLY_VERIFY");
201 const int skip_verify = (NULL == verify_env ? 0 : (atoi(verify_env) + 1));
202
203 if (0 == skip_verify) {
204 dbm_distribution_t *const dist_shared = matrix_c->dist;
205 dbm_create(&matrix_d, dist_shared, matrix_c->name, matrix_c->nrows,
206 matrix_c->ncols, matrix_c->row_sizes, matrix_c->col_sizes);
207 dbm_copy(matrix_d, matrix_c);
208 }
209
210 int64_t flop = 0;
211 const double time_start_multiply = omp_get_wtime();
212 dbm_multiply(false, false, 1.0, matrix_a, matrix_b, 1.0, matrix_c, false,
213 1e-8, &flop);
214 const double time_end_multiply = omp_get_wtime();
215
216 if (cp_mpi_comm_rank(comm) == 0) {
217 printf("%5i x %5i x %5i with %3i x %3i x %3i blocks: ", M, N, K, m, n, k);
218 }
219
220 if (NULL != matrix_d) { // Calculate result on the host for validation.
221 dbm_multiply(false, false, 1.0, matrix_a, matrix_b, 1.0, matrix_d, false,
222 1e-8, NULL);
223
224 const double maxeps = 1E-5, epsilon = dbm_maxeps(matrix_d, matrix_c);
225 if (maxeps < epsilon) {
226 printf("ERROR\n");
227 fprintf(stderr, "Failed validation (epsilon=%f).\n", epsilon);
228 exit(1);
229 }
230 dbm_release(matrix_d);
231 }
232
233 cp_mpi_sum_int64(&flop, 1, comm);
234 if (cp_mpi_comm_rank(comm) == 0) {
235 const double duration = time_end_multiply - time_start_multiply;
236 printf("%6.3f s => %6.1f GFLOP/s\n", duration, 1e-9 * flop / duration);
237 fflush(stdout);
238 }
239
240 dbm_release(matrix_a);
241 dbm_release(matrix_b);
242 dbm_release(matrix_c);
243}
244
245/*******************************************************************************
246 * \brief Stand-alone miniapp for smoke-testing and benchmarking dbm_multiply.
247 * \author Ole Schuett
248 ******************************************************************************/
249int main(int argc, char *argv[]) {
250 int result = EXIT_SUCCESS;
251
252 srand(25071975); // seed rng
253
254 cp_mpi_init(&argc, &argv);
255 dbm_library_init();
256
257 const cp_mpi_comm_t world_comm = cp_mpi_get_comm_world();
258 const int nranks = cp_mpi_comm_size(world_comm);
259 const int my_rank = cp_mpi_comm_rank(world_comm);
260
261 if (offload_get_device_count() > 0) {
262 offload_set_chosen_device(my_rank % offload_get_device_count());
263 }
264
265 // Create 2D cart.
266 int dims[2] = {0, 0};
267 cp_mpi_dims_create(nranks, 2, dims);
268 const int periods[2] = {true, true};
269 cp_mpi_comm_t comm = cp_mpi_cart_create(world_comm, 2, dims, periods, false);
270
271 if (my_rank == 0) {
272 printf("OpenMP-threads: %i GPUs: %i", omp_get_max_threads(),
273 imin(offload_get_device_count(), nranks));
274#if defined(__LIBXSMM)
275 printf(" Libxsmm: %s", LIBXSMM_VERSION);
276#else
277 printf(" Libxsmm: n/a");
278#endif
279#if defined(__parallel)
280 printf(" MPI-ranks: %i MPI-cart: %i x %i", nranks, dims[0], dims[1]);
281#else
282 printf(" MPI: n/a");
283#endif
284 printf("\n\n");
285 fflush(stdout);
286 }
287
288 if (1 >= argc) {
289 benchmark_multiply(16384, 128, 128, 4, 4, 4, comm);
290 benchmark_multiply(128, 16384, 128, 4, 4, 4, comm);
291 benchmark_multiply(128, 128, 16384, 4, 4, 4, comm);
292 benchmark_multiply(645, 645, 645, 4, 4, 4, comm);
293 if (my_rank == 0)
294 printf("\n");
295
296 benchmark_multiply(60, 500, 500, 128, 4, 4, comm);
297 benchmark_multiply(500, 60, 500, 4, 128, 4, comm);
298 benchmark_multiply(500, 500, 60, 4, 4, 128, comm);
299 if (my_rank == 0)
300 printf("\n");
301
302 benchmark_multiply(500, 60, 60, 4, 128, 128, comm);
303 benchmark_multiply(60, 500, 60, 128, 4, 128, comm);
304 benchmark_multiply(60, 60, 500, 128, 128, 4, comm);
305 if (my_rank == 0)
306 printf("\n");
307
308 benchmark_multiply(350, 350, 350, 23, 23, 23, comm);
309 benchmark_multiply(250, 250, 250, 32, 32, 32, comm);
310 benchmark_multiply(60, 60, 60, 128, 128, 128, comm);
311 } else { /* read triplet(s) from file or one triplet from command line */
312 FILE *const file = fopen(argv[1], "r"); /* try 1st arg as filename */
313 char buffer[1024];
314 const char delims[] = "x,;:|/\t ";
315 int mnk[] = {0, 0, 0}, i = 1, j = 0;
316 while (i < argc &&
317 (NULL == file || NULL != fgets(buffer, sizeof(buffer), file))) {
318 const char *arg = strtok(NULL != file ? buffer : argv[i], delims);
319 for (; NULL != arg && j < 3; arg = strtok(NULL, delims), ++j) {
320 mnk[j] = atoi(arg);
321 }
322 if (NULL != file) {
323 j = 0;
324 } else if (++i < argc) {
325 continue;
326 }
327 if (0 < mnk[0]) { /* valid MxNxK? */
328 const int m = mnk[0];
329 const int n = (0 < mnk[1] ? mnk[1] : m);
330 const int k = (0 < mnk[2] ? mnk[2] : m);
331 int M = (NULL == arg ? 0 : atoi(arg)), N, K;
332 if (0 < M) {
333 arg = strtok(NULL, delims);
334 N = (NULL == arg ? 1 : atoi(arg));
335 arg = strtok(NULL, delims);
336 K = (NULL == arg ? 1 : atoi(arg));
337 } else { /* default */
338 M = N = K = 128;
339 }
340 benchmark_multiply(M, N, K, m, n, k, comm);
341 mnk[0] = mnk[1] = mnk[2] = 0;
342 } else {
343 fprintf(stderr, "ERROR: invalid argument(s)\n");
344 result = EXIT_FAILURE;
345 i = argc; /* break */
346 }
347 }
348 if (NULL != file) {
349 fclose(file);
350 }
351 }
352
353 if (EXIT_SUCCESS == result) {
354 const int fortran_comm = cp_mpi_comm_c2f(comm);
355 dbm_library_print_stats(fortran_comm, &print_func, my_rank);
356 offload_mempool_stats_print(fortran_comm, &print_func, my_rank);
357 }
358 dbm_library_finalize();
359 cp_mpi_comm_free(&comm);
361 return result;
362}
363
364// EOF
int cp_mpi_comm_size(const cp_mpi_comm_t comm)
Wrapper around MPI_Comm_size.
Definition cp_mpi.c:110
cp_mpi_comm_t cp_mpi_cart_create(const cp_mpi_comm_t comm_old, const int ndims, const int dims[], const int periods[], const int reorder)
Wrapper around MPI_Cart_create.
Definition cp_mpi.c:140
void cp_mpi_init(int *argc, char ***argv)
Wrapper around MPI_Init.
Definition cp_mpi.c:34
int cp_mpi_comm_c2f(const cp_mpi_comm_t comm)
Wrapper around MPI_Comm_c2f.
Definition cp_mpi.c:82
cp_mpi_comm_t cp_mpi_get_comm_world(void)
Returns MPI_COMM_WORLD.
Definition cp_mpi.c:57
void cp_mpi_cart_get(const cp_mpi_comm_t comm, int maxdims, int dims[], int periods[], int coords[])
Wrapper around MPI_Cart_get.
Definition cp_mpi.c:161
void cp_mpi_sum_int64(int64_t *values, const int count, const cp_mpi_comm_t comm)
Wrapper around MPI_Allreduce for op MPI_SUM and datatype MPI_INT64_T.
Definition cp_mpi.c:325
void cp_mpi_comm_free(cp_mpi_comm_t *comm)
Wrapper around MPI_Comm_free.
Definition cp_mpi.c:212
void cp_mpi_dims_create(const int nnodes, const int ndims, int dims[])
Wrapper around MPI_Dims_create.
Definition cp_mpi.c:125
void cp_mpi_finalize(void)
Wrapper around MPI_Finalize.
Definition cp_mpi.c:47
int cp_mpi_comm_rank(const cp_mpi_comm_t comm)
Wrapper around MPI_Comm_rank.
Definition cp_mpi.c:95
int cp_mpi_comm_t
Definition cp_mpi.h:18
void dbm_get_col_sizes(const dbm_matrix_t *matrix, int *ncols, const int **col_sizes)
Returns the column block sizes of the given matrix.
Definition dbm_matrix.c:648
void dbm_get_row_sizes(const dbm_matrix_t *matrix, int *nrows, const int **row_sizes)
Returns the row block sizes of the given matrix.
Definition dbm_matrix.c:638
double dbm_maxeps(const dbm_matrix_t *matrix_a, const dbm_matrix_t *matrix_b)
Calculates maximum relative difference between matrix_a and matrix_b.
Definition dbm_matrix.c:566
static dbm_distribution_t * create_dist(const int nrows, const int ncols, const cp_mpi_comm_t comm)
Private routine for creating a distribution.
Definition dbm_miniapp.c:45
static int imin(int x, int y)
Returns the smaller of the two integers (missing from the C standard).
Definition dbm_miniapp.c:39
void benchmark_multiply(const int M, const int N, const int K, const int m, const int n, const int k, const cp_mpi_comm_t comm)
Run a benchmark of dbm_multiply with given block sizes.
static void set_all_blocks(dbm_matrix_t *matrix)
Private routine for setting all blocks.
static void print_func(const char *msg, int, int output_unit)
Wrapper for printf, passed to dbm_library_print_stats.
Definition dbm_miniapp.c:29
static dbm_matrix_t * create_some_matrix(const int nrows, const int ncols, const int nrows_min, const int nrows_max, const int ncols_min, const int ncols_max, const cp_mpi_comm_t comm)
Private routine for creating a distribution and an empty matrix.
Definition dbm_miniapp.c:74
static void reserve_all_blocks(dbm_matrix_t *matrix)
Private routine for reserving all blocks of the given matrix.
static void const int const int i
int main()
Unit test of the C-interface provided via libcp2k.h.
Internal struct for storing a two dimensional distribution.
Internal struct for storing a block iterator.
Definition dbm_matrix.h:34
Internal struct for storing a matrix.
Definition dbm_matrix.h:19
int * row_sizes
Definition dbm_matrix.h:24
char * name
Definition dbm_matrix.h:21
int * col_sizes
Definition dbm_matrix.h:25
dbm_distribution_t * dist
Definition dbm_matrix.h:20