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 : #include "cp_mpi.h"
8 :
9 : #include <assert.h>
10 : #include <omp.h>
11 : #include <stdio.h>
12 : #include <stdlib.h>
13 : #include <string.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(CMD) \
21 : do { \
22 : const int error = (CMD); \
23 : if (MPI_SUCCESS != error) { \
24 : fprintf(stderr, "MPI error #%i in %s:%i\n", error, __FILE__, __LINE__); \
25 : MPI_Abort(MPI_COMM_WORLD, EXIT_FAILURE); \
26 : } \
27 : } while (0)
28 : #endif
29 :
30 : /*******************************************************************************
31 : * \brief Wrapper around MPI_Init.
32 : * \author Ole Schuett
33 : ******************************************************************************/
34 0 : void cp_mpi_init(int *argc, char ***argv) {
35 : #if defined(__parallel)
36 0 : CHECK(MPI_Init(argc, argv));
37 : #else
38 0 : (void)argc; // mark used
39 0 : (void)argv;
40 : #endif
41 0 : }
42 :
43 : /*******************************************************************************
44 : * \brief Wrapper around MPI_Finalize.
45 : * \author Ole Schuett
46 : ******************************************************************************/
47 0 : void cp_mpi_finalize(void) {
48 : #if defined(__parallel)
49 0 : CHECK(MPI_Finalize());
50 : #endif
51 0 : }
52 :
53 : /*******************************************************************************
54 : * \brief Returns MPI_COMM_WORLD.
55 : * \author Ole Schuett
56 : ******************************************************************************/
57 0 : cp_mpi_comm_t cp_mpi_get_comm_world(void) {
58 : #if defined(__parallel)
59 0 : return MPI_COMM_WORLD;
60 : #else
61 0 : return -1;
62 : #endif
63 : }
64 :
65 : /*******************************************************************************
66 : * \brief Wrapper around MPI_Comm_f2c.
67 : * \author Ole Schuett
68 : ******************************************************************************/
69 859604 : cp_mpi_comm_t cp_mpi_comm_f2c(const int fortran_comm) {
70 : #if defined(__parallel)
71 : // Attempt to support "no MPI" if __parallel is defined (0 -> MPI_COMM_NULL).
72 859600 : return 0 != fortran_comm ? MPI_Comm_f2c(fortran_comm) : MPI_COMM_NULL;
73 : #else
74 4 : (void)fortran_comm; // mark used
75 4 : return -1;
76 : #endif
77 : }
78 :
79 : /*******************************************************************************
80 : * \brief Wrapper around MPI_Comm_c2f.
81 : * \author Ole Schuett
82 : ******************************************************************************/
83 0 : int cp_mpi_comm_c2f(const cp_mpi_comm_t comm) {
84 : #if defined(__parallel)
85 0 : return MPI_Comm_c2f(comm);
86 : #else
87 0 : (void)comm; // mark used
88 0 : return -1;
89 : #endif
90 : }
91 :
92 : /*******************************************************************************
93 : * \brief Wrapper around MPI_Comm_rank.
94 : * \author Ole Schuett
95 : ******************************************************************************/
96 2494152 : int cp_mpi_comm_rank(const cp_mpi_comm_t comm) {
97 : #if defined(__parallel)
98 2494152 : int rank = 0;
99 2494152 : if (MPI_COMM_NULL != comm) { // !MPI_Comm_compare
100 2494152 : CHECK(MPI_Comm_rank(comm, &rank));
101 : }
102 2494152 : return rank;
103 : #else
104 0 : (void)comm; // mark used
105 0 : return 0;
106 : #endif
107 : }
108 :
109 : /*******************************************************************************
110 : * \brief Wrapper around MPI_Comm_size.
111 : * \author Ole Schuett
112 : ******************************************************************************/
113 2494296 : int cp_mpi_comm_size(const cp_mpi_comm_t comm) {
114 : #if defined(__parallel)
115 2494296 : int nranks = 1;
116 2494296 : if (MPI_COMM_NULL != comm) { // !MPI_Comm_compare
117 2494296 : CHECK(MPI_Comm_size(comm, &nranks));
118 : }
119 2494296 : return nranks;
120 : #else
121 0 : (void)comm; // mark used
122 0 : return 1;
123 : #endif
124 : }
125 :
126 : /*******************************************************************************
127 : * \brief Wrapper around MPI_Dims_create.
128 : * \author Ole Schuett
129 : ******************************************************************************/
130 0 : void cp_mpi_dims_create(const int nnodes, const int ndims, int dims[]) {
131 : #if defined(__parallel)
132 0 : CHECK(MPI_Dims_create(nnodes, ndims, dims));
133 : #else
134 0 : dims[0] = nnodes;
135 0 : for (int i = 1; i < ndims; i++) {
136 0 : dims[i] = 1;
137 : }
138 : #endif
139 0 : }
140 :
141 : /*******************************************************************************
142 : * \brief Wrapper around MPI_Cart_create.
143 : * \author Ole Schuett
144 : ******************************************************************************/
145 0 : cp_mpi_comm_t cp_mpi_cart_create(const cp_mpi_comm_t comm_old, const int ndims,
146 : const int dims[], const int periods[],
147 : const int reorder) {
148 : #if defined(__parallel)
149 0 : cp_mpi_comm_t comm_cart;
150 0 : CHECK(MPI_Cart_create(comm_old, ndims, dims, periods, reorder, &comm_cart));
151 0 : return comm_cart;
152 : #else
153 0 : (void)comm_old; // mark used
154 0 : (void)ndims;
155 0 : (void)dims;
156 0 : (void)periods;
157 0 : (void)reorder;
158 0 : return -1;
159 : #endif
160 : }
161 :
162 : /*******************************************************************************
163 : * \brief Wrapper around MPI_Cart_get.
164 : * \author Ole Schuett
165 : ******************************************************************************/
166 1662768 : void cp_mpi_cart_get(const cp_mpi_comm_t comm, int maxdims, int dims[],
167 : int periods[], int coords[]) {
168 : #if defined(__parallel)
169 1662768 : CHECK(MPI_Cart_get(comm, maxdims, dims, periods, coords));
170 : #else
171 0 : (void)comm; // mark used
172 0 : for (int i = 0; i < maxdims; i++) {
173 0 : dims[i] = 1;
174 0 : periods[i] = 1;
175 0 : coords[i] = 0;
176 : }
177 : #endif
178 1662768 : }
179 :
180 : /*******************************************************************************
181 : * \brief Wrapper around MPI_Cart_rank.
182 : * \author Ole Schuett
183 : ******************************************************************************/
184 105044001 : int cp_mpi_cart_rank(const cp_mpi_comm_t comm, const int coords[]) {
185 : #if defined(__parallel)
186 105044001 : int rank;
187 105044001 : CHECK(MPI_Cart_rank(comm, coords, &rank));
188 105044001 : return rank;
189 : #else
190 0 : (void)comm; // mark used
191 0 : (void)coords;
192 0 : return 0;
193 : #endif
194 : }
195 :
196 : /*******************************************************************************
197 : * \brief Wrapper around MPI_Cart_sub.
198 : * \author Ole Schuett
199 : ******************************************************************************/
200 1662768 : cp_mpi_comm_t cp_mpi_cart_sub(const cp_mpi_comm_t comm,
201 : const int remain_dims[]) {
202 : #if defined(__parallel)
203 1662768 : cp_mpi_comm_t newcomm;
204 1662768 : CHECK(MPI_Cart_sub(comm, remain_dims, &newcomm));
205 1662768 : return newcomm;
206 : #else
207 0 : (void)comm; // mark used
208 0 : (void)remain_dims;
209 0 : return -1;
210 : #endif
211 : }
212 :
213 : /*******************************************************************************
214 : * \brief Wrapper around MPI_Comm_free.
215 : * \author Ole Schuett
216 : ******************************************************************************/
217 1662768 : void cp_mpi_comm_free(cp_mpi_comm_t *comm) {
218 : #if defined(__parallel)
219 1662768 : CHECK(MPI_Comm_free(comm));
220 : #else
221 0 : (void)comm; // mark used
222 : #endif
223 1662768 : }
224 :
225 : /*******************************************************************************
226 : * \brief Wrapper around MPI_Comm_compare.
227 : * \author Ole Schuett
228 : ******************************************************************************/
229 426648 : bool cp_mpi_comms_are_similar(const cp_mpi_comm_t comm1,
230 : const cp_mpi_comm_t comm2) {
231 : #if defined(__parallel)
232 426648 : int res;
233 426648 : CHECK(MPI_Comm_compare(comm1, comm2, &res));
234 426648 : return res == MPI_IDENT || res == MPI_CONGRUENT || res == MPI_SIMILAR;
235 : #else
236 0 : (void)comm1; // mark used
237 0 : (void)comm2;
238 0 : return true;
239 : #endif
240 : }
241 :
242 : /*******************************************************************************
243 : * \brief Wrapper around MPI_Allreduce for op MPI_MAX and datatype MPI_INT.
244 : * \author Ole Schuett
245 : ******************************************************************************/
246 853008 : void cp_mpi_max_int(int *values, const int count, const cp_mpi_comm_t comm) {
247 : #if defined(__parallel)
248 853008 : if (MPI_COMM_NULL != comm) { // !MPI_Comm_compare
249 853008 : int value = 0;
250 1706016 : void *recvbuf =
251 853008 : (1 < count ? cp_mpi_alloc_mem(count * sizeof(int)) : &value);
252 853008 : CHECK(MPI_Allreduce(values, recvbuf, count, MPI_INT, MPI_MAX, comm));
253 853008 : memcpy(values, recvbuf, count * sizeof(int));
254 853008 : if (1 < count) {
255 0 : cp_mpi_free_mem(recvbuf);
256 : }
257 : }
258 : #else
259 0 : (void)comm; // mark used
260 0 : (void)values;
261 0 : (void)count;
262 : #endif
263 853008 : }
264 :
265 : /*******************************************************************************
266 : * \brief Wrapper around MPI_Allreduce for op MPI_MAX and datatype MPI_UINT64_T.
267 : * \author Ole Schuett
268 : ******************************************************************************/
269 34870 : void cp_mpi_max_uint64(uint64_t *values, const int count,
270 : const cp_mpi_comm_t comm) {
271 : #if defined(__parallel)
272 34862 : if (MPI_COMM_NULL != comm) { // !MPI_Comm_compare
273 34862 : uint64_t value = 0;
274 69724 : void *recvbuf =
275 34862 : (1 < count ? cp_mpi_alloc_mem(count * sizeof(uint64_t)) : &value);
276 34862 : CHECK(MPI_Allreduce(values, recvbuf, count, MPI_UINT64_T, MPI_MAX, comm));
277 34862 : memcpy(values, recvbuf, count * sizeof(uint64_t));
278 34862 : if (1 < count) {
279 0 : cp_mpi_free_mem(recvbuf);
280 : }
281 : }
282 : #else
283 8 : (void)comm; // mark used
284 8 : (void)values;
285 8 : (void)count;
286 : #endif
287 34870 : }
288 :
289 : /*******************************************************************************
290 : * \brief Wrapper around MPI_Allreduce for op MPI_MAX and datatype MPI_DOUBLE.
291 : * \author Ole Schuett
292 : ******************************************************************************/
293 48 : void cp_mpi_max_double(double *values, const int count,
294 : const cp_mpi_comm_t comm) {
295 : #if defined(__parallel)
296 48 : if (MPI_COMM_NULL != comm) { // !MPI_Comm_compare
297 48 : double value = 0;
298 96 : void *recvbuf =
299 48 : (1 < count ? cp_mpi_alloc_mem(count * sizeof(double)) : &value);
300 48 : CHECK(MPI_Allreduce(values, recvbuf, count, MPI_DOUBLE, MPI_MAX, comm));
301 48 : memcpy(values, recvbuf, count * sizeof(double));
302 48 : if (1 < count) {
303 0 : cp_mpi_free_mem(recvbuf);
304 : }
305 : }
306 : #else
307 0 : (void)comm; // mark used
308 0 : (void)values;
309 0 : (void)count;
310 : #endif
311 48 : }
312 :
313 : /*******************************************************************************
314 : * \brief Wrapper around MPI_Allreduce for op MPI_SUM and datatype MPI_INT.
315 : * \author Ole Schuett
316 : ******************************************************************************/
317 213252 : void cp_mpi_sum_int(int *values, const int count, const cp_mpi_comm_t comm) {
318 : #if defined(__parallel)
319 213252 : if (MPI_COMM_NULL != comm) { // !MPI_Comm_compare
320 213252 : int value = 0;
321 426504 : void *recvbuf =
322 213252 : (1 < count ? cp_mpi_alloc_mem(count * sizeof(int)) : &value);
323 213252 : CHECK(MPI_Allreduce(values, recvbuf, count, MPI_INT, MPI_SUM, comm));
324 213252 : memcpy(values, recvbuf, count * sizeof(int));
325 213252 : if (1 < count) {
326 210629 : cp_mpi_free_mem(recvbuf);
327 : }
328 : }
329 : #else
330 0 : (void)comm; // mark used
331 0 : (void)values;
332 0 : (void)count;
333 : #endif
334 213252 : }
335 :
336 : /*******************************************************************************
337 : * \brief Wrapper around MPI_Allreduce for op MPI_SUM and datatype MPI_LONG.
338 : * \author Hans Pabst
339 : ******************************************************************************/
340 3762400 : void cp_mpi_sum_long(long *values, const int count, const cp_mpi_comm_t comm) {
341 : #if defined(__parallel)
342 3761600 : if (MPI_COMM_NULL != comm) { // !MPI_Comm_compare
343 3761600 : long value = 0;
344 7523200 : void *recvbuf =
345 3761600 : (1 < count ? cp_mpi_alloc_mem(count * sizeof(long)) : &value);
346 3761600 : CHECK(MPI_Allreduce(values, recvbuf, count, MPI_LONG, MPI_SUM, comm));
347 3761600 : memcpy(values, recvbuf, count * sizeof(long));
348 3761600 : if (1 < count) {
349 0 : cp_mpi_free_mem(recvbuf);
350 : }
351 : }
352 : #else
353 800 : (void)comm; // mark used
354 800 : (void)values;
355 800 : (void)count;
356 : #endif
357 3762400 : }
358 :
359 : /*******************************************************************************
360 : * \brief Wrapper around MPI_Allreduce for op MPI_SUM and datatype MPI_INT64_T.
361 : * \author Ole Schuett
362 : ******************************************************************************/
363 602112 : void cp_mpi_sum_int64(int64_t *values, const int count,
364 : const cp_mpi_comm_t comm) {
365 : #if defined(__parallel)
366 602112 : if (MPI_COMM_NULL != comm) { // !MPI_Comm_compare
367 602112 : int64_t value = 0;
368 1204224 : void *recvbuf =
369 602112 : (1 < count ? cp_mpi_alloc_mem(count * sizeof(int64_t)) : &value);
370 602112 : CHECK(MPI_Allreduce(values, recvbuf, count, MPI_INT64_T, MPI_SUM, comm));
371 602112 : memcpy(values, recvbuf, count * sizeof(int64_t));
372 602112 : if (1 < count) {
373 0 : cp_mpi_free_mem(recvbuf);
374 : }
375 : }
376 : #else
377 0 : (void)comm; // mark used
378 0 : (void)values;
379 0 : (void)count;
380 : #endif
381 602112 : }
382 :
383 : /*******************************************************************************
384 : * \brief Wrapper around MPI_Allreduce for op MPI_SUM and datatype MPI_DOUBLE.
385 : * \author Ole Schuett
386 : ******************************************************************************/
387 190 : void cp_mpi_sum_double(double *values, const int count,
388 : const cp_mpi_comm_t comm) {
389 : #if defined(__parallel)
390 190 : if (MPI_COMM_NULL != comm) { // !MPI_Comm_compare
391 190 : double value = 0;
392 380 : void *recvbuf =
393 190 : (1 < count ? cp_mpi_alloc_mem(count * sizeof(double)) : &value);
394 190 : CHECK(MPI_Allreduce(values, recvbuf, count, MPI_DOUBLE, MPI_SUM, comm));
395 190 : memcpy(values, recvbuf, count * sizeof(double));
396 190 : if (1 < count) {
397 0 : cp_mpi_free_mem(recvbuf);
398 : }
399 : }
400 : #else
401 0 : (void)comm; // mark used
402 0 : (void)values;
403 0 : (void)count;
404 : #endif
405 190 : }
406 :
407 : /*******************************************************************************
408 : * \brief Wrapper around MPI_Sendrecv for datatype MPI_BYTE.
409 : * \author Ole Schuett
410 : ******************************************************************************/
411 19134 : int cp_mpi_sendrecv_byte(const void *sendbuf, const int sendcount,
412 : const int dest, const int sendtag, void *recvbuf,
413 : const int recvcount, const int source,
414 : const int recvtag, const cp_mpi_comm_t comm) {
415 : #if defined(__parallel)
416 19134 : MPI_Status status;
417 19134 : CHECK(MPI_Sendrecv(sendbuf, sendcount, MPI_BYTE, dest, sendtag, recvbuf,
418 : recvcount, MPI_BYTE, source, recvtag, comm, &status));
419 19134 : int count_received = 0;
420 19134 : CHECK(MPI_Get_count(&status, MPI_BYTE, &count_received));
421 19134 : return count_received;
422 : #else
423 0 : (void)sendbuf; // mark used
424 0 : (void)sendcount;
425 0 : (void)dest;
426 0 : (void)sendtag;
427 0 : (void)recvbuf;
428 0 : (void)recvcount;
429 0 : (void)source;
430 0 : (void)recvtag;
431 0 : (void)comm;
432 0 : fprintf(stderr, "Error: cp_mpi_sendrecv_byte not available without MPI\n");
433 0 : abort();
434 : #endif
435 : }
436 :
437 : /*******************************************************************************
438 : * \brief Wrapper around MPI_Sendrecv for datatype MPI_DOUBLE.
439 : * \author Ole Schuett
440 : ******************************************************************************/
441 19134 : int cp_mpi_sendrecv_double(const double *sendbuf, const int sendcount,
442 : const int dest, const int sendtag, double *recvbuf,
443 : const int recvcount, const int source,
444 : const int recvtag, const cp_mpi_comm_t comm) {
445 : #if defined(__parallel)
446 19134 : MPI_Status status;
447 19134 : CHECK(MPI_Sendrecv(sendbuf, sendcount, MPI_DOUBLE, dest, sendtag, recvbuf,
448 : recvcount, MPI_DOUBLE, source, recvtag, comm, &status));
449 19134 : int count_received;
450 19134 : CHECK(MPI_Get_count(&status, MPI_DOUBLE, &count_received));
451 19134 : return count_received;
452 : #else
453 0 : (void)sendbuf; // mark used
454 0 : (void)sendcount;
455 0 : (void)dest;
456 0 : (void)sendtag;
457 0 : (void)recvbuf;
458 0 : (void)recvcount;
459 0 : (void)source;
460 0 : (void)recvtag;
461 0 : (void)comm;
462 0 : fprintf(stderr, "Error: cp_mpi_sendrecv_double not available without MPI\n");
463 0 : abort();
464 : #endif
465 : }
466 :
467 : /*******************************************************************************
468 : * \brief Wrapper around MPI_Alltoall for datatype MPI_INT.
469 : * \author Ole Schuett
470 : ******************************************************************************/
471 891420 : void cp_mpi_alltoall_int(const int *sendbuf, const int sendcount, int *recvbuf,
472 : const int recvcount, const cp_mpi_comm_t comm) {
473 : #if defined(__parallel)
474 891420 : CHECK(MPI_Alltoall(sendbuf, sendcount, MPI_INT, recvbuf, recvcount, MPI_INT,
475 : comm));
476 : #else
477 0 : (void)comm; // mark used
478 0 : assert(sendcount == recvcount);
479 0 : memcpy(recvbuf, sendbuf, sendcount * sizeof(int));
480 : #endif
481 891420 : }
482 :
483 : /*******************************************************************************
484 : * \brief Wrapper around MPI_Alltoallv for datatype MPI_BYTE.
485 : * \author Ole Schuett
486 : ******************************************************************************/
487 445638 : void cp_mpi_alltoallv_byte(const void *sendbuf, const int *sendcounts,
488 : const int *sdispls, void *recvbuf,
489 : const int *recvcounts, const int *rdispls,
490 : const cp_mpi_comm_t comm) {
491 : #if defined(__parallel)
492 445638 : CHECK(MPI_Alltoallv(sendbuf, sendcounts, sdispls, MPI_BYTE, recvbuf,
493 : recvcounts, rdispls, MPI_BYTE, comm));
494 : #else
495 0 : (void)comm; // mark used
496 0 : assert(sendcounts[0] == recvcounts[0]);
497 0 : assert(sdispls[0] == 0 && rdispls[0] == 0);
498 0 : memcpy(recvbuf, sendbuf, sendcounts[0]);
499 : #endif
500 445638 : }
501 :
502 : /*******************************************************************************
503 : * \brief Wrapper around MPI_Alltoallv for datatype MPI_DOUBLE.
504 : * \author Ole Schuett
505 : ******************************************************************************/
506 445782 : void cp_mpi_alltoallv_double(const double *sendbuf, const int *sendcounts,
507 : const int *sdispls, double *recvbuf,
508 : const int *recvcounts, const int *rdispls,
509 : const cp_mpi_comm_t comm) {
510 : #if defined(__parallel)
511 445782 : CHECK(MPI_Alltoallv(sendbuf, sendcounts, sdispls, MPI_DOUBLE, recvbuf,
512 : recvcounts, rdispls, MPI_DOUBLE, comm));
513 : #else
514 0 : (void)comm; // mark used
515 0 : assert(sendcounts[0] == recvcounts[0]);
516 0 : assert(sdispls[0] == 0 && rdispls[0] == 0);
517 0 : memcpy(recvbuf, sendbuf, sendcounts[0] * sizeof(double));
518 : #endif
519 445782 : }
520 :
521 : /*******************************************************************************
522 : * \brief Wrapper around MPI_Alloc_mem.
523 : * \author Hans Pabst
524 : ******************************************************************************/
525 1936067 : void *cp_mpi_alloc_mem(size_t size) {
526 1936067 : void *result = NULL;
527 : #if DBM_ALLOC_MPI && defined(__parallel)
528 : CHECK(MPI_Alloc_mem((MPI_Aint)size, MPI_INFO_NULL, &result));
529 : #elif DBM_ALLOC_OPENMP && (201811 /*v5.0*/ <= _OPENMP)
530 : result = omp_alloc(size, omp_null_allocator);
531 : #else
532 1936067 : result = malloc(size);
533 : #endif
534 1936067 : return result;
535 : }
536 :
537 : /*******************************************************************************
538 : * \brief Wrapper around MPI_Free_mem.
539 : * \author Hans Pabst
540 : ******************************************************************************/
541 1936067 : void cp_mpi_free_mem(void *mem) {
542 : #if DBM_ALLOC_MPI && defined(__parallel)
543 : CHECK(MPI_Free_mem(mem));
544 : #elif DBM_ALLOC_OPENMP && (201811 /*v5.0*/ <= _OPENMP)
545 : omp_free(mem, omp_null_allocator);
546 : #else
547 1936067 : free(mem);
548 : #endif
549 1936067 : }
550 :
551 : // EOF
|