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