Convert domdec code to C++
[alexxy/gromacs.git] / src / gromacs / mdlib / 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, 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 #include "gmxpre.h"
37
38 #include "gromacs/legacyheaders/domdec_network.h"
39
40 #include "config.h"
41
42 #include <string.h>
43
44 #include "gromacs/legacyheaders/types/commrec.h"
45 #include "gromacs/utility/gmxmpi.h"
46
47
48 #define DDMASTERRANK(dd)   (dd->masterrank)
49
50
51 void dd_sendrecv_int(const gmx_domdec_t gmx_unused *dd,
52                      int gmx_unused ddimind, int gmx_unused direction,
53                      int gmx_unused *buf_s, int gmx_unused n_s,
54                      int gmx_unused *buf_r, int gmx_unused n_r)
55 {
56 #ifdef GMX_MPI
57     int        rank_s, rank_r;
58     MPI_Status stat;
59
60     rank_s = dd->neighbor[ddimind][direction == dddirForward ? 0 : 1];
61     rank_r = dd->neighbor[ddimind][direction == dddirForward ? 1 : 0];
62
63     if (n_s && n_r)
64     {
65         MPI_Sendrecv(buf_s, n_s*sizeof(int), MPI_BYTE, rank_s, 0,
66                      buf_r, n_r*sizeof(int), MPI_BYTE, rank_r, 0,
67                      dd->mpi_comm_all, &stat);
68     }
69     else if (n_s)
70     {
71         MPI_Send(    buf_s, n_s*sizeof(int), MPI_BYTE, rank_s, 0,
72                      dd->mpi_comm_all);
73     }
74     else if (n_r)
75     {
76         MPI_Recv(    buf_r, n_r*sizeof(int), MPI_BYTE, rank_r, 0,
77                      dd->mpi_comm_all, &stat);
78     }
79
80 #endif
81 }
82
83 void dd_sendrecv_real(const gmx_domdec_t gmx_unused *dd,
84                       int gmx_unused ddimind, int gmx_unused direction,
85                       real gmx_unused *buf_s, int gmx_unused n_s,
86                       real gmx_unused *buf_r, int gmx_unused n_r)
87 {
88 #ifdef GMX_MPI
89     int        rank_s, rank_r;
90     MPI_Status stat;
91
92     rank_s = dd->neighbor[ddimind][direction == dddirForward ? 0 : 1];
93     rank_r = dd->neighbor[ddimind][direction == dddirForward ? 1 : 0];
94
95     if (n_s && n_r)
96     {
97         MPI_Sendrecv(buf_s, n_s*sizeof(real), MPI_BYTE, rank_s, 0,
98                      buf_r, n_r*sizeof(real), MPI_BYTE, rank_r, 0,
99                      dd->mpi_comm_all, &stat);
100     }
101     else if (n_s)
102     {
103         MPI_Send(    buf_s, n_s*sizeof(real), MPI_BYTE, rank_s, 0,
104                      dd->mpi_comm_all);
105     }
106     else if (n_r)
107     {
108         MPI_Recv(    buf_r, n_r*sizeof(real), MPI_BYTE, rank_r, 0,
109                      dd->mpi_comm_all, &stat);
110     }
111
112 #endif
113 }
114
115 void dd_sendrecv_rvec(const gmx_domdec_t gmx_unused *dd,
116                       int gmx_unused ddimind, int gmx_unused direction,
117                       rvec gmx_unused *buf_s, int gmx_unused n_s,
118                       rvec gmx_unused *buf_r, int gmx_unused n_r)
119 {
120 #ifdef GMX_MPI
121     int        rank_s, rank_r;
122     MPI_Status stat;
123
124     rank_s = dd->neighbor[ddimind][direction == dddirForward ? 0 : 1];
125     rank_r = dd->neighbor[ddimind][direction == dddirForward ? 1 : 0];
126
127     if (n_s && n_r)
128     {
129         MPI_Sendrecv(buf_s[0], n_s*sizeof(rvec), MPI_BYTE, rank_s, 0,
130                      buf_r[0], n_r*sizeof(rvec), MPI_BYTE, rank_r, 0,
131                      dd->mpi_comm_all, &stat);
132     }
133     else if (n_s)
134     {
135         MPI_Send(    buf_s[0], n_s*sizeof(rvec), MPI_BYTE, rank_s, 0,
136                      dd->mpi_comm_all);
137     }
138     else if (n_r)
139     {
140         MPI_Recv(    buf_r[0], n_r*sizeof(rvec), MPI_BYTE, rank_r, 0,
141                      dd->mpi_comm_all, &stat);
142     }
143
144 #endif
145 }
146
147 void dd_sendrecv2_rvec(const gmx_domdec_t gmx_unused *dd,
148                        int gmx_unused ddimind,
149                        rvec gmx_unused *buf_s_fw, int gmx_unused n_s_fw,
150                        rvec gmx_unused *buf_r_fw, int gmx_unused n_r_fw,
151                        rvec gmx_unused *buf_s_bw, int gmx_unused n_s_bw,
152                        rvec gmx_unused *buf_r_bw, int gmx_unused n_r_bw)
153 {
154 #ifdef GMX_MPI
155     int         rank_fw, rank_bw, nreq;
156     MPI_Request req[4];
157     MPI_Status  stat[4];
158
159     rank_fw = dd->neighbor[ddimind][0];
160     rank_bw = dd->neighbor[ddimind][1];
161
162     if (!dd->bSendRecv2)
163     {
164         /* Try to send and receive in two directions simultaneously.
165          * Should be faster, especially on machines
166          * with full 3D communication networks.
167          * However, it could be that communication libraries are
168          * optimized for MPI_Sendrecv and non-blocking MPI calls
169          * are slower.
170          * SendRecv2 can be turned on with the env.var. GMX_DD_SENDRECV2
171          */
172         nreq = 0;
173         if (n_r_fw)
174         {
175             MPI_Irecv(buf_r_fw[0], n_r_fw*sizeof(rvec), MPI_BYTE,
176                       rank_bw, 0, dd->mpi_comm_all, &req[nreq++]);
177         }
178         if (n_r_bw)
179         {
180             MPI_Irecv(buf_r_bw[0], n_r_bw*sizeof(rvec), MPI_BYTE,
181                       rank_fw, 1, dd->mpi_comm_all, &req[nreq++]);
182         }
183         if (n_s_fw)
184         {
185             MPI_Isend(buf_s_fw[0], n_s_fw*sizeof(rvec), MPI_BYTE,
186                       rank_fw, 0, dd->mpi_comm_all, &req[nreq++]);
187         }
188         if (n_s_bw)
189         {
190             MPI_Isend(buf_s_bw[0], n_s_bw*sizeof(rvec), MPI_BYTE,
191                       rank_bw, 1, dd->mpi_comm_all, &req[nreq++]);
192         }
193         if (nreq)
194         {
195             MPI_Waitall(nreq, req, stat);
196         }
197     }
198     else
199     {
200         /* Communicate in two ordered phases.
201          * This is slower, even on a dual-core Opteron cluster
202          * with a single full-duplex network connection per machine.
203          */
204         /* Forward */
205         MPI_Sendrecv(buf_s_fw[0], n_s_fw*sizeof(rvec), MPI_BYTE, rank_fw, 0,
206                      buf_r_fw[0], n_r_fw*sizeof(rvec), MPI_BYTE, rank_bw, 0,
207                      dd->mpi_comm_all, &stat[0]);
208         /* Backward */
209         MPI_Sendrecv(buf_s_bw[0], n_s_bw*sizeof(rvec), MPI_BYTE, rank_bw, 0,
210                      buf_r_bw[0], n_r_bw*sizeof(rvec), MPI_BYTE, rank_fw, 0,
211                      dd->mpi_comm_all, &stat[0]);
212     }
213 #endif
214 }
215
216 /* IBM's BlueGene(/L) MPI_Bcast dereferences the data pointer
217  * even when 0 == nbytes, so we protect calls to it on BlueGene.
218  * Fortunately dd_bcast() and dd_bcastc() are only
219  * called during DD setup and partition.
220  */
221
222 void dd_bcast(gmx_domdec_t gmx_unused *dd, int gmx_unused nbytes, void gmx_unused *data)
223 {
224 #ifdef GMX_MPI
225 #ifdef GMX_BLUEGENE
226     if (nbytes > 0)
227     {
228 #endif
229     MPI_Bcast(data, nbytes, MPI_BYTE,
230               DDMASTERRANK(dd), dd->mpi_comm_all);
231 #ifdef GMX_BLUEGENE
232 }
233 #endif
234 #endif
235 }
236
237 void dd_bcastc(gmx_domdec_t *dd, int nbytes, void *src, void *dest)
238 {
239     if (DDMASTER(dd))
240     {
241         memcpy(dest, src, nbytes);
242     }
243 #ifdef GMX_MPI
244 #ifdef GMX_BLUEGENE
245     if (nbytes > 0)
246     {
247 #endif
248     MPI_Bcast(dest, nbytes, MPI_BYTE,
249               DDMASTERRANK(dd), dd->mpi_comm_all);
250 #ifdef GMX_BLUEGENE
251 }
252 #endif
253 #endif
254 }
255
256 void dd_scatter(gmx_domdec_t gmx_unused *dd, int gmx_unused nbytes, void gmx_unused *src, void gmx_unused *dest)
257 {
258 #ifdef GMX_MPI
259     MPI_Scatter(src, nbytes, MPI_BYTE,
260                 dest, nbytes, MPI_BYTE,
261                 DDMASTERRANK(dd), dd->mpi_comm_all);
262 #endif
263 }
264
265 void dd_gather(gmx_domdec_t gmx_unused *dd, int gmx_unused nbytes, void gmx_unused *src, void gmx_unused *dest)
266 {
267 #ifdef GMX_MPI
268     MPI_Gather(src, nbytes, MPI_BYTE,
269                dest, nbytes, MPI_BYTE,
270                DDMASTERRANK(dd), dd->mpi_comm_all);
271 #endif
272 }
273
274 void dd_scatterv(gmx_domdec_t gmx_unused *dd,
275                  int gmx_unused *scounts, int gmx_unused *disps, void gmx_unused *sbuf,
276                  int gmx_unused rcount, void gmx_unused *rbuf)
277 {
278 #ifdef GMX_MPI
279     int dum;
280
281     if (rcount == 0)
282     {
283         /* MPI does not allow NULL pointers */
284         rbuf = &dum;
285     }
286     MPI_Scatterv(sbuf, scounts, disps, MPI_BYTE,
287                  rbuf, rcount, MPI_BYTE,
288                  DDMASTERRANK(dd), dd->mpi_comm_all);
289 #endif
290 }
291
292 void dd_gatherv(gmx_domdec_t gmx_unused *dd,
293                 int gmx_unused scount, void gmx_unused *sbuf,
294                 int gmx_unused *rcounts, int gmx_unused *disps, void gmx_unused *rbuf)
295 {
296 #ifdef GMX_MPI
297     int dum;
298
299     if (scount == 0)
300     {
301         /* MPI does not allow NULL pointers */
302         sbuf = &dum;
303     }
304     MPI_Gatherv(sbuf, scount, MPI_BYTE,
305                 rbuf, rcounts, disps, MPI_BYTE,
306                 DDMASTERRANK(dd), dd->mpi_comm_all);
307 #endif
308 }