Line data Source code
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 <stdio.h>
10 : #include <stdlib.h>
11 : #include <string.h>
12 :
13 : #include "dbm_mpi.h"
14 :
15 : #if defined(__parallel)
16 : /*******************************************************************************
17 : * \brief Check given MPI status and upon failure abort with a nice message.
18 : * \author Ole Schuett
19 : ******************************************************************************/
20 : #define CHECK(status) \
21 : if (status != MPI_SUCCESS) { \
22 : fprintf(stderr, "MPI error in %s:%i\n", __FILE__, __LINE__); \
23 : MPI_Abort(MPI_COMM_WORLD, EXIT_FAILURE); \
24 : }
25 : #endif
26 :
27 : /*******************************************************************************
28 : * \brief Wrapper around MPI_Init.
29 : * \author Ole Schuett
30 : ******************************************************************************/
31 0 : void dbm_mpi_init(int *argc, char ***argv) {
32 : #if defined(__parallel)
33 0 : CHECK(MPI_Init(argc, argv));
34 : #else
35 : (void)argc; // mark used
36 : (void)argv;
37 : #endif
38 0 : }
39 :
40 : /*******************************************************************************
41 : * \brief Wrapper around MPI_Finalize.
42 : * \author Ole Schuett
43 : ******************************************************************************/
44 0 : void dbm_mpi_finalize() {
45 : #if defined(__parallel)
46 0 : CHECK(MPI_Finalize());
47 : #endif
48 0 : }
49 :
50 : /*******************************************************************************
51 : * \brief Returns MPI_COMM_WORLD.
52 : * \author Ole Schuett
53 : ******************************************************************************/
54 0 : dbm_mpi_comm_t dbm_mpi_get_comm_world() {
55 : #if defined(__parallel)
56 0 : return MPI_COMM_WORLD;
57 : #else
58 : return -1;
59 : #endif
60 : }
61 :
62 : /*******************************************************************************
63 : * \brief Wrapper around MPI_Comm_f2c.
64 : * \author Ole Schuett
65 : ******************************************************************************/
66 816488 : dbm_mpi_comm_t dbm_mpi_comm_f2c(const int fortran_comm) {
67 : #if defined(__parallel)
68 816488 : return MPI_Comm_f2c(fortran_comm);
69 : #else
70 : (void)fortran_comm; // mark used
71 : return -1;
72 : #endif
73 : }
74 :
75 : /*******************************************************************************
76 : * \brief Wrapper around MPI_Comm_c2f.
77 : * \author Ole Schuett
78 : ******************************************************************************/
79 0 : int dbm_mpi_comm_c2f(const dbm_mpi_comm_t comm) {
80 : #if defined(__parallel)
81 0 : return MPI_Comm_c2f(comm);
82 : #else
83 : (void)comm; // mark used
84 : return -1;
85 : #endif
86 : }
87 :
88 : /*******************************************************************************
89 : * \brief Wrapper around MPI_Comm_rank.
90 : * \author Ole Schuett
91 : ******************************************************************************/
92 2421486 : int dbm_mpi_comm_rank(const dbm_mpi_comm_t comm) {
93 : #if defined(__parallel)
94 2421486 : int rank;
95 2421486 : CHECK(MPI_Comm_rank(comm, &rank));
96 2421486 : return rank;
97 : #else
98 : (void)comm; // mark used
99 : return 0;
100 : #endif
101 : }
102 :
103 : /*******************************************************************************
104 : * \brief Wrapper around MPI_Comm_size.
105 : * \author Ole Schuett
106 : ******************************************************************************/
107 2421630 : int dbm_mpi_comm_size(const dbm_mpi_comm_t comm) {
108 : #if defined(__parallel)
109 2421630 : int nranks;
110 2421630 : CHECK(MPI_Comm_size(comm, &nranks));
111 2421630 : return nranks;
112 : #else
113 : (void)comm; // mark used
114 : return 1;
115 : #endif
116 : }
117 :
118 : /*******************************************************************************
119 : * \brief Wrapper around MPI_Dims_create.
120 : * \author Ole Schuett
121 : ******************************************************************************/
122 0 : void dbm_mpi_dims_create(const int nnodes, const int ndims, int dims[]) {
123 : #if defined(__parallel)
124 0 : CHECK(MPI_Dims_create(nnodes, ndims, dims));
125 : #else
126 : dims[0] = nnodes;
127 : for (int i = 1; i < ndims; i++) {
128 : dims[i] = 1;
129 : }
130 : #endif
131 0 : }
132 :
133 : /*******************************************************************************
134 : * \brief Wrapper around MPI_Cart_create.
135 : * \author Ole Schuett
136 : ******************************************************************************/
137 0 : dbm_mpi_comm_t dbm_mpi_cart_create(const dbm_mpi_comm_t comm_old,
138 : const int ndims, const int dims[],
139 : const int periods[], const int reorder) {
140 : #if defined(__parallel)
141 0 : dbm_mpi_comm_t comm_cart;
142 0 : CHECK(MPI_Cart_create(comm_old, ndims, dims, periods, reorder, &comm_cart));
143 0 : return comm_cart;
144 : #else
145 : (void)comm_old; // mark used
146 : (void)ndims;
147 : (void)dims;
148 : (void)periods;
149 : (void)reorder;
150 : return -1;
151 : #endif
152 : }
153 :
154 : /*******************************************************************************
155 : * \brief Wrapper around MPI_Cart_get.
156 : * \author Ole Schuett
157 : ******************************************************************************/
158 1614324 : void dbm_mpi_cart_get(const dbm_mpi_comm_t comm, int maxdims, int dims[],
159 : int periods[], int coords[]) {
160 : #if defined(__parallel)
161 1614324 : CHECK(MPI_Cart_get(comm, maxdims, dims, periods, coords));
162 : #else
163 : (void)comm; // mark used
164 : for (int i = 0; i < maxdims; i++) {
165 : dims[i] = 1;
166 : periods[i] = 1;
167 : coords[i] = 0;
168 : }
169 : #endif
170 1614324 : }
171 :
172 : /*******************************************************************************
173 : * \brief Wrapper around MPI_Cart_rank.
174 : * \author Ole Schuett
175 : ******************************************************************************/
176 105003511 : int dbm_mpi_cart_rank(const dbm_mpi_comm_t comm, const int coords[]) {
177 : #if defined(__parallel)
178 105003511 : int rank;
179 105003511 : CHECK(MPI_Cart_rank(comm, coords, &rank));
180 105003511 : return rank;
181 : #else
182 : (void)comm; // mark used
183 : (void)coords;
184 : return 0;
185 : #endif
186 : }
187 :
188 : /*******************************************************************************
189 : * \brief Wrapper around MPI_Cart_sub.
190 : * \author Ole Schuett
191 : ******************************************************************************/
192 1614324 : dbm_mpi_comm_t dbm_mpi_cart_sub(const dbm_mpi_comm_t comm,
193 : const int remain_dims[]) {
194 : #if defined(__parallel)
195 1614324 : dbm_mpi_comm_t newcomm;
196 1614324 : CHECK(MPI_Cart_sub(comm, remain_dims, &newcomm));
197 1614324 : return newcomm;
198 : #else
199 : (void)comm; // mark used
200 : (void)remain_dims;
201 : return -1;
202 : #endif
203 : }
204 :
205 : /*******************************************************************************
206 : * \brief Wrapper around MPI_Comm_free.
207 : * \author Ole Schuett
208 : ******************************************************************************/
209 1614324 : void dbm_mpi_comm_free(dbm_mpi_comm_t *comm) {
210 : #if defined(__parallel)
211 1614324 : CHECK(MPI_Comm_free(comm));
212 : #else
213 : (void)comm; // mark used
214 : #endif
215 1614324 : }
216 :
217 : /*******************************************************************************
218 : * \brief Wrapper around MPI_Comm_compare.
219 : * \author Ole Schuett
220 : ******************************************************************************/
221 404070 : bool dbm_mpi_comms_are_similar(const dbm_mpi_comm_t comm1,
222 : const dbm_mpi_comm_t comm2) {
223 : #if defined(__parallel)
224 404070 : int res;
225 404070 : CHECK(MPI_Comm_compare(comm1, comm2, &res));
226 404070 : return res == MPI_IDENT || res == MPI_CONGRUENT || res == MPI_SIMILAR;
227 : #else
228 : (void)comm1; // mark used
229 : (void)comm2;
230 : return true;
231 : #endif
232 : }
233 :
234 : /*******************************************************************************
235 : * \brief Wrapper around MPI_Allreduce for op MPI_MAX and datatype MPI_INT.
236 : * \author Ole Schuett
237 : ******************************************************************************/
238 807852 : void dbm_mpi_max_int(int *values, const int count, const dbm_mpi_comm_t comm) {
239 : #if defined(__parallel)
240 807852 : int value = 0;
241 807852 : void *recvbuf = (0 < count ? dbm_mpi_alloc_mem(count * sizeof(int)) : &value);
242 807852 : CHECK(MPI_Allreduce(values, recvbuf, count, MPI_INT, MPI_MAX, comm));
243 807852 : memcpy(values, recvbuf, count * sizeof(int));
244 807852 : if (0 < count) {
245 807852 : dbm_mpi_free_mem(recvbuf);
246 : }
247 : #else
248 : (void)comm; // mark used
249 : (void)values;
250 : (void)count;
251 : #endif
252 807852 : }
253 :
254 : /*******************************************************************************
255 : * \brief Wrapper around MPI_Allreduce for op MPI_MAX and datatype MPI_UINT64_T.
256 : * \author Ole Schuett
257 : ******************************************************************************/
258 19012 : void dbm_mpi_max_uint64(uint64_t *values, const int count,
259 : const dbm_mpi_comm_t comm) {
260 : #if defined(__parallel)
261 19012 : uint64_t value = 0;
262 38024 : void *recvbuf =
263 19012 : (0 < count ? dbm_mpi_alloc_mem(count * sizeof(uint64_t)) : &value);
264 19012 : CHECK(MPI_Allreduce(values, recvbuf, count, MPI_UINT64_T, MPI_MAX, comm));
265 19012 : memcpy(values, recvbuf, count * sizeof(uint64_t));
266 19012 : if (0 < count) {
267 19012 : dbm_mpi_free_mem(recvbuf);
268 : }
269 : #else
270 : (void)comm; // mark used
271 : (void)values;
272 : (void)count;
273 : #endif
274 19012 : }
275 :
276 : /*******************************************************************************
277 : * \brief Wrapper around MPI_Allreduce for op MPI_MAX and datatype MPI_DOUBLE.
278 : * \author Ole Schuett
279 : ******************************************************************************/
280 48 : void dbm_mpi_max_double(double *values, const int count,
281 : const dbm_mpi_comm_t comm) {
282 : #if defined(__parallel)
283 48 : double value = 0;
284 96 : void *recvbuf =
285 48 : (0 < count ? dbm_mpi_alloc_mem(count * sizeof(double)) : &value);
286 48 : CHECK(MPI_Allreduce(values, recvbuf, count, MPI_DOUBLE, MPI_MAX, comm));
287 48 : memcpy(values, recvbuf, count * sizeof(double));
288 48 : if (0 < count) {
289 48 : dbm_mpi_free_mem(recvbuf);
290 : }
291 : #else
292 : (void)comm; // mark used
293 : (void)values;
294 : (void)count;
295 : #endif
296 48 : }
297 :
298 : /*******************************************************************************
299 : * \brief Wrapper around MPI_Allreduce for op MPI_SUM and datatype MPI_INT.
300 : * \author Ole Schuett
301 : ******************************************************************************/
302 201963 : void dbm_mpi_sum_int(int *values, const int count, const dbm_mpi_comm_t comm) {
303 : #if defined(__parallel)
304 201963 : int value = 0;
305 201963 : void *recvbuf = (0 < count ? dbm_mpi_alloc_mem(count * sizeof(int)) : &value);
306 201963 : CHECK(MPI_Allreduce(values, recvbuf, count, MPI_INT, MPI_SUM, comm));
307 201963 : memcpy(values, recvbuf, count * sizeof(int));
308 201963 : if (0 < count) {
309 201963 : dbm_mpi_free_mem(recvbuf);
310 : }
311 : #else
312 : (void)comm; // mark used
313 : (void)values;
314 : (void)count;
315 : #endif
316 201963 : }
317 :
318 : /*******************************************************************************
319 : * \brief Wrapper around MPI_Allreduce for op MPI_SUM and datatype MPI_INT64_T.
320 : * \author Ole Schuett
321 : ******************************************************************************/
322 798827 : void dbm_mpi_sum_int64(int64_t *values, const int count,
323 : const dbm_mpi_comm_t comm) {
324 : #if defined(__parallel)
325 798827 : int64_t value = 0;
326 1597654 : void *recvbuf =
327 798827 : (0 < count ? dbm_mpi_alloc_mem(count * sizeof(int64_t)) : &value);
328 798827 : CHECK(MPI_Allreduce(values, recvbuf, count, MPI_INT64_T, MPI_SUM, comm));
329 798827 : memcpy(values, recvbuf, count * sizeof(int64_t));
330 798827 : if (0 < count) {
331 798827 : dbm_mpi_free_mem(recvbuf);
332 : }
333 : #else
334 : (void)comm; // mark used
335 : (void)values;
336 : (void)count;
337 : #endif
338 798827 : }
339 :
340 : /*******************************************************************************
341 : * \brief Wrapper around MPI_Allreduce for op MPI_SUM and datatype MPI_UINT64_T.
342 : * \author Hans Pabst
343 : ******************************************************************************/
344 0 : void dbm_mpi_sum_uint64(uint64_t *values, const int count,
345 : const dbm_mpi_comm_t comm) {
346 : #if defined(__parallel)
347 0 : uint64_t value = 0;
348 0 : void *recvbuf =
349 0 : (0 < count ? dbm_mpi_alloc_mem(count * sizeof(uint64_t)) : &value);
350 0 : CHECK(MPI_Allreduce(values, recvbuf, count, MPI_UINT64_T, MPI_SUM, comm));
351 0 : memcpy(values, recvbuf, count * sizeof(uint64_t));
352 0 : if (0 < count) {
353 0 : dbm_mpi_free_mem(recvbuf);
354 : }
355 : #else
356 : (void)comm; // mark used
357 : (void)values;
358 : (void)count;
359 : #endif
360 0 : }
361 :
362 : /*******************************************************************************
363 : * \brief Wrapper around MPI_Allreduce for op MPI_SUM and datatype MPI_DOUBLE.
364 : * \author Ole Schuett
365 : ******************************************************************************/
366 190 : void dbm_mpi_sum_double(double *values, const int count,
367 : const dbm_mpi_comm_t comm) {
368 : #if defined(__parallel)
369 190 : double value = 0;
370 380 : void *recvbuf =
371 190 : (0 < count ? dbm_mpi_alloc_mem(count * sizeof(double)) : &value);
372 190 : CHECK(MPI_Allreduce(values, recvbuf, count, MPI_DOUBLE, MPI_SUM, comm));
373 190 : memcpy(values, recvbuf, count * sizeof(double));
374 190 : if (0 < count) {
375 190 : dbm_mpi_free_mem(recvbuf);
376 : }
377 : #else
378 : (void)comm; // mark used
379 : (void)values;
380 : (void)count;
381 : #endif
382 190 : }
383 :
384 : /*******************************************************************************
385 : * \brief Wrapper around MPI_Sendrecv for datatype MPI_BYTE.
386 : * \author Ole Schuett
387 : ******************************************************************************/
388 21800 : int dbm_mpi_sendrecv_byte(const void *sendbuf, const int sendcount,
389 : const int dest, const int sendtag, void *recvbuf,
390 : const int recvcount, const int source,
391 : const int recvtag, const dbm_mpi_comm_t comm) {
392 : #if defined(__parallel)
393 21800 : MPI_Status status;
394 21800 : CHECK(MPI_Sendrecv(sendbuf, sendcount, MPI_BYTE, dest, sendtag, recvbuf,
395 : recvcount, MPI_BYTE, source, recvtag, comm, &status))
396 21800 : int count_received;
397 21800 : CHECK(MPI_Get_count(&status, MPI_BYTE, &count_received));
398 21800 : return count_received;
399 : #else
400 : (void)sendbuf; // mark used
401 : (void)sendcount;
402 : (void)dest;
403 : (void)sendtag;
404 : (void)recvbuf;
405 : (void)recvcount;
406 : (void)source;
407 : (void)recvtag;
408 : (void)comm;
409 : fprintf(stderr, "Error: dbm_mpi_sendrecv_byte not available without MPI\n");
410 : abort();
411 : #endif
412 : }
413 :
414 : /*******************************************************************************
415 : * \brief Wrapper around MPI_Sendrecv for datatype MPI_DOUBLE.
416 : * \author Ole Schuett
417 : ******************************************************************************/
418 21800 : int dbm_mpi_sendrecv_double(const double *sendbuf, const int sendcount,
419 : const int dest, const int sendtag, double *recvbuf,
420 : const int recvcount, const int source,
421 : const int recvtag, const dbm_mpi_comm_t comm) {
422 : #if defined(__parallel)
423 21800 : MPI_Status status;
424 21800 : CHECK(MPI_Sendrecv(sendbuf, sendcount, MPI_DOUBLE, dest, sendtag, recvbuf,
425 : recvcount, MPI_DOUBLE, source, recvtag, comm, &status))
426 21800 : int count_received;
427 21800 : CHECK(MPI_Get_count(&status, MPI_DOUBLE, &count_received));
428 21800 : return count_received;
429 : #else
430 : (void)sendbuf; // mark used
431 : (void)sendcount;
432 : (void)dest;
433 : (void)sendtag;
434 : (void)recvbuf;
435 : (void)recvcount;
436 : (void)source;
437 : (void)recvtag;
438 : (void)comm;
439 : fprintf(stderr, "Error: dbm_mpi_sendrecv_double not available without MPI\n");
440 : abort();
441 : #endif
442 : }
443 :
444 : /*******************************************************************************
445 : * \brief Wrapper around MPI_Alltoall for datatype MPI_INT.
446 : * \author Ole Schuett
447 : ******************************************************************************/
448 851596 : void dbm_mpi_alltoall_int(const int *sendbuf, const int sendcount, int *recvbuf,
449 : const int recvcount, const dbm_mpi_comm_t comm) {
450 : #if defined(__parallel)
451 851596 : CHECK(MPI_Alltoall(sendbuf, sendcount, MPI_INT, recvbuf, recvcount, MPI_INT,
452 851596 : comm));
453 : #else
454 : (void)comm; // mark used
455 : assert(sendcount == recvcount);
456 : memcpy(recvbuf, sendbuf, sendcount * sizeof(int));
457 : #endif
458 851596 : }
459 :
460 : /*******************************************************************************
461 : * \brief Wrapper around MPI_Alltoallv for datatype MPI_BYTE.
462 : * \author Ole Schuett
463 : ******************************************************************************/
464 425726 : void dbm_mpi_alltoallv_byte(const void *sendbuf, const int *sendcounts,
465 : const int *sdispls, void *recvbuf,
466 : const int *recvcounts, const int *rdispls,
467 : const dbm_mpi_comm_t comm) {
468 : #if defined(__parallel)
469 425726 : CHECK(MPI_Alltoallv(sendbuf, sendcounts, sdispls, MPI_BYTE, recvbuf,
470 425726 : recvcounts, rdispls, MPI_BYTE, comm));
471 : #else
472 : (void)comm; // mark used
473 : assert(sendcounts[0] == recvcounts[0]);
474 : assert(sdispls[0] == 0 && rdispls[0] == 0);
475 : memcpy(recvbuf, sendbuf, sendcounts[0]);
476 : #endif
477 425726 : }
478 :
479 : /*******************************************************************************
480 : * \brief Wrapper around MPI_Alltoallv for datatype MPI_DOUBLE.
481 : * \author Ole Schuett
482 : ******************************************************************************/
483 425870 : void dbm_mpi_alltoallv_double(const double *sendbuf, const int *sendcounts,
484 : const int *sdispls, double *recvbuf,
485 : const int *recvcounts, const int *rdispls,
486 : const dbm_mpi_comm_t comm) {
487 : #if defined(__parallel)
488 425870 : CHECK(MPI_Alltoallv(sendbuf, sendcounts, sdispls, MPI_DOUBLE, recvbuf,
489 425870 : recvcounts, rdispls, MPI_DOUBLE, comm));
490 : #else
491 : (void)comm; // mark used
492 : assert(sendcounts[0] == recvcounts[0]);
493 : assert(sdispls[0] == 0 && rdispls[0] == 0);
494 : memcpy(recvbuf, sendbuf, sendcounts[0] * sizeof(double));
495 : #endif
496 425870 : }
497 :
498 : /*******************************************************************************
499 : * \brief Wrapper around MPI_Alloc_mem.
500 : * \author Hans Pabst
501 : ******************************************************************************/
502 3067644 : void *dbm_mpi_alloc_mem(size_t size) {
503 3067644 : void *result = NULL;
504 : #if defined(__parallel)
505 3067644 : CHECK(MPI_Alloc_mem((MPI_Aint)size, MPI_INFO_NULL, &result));
506 : #else
507 : result = malloc(size);
508 : #endif
509 3067644 : return result;
510 : }
511 :
512 : /*******************************************************************************
513 : * \brief Wrapper around MPI_Free_mem.
514 : * \author Hans Pabst
515 : ******************************************************************************/
516 3067644 : void dbm_mpi_free_mem(void *mem) {
517 : #if defined(__parallel)
518 3067644 : CHECK(MPI_Free_mem(mem));
519 : #else
520 : free(mem);
521 : #endif
522 3067644 : }
523 :
524 : // EOF
|