(git:e966546)
Loading...
Searching...
No Matches
cp_mpi.h
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: GPL-2.0-or-later */
6/*----------------------------------------------------------------------------*/
7#ifndef CP_MPI_H
8#define CP_MPI_H
9
10#include <stdbool.h>
11#include <stddef.h>
12#include <stdint.h>
13
14#if defined(__parallel)
15#include <mpi.h>
16typedef MPI_Comm cp_mpi_comm_t;
17#else
18typedef int cp_mpi_comm_t;
19#endif
20
21/*******************************************************************************
22 * \brief Wrapper around MPI_Init.
23 * \author Ole Schuett
24 ******************************************************************************/
25void cp_mpi_init(int *argc, char ***argv);
26
27/*******************************************************************************
28 * \brief Wrapper around MPI_Finalize.
29 * \author Ole Schuett
30 ******************************************************************************/
31void cp_mpi_finalize(void);
32
33/*******************************************************************************
34 * \brief Returns MPI_COMM_WORLD.
35 * \author Ole Schuett
36 ******************************************************************************/
38
39/*******************************************************************************
40 * \brief Wrapper around MPI_Comm_f2c.
41 * \author Ole Schuett
42 ******************************************************************************/
43cp_mpi_comm_t cp_mpi_comm_f2c(const int fortran_comm);
44
45/*******************************************************************************
46 * \brief Wrapper around MPI_Comm_c2f.
47 * \author Ole Schuett
48 ******************************************************************************/
49int cp_mpi_comm_c2f(const cp_mpi_comm_t comm);
50
51/*******************************************************************************
52 * \brief Wrapper around MPI_Comm_rank.
53 * \author Ole Schuett
54 ******************************************************************************/
55int cp_mpi_comm_rank(const cp_mpi_comm_t comm);
56
57/*******************************************************************************
58 * \brief Wrapper around MPI_Comm_size.
59 * \author Ole Schuett
60 ******************************************************************************/
61int cp_mpi_comm_size(const cp_mpi_comm_t comm);
62
63/*******************************************************************************
64 * \brief Wrapper around MPI_Dims_create.
65 * \author Ole Schuett
66 ******************************************************************************/
67void cp_mpi_dims_create(const int nnodes, const int ndims, int dims[]);
68
69/*******************************************************************************
70 * \brief Wrapper around MPI_Cart_create.
71 * \author Ole Schuett
72 ******************************************************************************/
73cp_mpi_comm_t cp_mpi_cart_create(const cp_mpi_comm_t comm_old, const int ndims,
74 const int dims[], const int periods[],
75 const int reorder);
76
77/*******************************************************************************
78 * \brief Wrapper around MPI_Cart_get.
79 * \author Ole Schuett
80 ******************************************************************************/
81void cp_mpi_cart_get(const cp_mpi_comm_t comm, int maxdims, int dims[],
82 int periods[], int coords[]);
83
84/*******************************************************************************
85 * \brief Wrapper around MPI_Cart_rank.
86 * \author Ole Schuett
87 ******************************************************************************/
88int cp_mpi_cart_rank(const cp_mpi_comm_t comm, const int coords[]);
89
90/*******************************************************************************
91 * \brief Wrapper around MPI_Cart_sub.
92 * \author Ole Schuett
93 ******************************************************************************/
95 const int remain_dims[]);
96
97/*******************************************************************************
98 * \brief Wrapper around MPI_Comm_free.
99 * \author Ole Schuett
100 ******************************************************************************/
102
103/*******************************************************************************
104 * \brief Wrapper around MPI_Comm_compare.
105 * \author Ole Schuett
106 ******************************************************************************/
108 const cp_mpi_comm_t comm2);
109
110/*******************************************************************************
111 * \brief Wrapper around MPI_Allreduce for op MPI_MAX and datatype MPI_INT.
112 * \author Ole Schuett
113 ******************************************************************************/
114void cp_mpi_max_int(int *values, const int count, const cp_mpi_comm_t comm);
115
116/*******************************************************************************
117 * \brief Wrapper around MPI_Allreduce for op MPI_MAX and datatype MPI_UINT64_T.
118 * \author Hans Pabst
119 ******************************************************************************/
120void cp_mpi_max_uint64(uint64_t *values, const int count,
121 const cp_mpi_comm_t comm);
122
123/*******************************************************************************
124 * \brief Wrapper around MPI_Allreduce for op MPI_MAX and datatype MPI_DOUBLE.
125 * \author Ole Schuett
126 ******************************************************************************/
127void cp_mpi_max_double(double *values, const int count,
128 const cp_mpi_comm_t comm);
129
130/*******************************************************************************
131 * \brief Wrapper around MPI_Allreduce for op MPI_SUM and datatype MPI_INT.
132 * \author Ole Schuett
133 ******************************************************************************/
134void cp_mpi_sum_int(int *values, const int count, const cp_mpi_comm_t comm);
135
136/*******************************************************************************
137 * \brief Wrapper around MPI_Allreduce for op MPI_SUM and datatype MPI_INT64_T.
138 * \author Ole Schuett
139 ******************************************************************************/
140void cp_mpi_sum_int64(int64_t *values, const int count,
141 const cp_mpi_comm_t comm);
142
143/*******************************************************************************
144 * \brief Wrapper around MPI_Allreduce for op MPI_SUM and datatype MPI_DOUBLE.
145 * \author Ole Schuett
146 ******************************************************************************/
147void cp_mpi_sum_double(double *values, const int count,
148 const cp_mpi_comm_t comm);
149
150/*******************************************************************************
151 * \brief Wrapper around MPI_Sendrecv for datatype MPI_BYTE.
152 * \author Ole Schuett
153 ******************************************************************************/
154int cp_mpi_sendrecv_byte(const void *sendbuf, const int sendcount,
155 const int dest, const int sendtag, void *recvbuf,
156 const int recvcount, const int source,
157 const int recvtag, const cp_mpi_comm_t comm);
158
159/*******************************************************************************
160 * \brief Wrapper around MPI_Sendrecv for datatype MPI_DOUBLE.
161 * \author Ole Schuett
162 ******************************************************************************/
163int cp_mpi_sendrecv_double(const double *sendbuf, const int sendcount,
164 const int dest, const int sendtag, double *recvbuf,
165 const int recvcount, const int source,
166 const int recvtag, const cp_mpi_comm_t comm);
167
168/*******************************************************************************
169 * \brief Wrapper around MPI_Alltoall for datatype MPI_INT.
170 * \author Ole Schuett
171 ******************************************************************************/
172void cp_mpi_alltoall_int(const int *sendbuf, const int sendcount, int *recvbuf,
173 const int recvcount, const cp_mpi_comm_t comm);
174
175/*******************************************************************************
176 * \brief Wrapper around MPI_Alltoallv for datatype MPI_BYTE.
177 * \author Ole Schuett
178 ******************************************************************************/
179void cp_mpi_alltoallv_byte(const void *sendbuf, const int *sendcounts,
180 const int *sdispls, void *recvbuf,
181 const int *recvcounts, const int *rdispls,
182 const cp_mpi_comm_t comm);
183
184/*******************************************************************************
185 * \brief Wrapper around MPI_Alltoallv for datatype MPI_DOUBLE.
186 * \author Ole Schuett
187 ******************************************************************************/
188void cp_mpi_alltoallv_double(const double *sendbuf, const int *sendcounts,
189 const int *sdispls, double *recvbuf,
190 const int *recvcounts, const int *rdispls,
191 const cp_mpi_comm_t comm);
192
193/*******************************************************************************
194 * \brief Wrapper around MPI_Alloc_mem.
195 * \author Hans Pabst
196 ******************************************************************************/
197void *cp_mpi_alloc_mem(size_t size);
198
199/*******************************************************************************
200 * \brief Wrapper around MPI_Free_mem.
201 * \author Hans Pabst
202 ******************************************************************************/
203void cp_mpi_free_mem(void *ptr);
204
205#endif
206
207// EOF
int cp_mpi_comm_size(const cp_mpi_comm_t comm)
Wrapper around MPI_Comm_size.
Definition cp_mpi.c:110
void cp_mpi_max_int(int *values, const int count, const cp_mpi_comm_t comm)
Wrapper around MPI_Allreduce for op MPI_MAX and datatype MPI_INT.
Definition cp_mpi.c:241
void cp_mpi_sum_double(double *values, const int count, const cp_mpi_comm_t comm)
Wrapper around MPI_Allreduce for op MPI_SUM and datatype MPI_DOUBLE.
Definition cp_mpi.c:347
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
int cp_mpi_sendrecv_byte(const void *sendbuf, const int sendcount, const int dest, const int sendtag, void *recvbuf, const int recvcount, const int source, const int recvtag, const cp_mpi_comm_t comm)
Wrapper around MPI_Sendrecv for datatype MPI_BYTE.
Definition cp_mpi.c:369
void cp_mpi_alltoallv_double(const double *sendbuf, const int *sendcounts, const int *sdispls, double *recvbuf, const int *recvcounts, const int *rdispls, const cp_mpi_comm_t comm)
Wrapper around MPI_Alltoallv for datatype MPI_DOUBLE.
Definition cp_mpi.c:464
int cp_mpi_cart_rank(const cp_mpi_comm_t comm, const int coords[])
Wrapper around MPI_Cart_rank.
Definition cp_mpi.c:179
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
void * cp_mpi_alloc_mem(size_t size)
Wrapper around MPI_Alloc_mem.
Definition cp_mpi.c:483
void cp_mpi_max_uint64(uint64_t *values, const int count, const cp_mpi_comm_t comm)
Wrapper around MPI_Allreduce for op MPI_MAX and datatype MPI_UINT64_T.
Definition cp_mpi.c:261
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_max_double(double *values, const int count, const cp_mpi_comm_t comm)
Wrapper around MPI_Allreduce for op MPI_MAX and datatype MPI_DOUBLE.
Definition cp_mpi.c:283
int cp_mpi_comm_t
Definition cp_mpi.h:18
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_alltoall_int(const int *sendbuf, const int sendcount, int *recvbuf, const int recvcount, const cp_mpi_comm_t comm)
Wrapper around MPI_Alltoall for datatype MPI_INT.
Definition cp_mpi.c:429
void cp_mpi_dims_create(const int nnodes, const int ndims, int dims[])
Wrapper around MPI_Dims_create.
Definition cp_mpi.c:125
bool cp_mpi_comms_are_similar(const cp_mpi_comm_t comm1, const cp_mpi_comm_t comm2)
Wrapper around MPI_Comm_compare.
Definition cp_mpi.c:224
void cp_mpi_sum_int(int *values, const int count, const cp_mpi_comm_t comm)
Wrapper around MPI_Allreduce for op MPI_SUM and datatype MPI_INT.
Definition cp_mpi.c:305
void cp_mpi_free_mem(void *ptr)
Wrapper around MPI_Free_mem.
Definition cp_mpi.c:499
void cp_mpi_alltoallv_byte(const void *sendbuf, const int *sendcounts, const int *sdispls, void *recvbuf, const int *recvcounts, const int *rdispls, const cp_mpi_comm_t comm)
Wrapper around MPI_Alltoallv for datatype MPI_BYTE.
Definition cp_mpi.c:445
void cp_mpi_finalize(void)
Wrapper around MPI_Finalize.
Definition cp_mpi.c:47
cp_mpi_comm_t cp_mpi_comm_f2c(const int fortran_comm)
Wrapper around MPI_Comm_f2c.
Definition cp_mpi.c:69
cp_mpi_comm_t cp_mpi_cart_sub(const cp_mpi_comm_t comm, const int remain_dims[])
Wrapper around MPI_Cart_sub.
Definition cp_mpi.c:195
int cp_mpi_sendrecv_double(const double *sendbuf, const int sendcount, const int dest, const int sendtag, double *recvbuf, const int recvcount, const int source, const int recvtag, const cp_mpi_comm_t comm)
Wrapper around MPI_Sendrecv for datatype MPI_DOUBLE.
Definition cp_mpi.c:399
int cp_mpi_comm_rank(const cp_mpi_comm_t comm)
Wrapper around MPI_Comm_rank.
Definition cp_mpi.c:95