2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2008,2009,2010,2012,2013,2014,2015,2016,2018, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
38 * \brief This file defines functions for (mostly) the domdec module
39 * to use MPI functionality
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_domdec
47 #include "domdec_network.h"
53 #include "gromacs/domdec/domdec_struct.h"
54 #include "gromacs/utility/gmxmpi.h"
57 /*! \brief Returns the MPI rank of the domain decomposition master rank */
58 #define DDMASTERRANK(dd) ((dd)->masterrank)
61 /*! \brief Move data of type \p T in the communication region one cell along
62 * the domain decomposition
64 * Moves in the dimension indexed by ddDimensionIndex, either forward
65 * (direction=dddirFoward) or backward (direction=dddirBackward).
69 ddSendrecv(const struct gmx_domdec_t
*dd
,
73 int numElementsToSend
,
75 int numElementsToReceive
)
78 int sendRank
= dd
->neighbor
[ddDimensionIndex
][direction
== dddirForward
? 0 : 1];
79 int receiveRank
= dd
->neighbor
[ddDimensionIndex
][direction
== dddirForward
? 1 : 0];
81 constexpr int mpiTag
= 0;
83 if (numElementsToSend
> 0 && numElementsToReceive
> 0)
85 MPI_Sendrecv(sendBuffer
, numElementsToSend
*sizeof(T
), MPI_BYTE
,
87 receiveBuffer
, numElementsToReceive
*sizeof(T
), MPI_BYTE
,
92 else if (numElementsToSend
> 0)
94 MPI_Send( sendBuffer
, numElementsToSend
*sizeof(T
), MPI_BYTE
,
98 else if (numElementsToReceive
> 0)
100 MPI_Recv( receiveBuffer
, numElementsToReceive
*sizeof(T
), MPI_BYTE
,
106 GMX_UNUSED_VALUE(dd
);
107 GMX_UNUSED_VALUE(ddDimensionIndex
);
108 GMX_UNUSED_VALUE(direction
);
109 GMX_UNUSED_VALUE(sendBuffer
);
110 GMX_UNUSED_VALUE(numElementsToSend
);
111 GMX_UNUSED_VALUE(receiveBuffer
);
112 GMX_UNUSED_VALUE(numElementsToReceive
);
116 //! Specialization of extern template for int
117 template void ddSendrecv(const gmx_domdec_t
*, int, int,
118 int *, int, int *, int);
119 //! Specialization of extern template for real
120 template void ddSendrecv(const gmx_domdec_t
*, int, int,
121 real
*, int, real
*, int);
122 //! Specialization of extern template for gmx::RVec
123 template void ddSendrecv(const gmx_domdec_t
*, int, int,
124 rvec
*, int, rvec
*, int);
126 template <typename T
>
128 ddSendrecv(const gmx_domdec_t
*dd
,
129 int ddDimensionIndex
,
131 gmx::ArrayRef
<T
> sendBuffer
,
132 gmx::ArrayRef
<T
> receiveBuffer
)
134 ddSendrecv(dd
, ddDimensionIndex
, direction
,
135 sendBuffer
.data(), sendBuffer
.size(),
136 receiveBuffer
.data(), receiveBuffer
.size());
139 //! Specialization of extern template for int
140 template void ddSendrecv(const gmx_domdec_t
*, int, int,
141 gmx::ArrayRef
<int>, gmx::ArrayRef
<int>);
142 //! Specialization of extern template for real
143 template void ddSendrecv(const gmx_domdec_t
*, int, int,
144 gmx::ArrayRef
<real
>, gmx::ArrayRef
<real
>);
145 //! Specialization of extern template for gmx::RVec
146 template void ddSendrecv(const gmx_domdec_t
*, int, int,
147 gmx::ArrayRef
<gmx::RVec
>, gmx::ArrayRef
<gmx::RVec
>);
149 void dd_sendrecv2_rvec(const struct gmx_domdec_t gmx_unused
*dd
,
150 int gmx_unused ddimind
,
151 rvec gmx_unused
*buf_s_fw
, int gmx_unused n_s_fw
,
152 rvec gmx_unused
*buf_r_fw
, int gmx_unused n_r_fw
,
153 rvec gmx_unused
*buf_s_bw
, int gmx_unused n_s_bw
,
154 rvec gmx_unused
*buf_r_bw
, int gmx_unused n_r_bw
)
157 int rank_fw
, rank_bw
, nreq
;
161 rank_fw
= dd
->neighbor
[ddimind
][0];
162 rank_bw
= dd
->neighbor
[ddimind
][1];
166 /* Try to send and receive in two directions simultaneously.
167 * Should be faster, especially on machines
168 * with full 3D communication networks.
169 * However, it could be that communication libraries are
170 * optimized for MPI_Sendrecv and non-blocking MPI calls
172 * SendRecv2 can be turned on with the env.var. GMX_DD_SENDRECV2
177 MPI_Irecv(buf_r_fw
[0], n_r_fw
*sizeof(rvec
), MPI_BYTE
,
178 rank_bw
, 0, dd
->mpi_comm_all
, &req
[nreq
++]);
182 MPI_Irecv(buf_r_bw
[0], n_r_bw
*sizeof(rvec
), MPI_BYTE
,
183 rank_fw
, 1, dd
->mpi_comm_all
, &req
[nreq
++]);
187 MPI_Isend(buf_s_fw
[0], n_s_fw
*sizeof(rvec
), MPI_BYTE
,
188 rank_fw
, 0, dd
->mpi_comm_all
, &req
[nreq
++]);
192 MPI_Isend(buf_s_bw
[0], n_s_bw
*sizeof(rvec
), MPI_BYTE
,
193 rank_bw
, 1, dd
->mpi_comm_all
, &req
[nreq
++]);
197 MPI_Waitall(nreq
, req
, stat
); //NOLINT(clang-analyzer-optin.mpi.MPI-Checker)
202 /* Communicate in two ordered phases.
203 * This is slower, even on a dual-core Opteron cluster
204 * with a single full-duplex network connection per machine.
207 MPI_Sendrecv(buf_s_fw
[0], n_s_fw
*sizeof(rvec
), MPI_BYTE
, rank_fw
, 0,
208 buf_r_fw
[0], n_r_fw
*sizeof(rvec
), MPI_BYTE
, rank_bw
, 0,
209 dd
->mpi_comm_all
, &stat
[0]);
211 MPI_Sendrecv(buf_s_bw
[0], n_s_bw
*sizeof(rvec
), MPI_BYTE
, rank_bw
, 0,
212 buf_r_bw
[0], n_r_bw
*sizeof(rvec
), MPI_BYTE
, rank_fw
, 0,
213 dd
->mpi_comm_all
, &stat
[0]);
218 void dd_bcast(const gmx_domdec_t gmx_unused
*dd
, int gmx_unused nbytes
, void gmx_unused
*data
)
223 MPI_Bcast(data
, nbytes
, MPI_BYTE
,
224 DDMASTERRANK(dd
), dd
->mpi_comm_all
);
229 void dd_bcastc(const gmx_domdec_t
*dd
, int nbytes
, void *src
, void *dest
)
231 if (DDMASTER(dd
) || dd
->nnodes
== 1)
233 memcpy(dest
, src
, nbytes
);
238 MPI_Bcast(dest
, nbytes
, MPI_BYTE
,
239 DDMASTERRANK(dd
), dd
->mpi_comm_all
);
244 void dd_scatter(const gmx_domdec_t gmx_unused
*dd
, int gmx_unused nbytes
, const void gmx_unused
*src
, void *dest
)
249 /* Some MPI implementions don't specify const */
250 MPI_Scatter(const_cast<void *>(src
), nbytes
, MPI_BYTE
,
251 dest
, nbytes
, MPI_BYTE
,
252 DDMASTERRANK(dd
), dd
->mpi_comm_all
);
257 /* 1 rank, either we copy everything, or dest=src: nothing to do */
260 memcpy(dest
, src
, nbytes
);
265 void dd_gather(const gmx_domdec_t gmx_unused
*dd
, int gmx_unused nbytes
, const void gmx_unused
*src
, void gmx_unused
*dest
)
268 /* Some MPI implementions don't specify const */
269 MPI_Gather(const_cast<void *>(src
), nbytes
, MPI_BYTE
,
270 dest
, nbytes
, MPI_BYTE
,
271 DDMASTERRANK(dd
), dd
->mpi_comm_all
);
275 void dd_scatterv(const gmx_domdec_t gmx_unused
*dd
,
276 int gmx_unused
*scounts
, int gmx_unused
*disps
, const void *sbuf
,
277 int rcount
, void *rbuf
)
286 /* MPI does not allow NULL pointers */
289 /* Some MPI implementions don't specify const */
290 MPI_Scatterv(const_cast<void *>(sbuf
), scounts
, disps
, MPI_BYTE
,
291 rbuf
, rcount
, MPI_BYTE
,
292 DDMASTERRANK(dd
), dd
->mpi_comm_all
);
297 /* 1 rank, either we copy everything, or rbuf=sbuf: nothing to do */
300 memcpy(rbuf
, sbuf
, rcount
);
305 void dd_gatherv(const gmx_domdec_t gmx_unused
*dd
,
306 int gmx_unused scount
, const void gmx_unused
*sbuf
,
307 int gmx_unused
*rcounts
, int gmx_unused
*disps
, void gmx_unused
*rbuf
)
314 /* MPI does not allow NULL pointers */
317 /* Some MPI implementions don't specify const */
318 MPI_Gatherv(const_cast<void *>(sbuf
), scount
, MPI_BYTE
,
319 rbuf
, rcounts
, disps
, MPI_BYTE
,
320 DDMASTERRANK(dd
), dd
->mpi_comm_all
);