Collect settings for DD in DDSettings
[alexxy/gromacs.git] / src / gromacs / domdec / domdec_network.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2008,2009,2010,2012,2013,2014,2015,2016,2018,2019, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 /*! \internal \file
37  *
38  * \brief This file defines functions for (mostly) the domdec module
39  * to use MPI functionality
40  *
41  * \author Berk Hess <hess@kth.se>
42  * \ingroup module_domdec
43  */
44
45 #include "gmxpre.h"
46
47 #include "domdec_network.h"
48
49 #include "config.h"
50
51 #include <cstring>
52
53 #include "gromacs/domdec/domdec_struct.h"
54 #include "gromacs/utility/gmxmpi.h"
55
56 #include "domdec_internal.h"
57
58
59 /*! \brief Returns the MPI rank of the domain decomposition master rank */
60 #define DDMASTERRANK(dd)   ((dd)->masterrank)
61
62
63 /*! \brief Move data of type \p T in the communication region one cell along
64  * the domain decomposition
65  *
66  * Moves in the dimension indexed by ddDimensionIndex, either forward
67  * (direction=dddirFoward) or backward (direction=dddirBackward).
68  */
69 template <typename T>
70 void
71 ddSendrecv(const struct gmx_domdec_t *dd,
72            int                        ddDimensionIndex,
73            int                        direction,
74            T                         *sendBuffer,
75            int                        numElementsToSend,
76            T                         *receiveBuffer,
77            int                        numElementsToReceive)
78 {
79 #if GMX_MPI
80     int           sendRank    = dd->neighbor[ddDimensionIndex][direction == dddirForward ? 0 : 1];
81     int           receiveRank = dd->neighbor[ddDimensionIndex][direction == dddirForward ? 1 : 0];
82
83     constexpr int mpiTag      = 0;
84     MPI_Status    mpiStatus;
85     if (numElementsToSend > 0 && numElementsToReceive > 0)
86     {
87         MPI_Sendrecv(sendBuffer,    numElementsToSend*sizeof(T),    MPI_BYTE,
88                      sendRank,      mpiTag,
89                      receiveBuffer, numElementsToReceive*sizeof(T), MPI_BYTE,
90                      receiveRank,   mpiTag,
91                      dd->mpi_comm_all,
92                      &mpiStatus);
93     }
94     else if (numElementsToSend > 0)
95     {
96         MPI_Send(    sendBuffer,    numElementsToSend*sizeof(T),    MPI_BYTE,
97                      sendRank,      mpiTag,
98                      dd->mpi_comm_all);
99     }
100     else if (numElementsToReceive > 0)
101     {
102         MPI_Recv(    receiveBuffer, numElementsToReceive*sizeof(T), MPI_BYTE,
103                      receiveRank,   mpiTag,
104                      dd->mpi_comm_all,
105                      &mpiStatus);
106     }
107 #else // GMX_MPI
108     GMX_UNUSED_VALUE(dd);
109     GMX_UNUSED_VALUE(ddDimensionIndex);
110     GMX_UNUSED_VALUE(direction);
111     GMX_UNUSED_VALUE(sendBuffer);
112     GMX_UNUSED_VALUE(numElementsToSend);
113     GMX_UNUSED_VALUE(receiveBuffer);
114     GMX_UNUSED_VALUE(numElementsToReceive);
115 #endif // GMX_MPI
116 }
117
118 //! Specialization of extern template for int
119 template void ddSendrecv(const gmx_domdec_t *, int, int,
120                          int *, int, int *, int);
121 //! Specialization of extern template for real
122 template void ddSendrecv(const gmx_domdec_t *, int, int,
123                          real *, int, real *, int);
124 //! Specialization of extern template for gmx::RVec
125 template void ddSendrecv(const gmx_domdec_t *, int, int,
126                          rvec *, int, rvec *, int);
127
128 template <typename T>
129 void
130 ddSendrecv(const gmx_domdec_t *dd,
131            int                 ddDimensionIndex,
132            int                 direction,
133            gmx::ArrayRef<T>    sendBuffer,
134            gmx::ArrayRef<T>    receiveBuffer)
135 {
136     ddSendrecv(dd, ddDimensionIndex, direction,
137                sendBuffer.data(), sendBuffer.size(),
138                receiveBuffer.data(), receiveBuffer.size());
139 }
140
141 //! Specialization of extern template for int
142 template void ddSendrecv(const gmx_domdec_t *, int, int,
143                          gmx::ArrayRef<int>, gmx::ArrayRef<int>);
144 //! Specialization of extern template for real
145 template void ddSendrecv(const gmx_domdec_t *, int, int,
146                          gmx::ArrayRef<real>, gmx::ArrayRef<real>);
147 //! Specialization of extern template for gmx::RVec
148 template void ddSendrecv(const gmx_domdec_t *, int, int,
149                          gmx::ArrayRef<gmx::RVec>, gmx::ArrayRef<gmx::RVec>);
150
151 void dd_sendrecv2_rvec(const struct gmx_domdec_t gmx_unused *dd,
152                        int gmx_unused ddimind,
153                        rvec gmx_unused *buf_s_fw, int gmx_unused n_s_fw,
154                        rvec gmx_unused *buf_r_fw, int gmx_unused n_r_fw,
155                        rvec gmx_unused *buf_s_bw, int gmx_unused n_s_bw,
156                        rvec gmx_unused *buf_r_bw, int gmx_unused n_r_bw)
157 {
158 #if GMX_MPI
159     int         rank_fw, rank_bw, nreq;
160     MPI_Request req[4];
161     MPI_Status  stat[4];
162
163     rank_fw = dd->neighbor[ddimind][0];
164     rank_bw = dd->neighbor[ddimind][1];
165
166     if (!dd->comm->ddSettings.useSendRecv2)
167     {
168         /* Try to send and receive in two directions simultaneously.
169          * Should be faster, especially on machines
170          * with full 3D communication networks.
171          * However, it could be that communication libraries are
172          * optimized for MPI_Sendrecv and non-blocking MPI calls
173          * are slower.
174          * SendRecv2 can be turned on with the env.var. GMX_DD_SENDRECV2
175          */
176         nreq = 0;
177         if (n_r_fw)
178         {
179             MPI_Irecv(buf_r_fw[0], n_r_fw*sizeof(rvec), MPI_BYTE,
180                       rank_bw, 0, dd->mpi_comm_all, &req[nreq++]);
181         }
182         if (n_r_bw)
183         {
184             MPI_Irecv(buf_r_bw[0], n_r_bw*sizeof(rvec), MPI_BYTE,
185                       rank_fw, 1, dd->mpi_comm_all, &req[nreq++]);
186         }
187         if (n_s_fw)
188         {
189             MPI_Isend(buf_s_fw[0], n_s_fw*sizeof(rvec), MPI_BYTE,
190                       rank_fw, 0, dd->mpi_comm_all, &req[nreq++]);
191         }
192         if (n_s_bw)
193         {
194             MPI_Isend(buf_s_bw[0], n_s_bw*sizeof(rvec), MPI_BYTE,
195                       rank_bw, 1, dd->mpi_comm_all, &req[nreq++]);
196         }
197         if (nreq)
198         {
199             MPI_Waitall(nreq, req, stat); //NOLINT(clang-analyzer-optin.mpi.MPI-Checker)
200         }
201     }
202     else
203     {
204         /* Communicate in two ordered phases.
205          * This is slower, even on a dual-core Opteron cluster
206          * with a single full-duplex network connection per machine.
207          */
208         /* Forward */
209         MPI_Sendrecv(buf_s_fw[0], n_s_fw*sizeof(rvec), MPI_BYTE, rank_fw, 0,
210                      buf_r_fw[0], n_r_fw*sizeof(rvec), MPI_BYTE, rank_bw, 0,
211                      dd->mpi_comm_all, &stat[0]);
212         /* Backward */
213         MPI_Sendrecv(buf_s_bw[0], n_s_bw*sizeof(rvec), MPI_BYTE, rank_bw, 0,
214                      buf_r_bw[0], n_r_bw*sizeof(rvec), MPI_BYTE, rank_fw, 0,
215                      dd->mpi_comm_all, &stat[0]);
216     }
217 #endif
218 }
219
220 void dd_bcast(const gmx_domdec_t gmx_unused *dd, int gmx_unused nbytes, void gmx_unused *data)
221 {
222 #if GMX_MPI
223     if (dd->nnodes > 1)
224     {
225         MPI_Bcast(data, nbytes, MPI_BYTE,
226                   DDMASTERRANK(dd), dd->mpi_comm_all);
227     }
228 #endif
229 }
230
231 void dd_bcastc(const gmx_domdec_t *dd, int nbytes, void *src, void *dest)
232 {
233     if (DDMASTER(dd) || dd->nnodes == 1)
234     {
235         memcpy(dest, src, nbytes);
236     }
237 #if GMX_MPI
238     if (dd->nnodes > 1)
239     {
240         MPI_Bcast(dest, nbytes, MPI_BYTE,
241                   DDMASTERRANK(dd), dd->mpi_comm_all);
242     }
243 #endif
244 }
245
246 void dd_scatter(const gmx_domdec_t gmx_unused *dd, int gmx_unused nbytes, const void gmx_unused *src, void *dest)
247 {
248 #if GMX_MPI
249     if (dd->nnodes > 1)
250     {
251         /* Some MPI implementions don't specify const */
252         MPI_Scatter(const_cast<void *>(src), nbytes, MPI_BYTE,
253                     dest, nbytes, MPI_BYTE,
254                     DDMASTERRANK(dd), dd->mpi_comm_all);
255     }
256     else
257 #endif
258     {
259         /* 1 rank, either we copy everything, or dest=src: nothing to do */
260         if (dest != src)
261         {
262             memcpy(dest, src, nbytes);
263         }
264     }
265 }
266
267 void dd_gather(const gmx_domdec_t gmx_unused *dd, int gmx_unused nbytes, const void gmx_unused *src, void gmx_unused *dest)
268 {
269 #if GMX_MPI
270     /* Some MPI implementions don't specify const */
271     MPI_Gather(const_cast<void *>(src), nbytes, MPI_BYTE,
272                dest, nbytes, MPI_BYTE,
273                DDMASTERRANK(dd), dd->mpi_comm_all);
274 #endif
275 }
276
277 void dd_scatterv(const gmx_domdec_t gmx_unused *dd,
278                  int gmx_unused *scounts, int gmx_unused *disps, const void *sbuf,
279                  int rcount, void *rbuf)
280 {
281 #if GMX_MPI
282     int dum;
283
284     if (dd->nnodes > 1)
285     {
286         if (rcount == 0)
287         {
288             /* MPI does not allow NULL pointers */
289             rbuf = &dum;
290         }
291         /* Some MPI implementions don't specify const */
292         MPI_Scatterv(const_cast<void *>(sbuf), scounts, disps, MPI_BYTE,
293                      rbuf, rcount, MPI_BYTE,
294                      DDMASTERRANK(dd), dd->mpi_comm_all);
295     }
296     else
297 #endif
298     {
299         /* 1 rank, either we copy everything, or rbuf=sbuf: nothing to do */
300         if (rbuf != sbuf)
301         {
302             memcpy(rbuf, sbuf, rcount);
303         }
304     }
305 }
306
307 void dd_gatherv(const gmx_domdec_t gmx_unused *dd,
308                 int gmx_unused scount, const void gmx_unused *sbuf,
309                 int gmx_unused *rcounts, int gmx_unused *disps, void gmx_unused *rbuf)
310 {
311 #if GMX_MPI
312     int dum;
313
314     if (scount == 0)
315     {
316         /* MPI does not allow NULL pointers */
317         sbuf = &dum;
318     }
319     /* Some MPI implementions don't specify const */
320     MPI_Gatherv(const_cast<void *>(sbuf), scount, MPI_BYTE,
321                 rbuf, rcounts, disps, MPI_BYTE,
322                 DDMASTERRANK(dd), dd->mpi_comm_all);
323 #endif
324 }