Merge branch 'release-4-6'
[alexxy/gromacs.git] / src / gromacs / mdlib / domdec_network.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *
4  * This file is part of Gromacs        Copyright (c) 1991-2008
5  * David van der Spoel, Erik Lindahl, Berk Hess, University of Groningen.
6  *
7  * This program is free software; you can redistribute it and/or
8  * modify it under the terms of the GNU General Public License
9  * as published by the Free Software Foundation; either version 2
10  * of the License, or (at your option) any later version.
11  *
12  * To help us fund GROMACS development, we humbly ask that you cite
13  * the research papers on the package. Check out http://www.gromacs.org
14  *
15  * And Hey:
16  * Gnomes, ROck Monsters And Chili Sauce
17  */
18
19 #ifdef HAVE_CONFIG_H
20 #include <config.h>
21 #endif
22
23 #include <string.h>
24 #include "domdec_network.h"
25
26 #include "gromacs/utility/gmxmpi.h"
27
28
29 #define DDMASTERRANK(dd)   (dd->masterrank)
30
31
32 void dd_sendrecv_int(const gmx_domdec_t *dd,
33                      int ddimind, int direction,
34                      int *buf_s, int n_s,
35                      int *buf_r, int n_r)
36 {
37 #ifdef GMX_MPI
38     int        rank_s, rank_r;
39     MPI_Status stat;
40
41     rank_s = dd->neighbor[ddimind][direction == dddirForward ? 0 : 1];
42     rank_r = dd->neighbor[ddimind][direction == dddirForward ? 1 : 0];
43
44     if (n_s && n_r)
45     {
46         MPI_Sendrecv(buf_s, n_s*sizeof(int), MPI_BYTE, rank_s, 0,
47                      buf_r, n_r*sizeof(int), MPI_BYTE, rank_r, 0,
48                      dd->mpi_comm_all, &stat);
49     }
50     else if (n_s)
51     {
52         MPI_Send(    buf_s, n_s*sizeof(int), MPI_BYTE, rank_s, 0,
53                      dd->mpi_comm_all);
54     }
55     else if (n_r)
56     {
57         MPI_Recv(    buf_r, n_r*sizeof(int), MPI_BYTE, rank_r, 0,
58                      dd->mpi_comm_all, &stat);
59     }
60
61 #endif
62 }
63
64 void dd_sendrecv_real(const gmx_domdec_t *dd,
65                       int ddimind, int direction,
66                       real *buf_s, int n_s,
67                       real *buf_r, int n_r)
68 {
69 #ifdef GMX_MPI
70     int        rank_s, rank_r;
71     MPI_Status stat;
72
73     rank_s = dd->neighbor[ddimind][direction == dddirForward ? 0 : 1];
74     rank_r = dd->neighbor[ddimind][direction == dddirForward ? 1 : 0];
75
76     if (n_s && n_r)
77     {
78         MPI_Sendrecv(buf_s, n_s*sizeof(real), MPI_BYTE, rank_s, 0,
79                      buf_r, n_r*sizeof(real), MPI_BYTE, rank_r, 0,
80                      dd->mpi_comm_all, &stat);
81     }
82     else if (n_s)
83     {
84         MPI_Send(    buf_s, n_s*sizeof(real), MPI_BYTE, rank_s, 0,
85                      dd->mpi_comm_all);
86     }
87     else if (n_r)
88     {
89         MPI_Recv(    buf_r, n_r*sizeof(real), MPI_BYTE, rank_r, 0,
90                      dd->mpi_comm_all, &stat);
91     }
92
93 #endif
94 }
95
96 void dd_sendrecv_rvec(const gmx_domdec_t *dd,
97                       int ddimind, int direction,
98                       rvec *buf_s, int n_s,
99                       rvec *buf_r, int n_r)
100 {
101 #ifdef GMX_MPI
102     int        rank_s, rank_r;
103     MPI_Status stat;
104
105     rank_s = dd->neighbor[ddimind][direction == dddirForward ? 0 : 1];
106     rank_r = dd->neighbor[ddimind][direction == dddirForward ? 1 : 0];
107
108     if (n_s && n_r)
109     {
110         MPI_Sendrecv(buf_s[0], n_s*sizeof(rvec), MPI_BYTE, rank_s, 0,
111                      buf_r[0], n_r*sizeof(rvec), MPI_BYTE, rank_r, 0,
112                      dd->mpi_comm_all, &stat);
113     }
114     else if (n_s)
115     {
116         MPI_Send(    buf_s[0], n_s*sizeof(rvec), MPI_BYTE, rank_s, 0,
117                      dd->mpi_comm_all);
118     }
119     else if (n_r)
120     {
121         MPI_Recv(    buf_r[0], n_r*sizeof(rvec), MPI_BYTE, rank_r, 0,
122                      dd->mpi_comm_all, &stat);
123     }
124
125 #endif
126 }
127
128 void dd_sendrecv2_rvec(const gmx_domdec_t *dd,
129                        int ddimind,
130                        rvec *buf_s_fw, int n_s_fw,
131                        rvec *buf_r_fw, int n_r_fw,
132                        rvec *buf_s_bw, int n_s_bw,
133                        rvec *buf_r_bw, int n_r_bw)
134 {
135 #ifdef GMX_MPI
136     int         rank_fw, rank_bw, nreq;
137     MPI_Request req[4];
138     MPI_Status  stat[4];
139
140     rank_fw = dd->neighbor[ddimind][0];
141     rank_bw = dd->neighbor[ddimind][1];
142
143     if (!dd->bSendRecv2)
144     {
145         /* Try to send and receive in two directions simultaneously.
146          * Should be faster, especially on machines
147          * with full 3D communication networks.
148          * However, it could be that communication libraries are
149          * optimized for MPI_Sendrecv and non-blocking MPI calls
150          * are slower.
151          * SendRecv2 can be turned on with the env.var. GMX_DD_SENDRECV2
152          */
153         nreq = 0;
154         if (n_r_fw)
155         {
156             MPI_Irecv(buf_r_fw[0], n_r_fw*sizeof(rvec), MPI_BYTE,
157                       rank_bw, 0, dd->mpi_comm_all, &req[nreq++]);
158         }
159         if (n_r_bw)
160         {
161             MPI_Irecv(buf_r_bw[0], n_r_bw*sizeof(rvec), MPI_BYTE,
162                       rank_fw, 1, dd->mpi_comm_all, &req[nreq++]);
163         }
164         if (n_s_fw)
165         {
166             MPI_Isend(buf_s_fw[0], n_s_fw*sizeof(rvec), MPI_BYTE,
167                       rank_fw, 0, dd->mpi_comm_all, &req[nreq++]);
168         }
169         if (n_s_bw)
170         {
171             MPI_Isend(buf_s_bw[0], n_s_bw*sizeof(rvec), MPI_BYTE,
172                       rank_bw, 1, dd->mpi_comm_all, &req[nreq++]);
173         }
174         if (nreq)
175         {
176             MPI_Waitall(nreq, req, stat);
177         }
178     }
179     else
180     {
181         /* Communicate in two ordered phases.
182          * This is slower, even on a dual-core Opteron cluster
183          * with a single full-duplex network connection per machine.
184          */
185         /* Forward */
186         MPI_Sendrecv(buf_s_fw[0], n_s_fw*sizeof(rvec), MPI_BYTE, rank_fw, 0,
187                      buf_r_fw[0], n_r_fw*sizeof(rvec), MPI_BYTE, rank_bw, 0,
188                      dd->mpi_comm_all, &stat[0]);
189         /* Backward */
190         MPI_Sendrecv(buf_s_bw[0], n_s_bw*sizeof(rvec), MPI_BYTE, rank_bw, 0,
191                      buf_r_bw[0], n_r_bw*sizeof(rvec), MPI_BYTE, rank_fw, 0,
192                      dd->mpi_comm_all, &stat[0]);
193     }
194 #endif
195 }
196
197 /* IBM's BlueGene(/L) MPI_Bcast dereferences the data pointer
198  * even when 0 == nbytes, so we protect calls to it on BlueGene.
199  * Fortunately dd_bcast() and dd_bcastc() are only
200  * called during DD setup and partition.
201  */
202
203 void dd_bcast(gmx_domdec_t *dd, int nbytes, void *data)
204 {
205 #ifdef GMX_MPI
206 #ifdef GMX_BLUEGENE
207     if (nbytes > 0)
208     {
209 #endif
210     MPI_Bcast(data, nbytes, MPI_BYTE,
211               DDMASTERRANK(dd), dd->mpi_comm_all);
212 #ifdef GMX_BLUEGENE
213 }
214 #endif
215 #endif
216 }
217
218 void dd_bcastc(gmx_domdec_t *dd, int nbytes, void *src, void *dest)
219 {
220     if (DDMASTER(dd))
221     {
222         memcpy(dest, src, nbytes);
223     }
224 #ifdef GMX_MPI
225 #ifdef GMX_BLUEGENE
226     if (nbytes > 0)
227     {
228 #endif
229     MPI_Bcast(dest, nbytes, MPI_BYTE,
230               DDMASTERRANK(dd), dd->mpi_comm_all);
231 #ifdef GMX_BLUEGENE
232 }
233 #endif
234 #endif
235 }
236
237 void dd_scatter(gmx_domdec_t *dd, int nbytes, void *src, void *dest)
238 {
239 #ifdef GMX_MPI
240     MPI_Scatter(src, nbytes, MPI_BYTE,
241                 dest, nbytes, MPI_BYTE,
242                 DDMASTERRANK(dd), dd->mpi_comm_all);
243 #endif
244 }
245
246 void dd_gather(gmx_domdec_t *dd, int nbytes, void *src, void *dest)
247 {
248 #ifdef GMX_MPI
249     MPI_Gather(src, nbytes, MPI_BYTE,
250                dest, nbytes, MPI_BYTE,
251                DDMASTERRANK(dd), dd->mpi_comm_all);
252 #endif
253 }
254
255 void dd_scatterv(gmx_domdec_t *dd,
256                  int *scounts, int *disps, void *sbuf,
257                  int rcount, void *rbuf)
258 {
259 #ifdef GMX_MPI
260     int dum;
261
262     if (rcount == 0)
263     {
264         /* MPI does not allow NULL pointers */
265         rbuf = &dum;
266     }
267     MPI_Scatterv(sbuf, scounts, disps, MPI_BYTE,
268                  rbuf, rcount, MPI_BYTE,
269                  DDMASTERRANK(dd), dd->mpi_comm_all);
270 #endif
271 }
272
273 void dd_gatherv(gmx_domdec_t *dd,
274                 int scount, void *sbuf,
275                 int *rcounts, int *disps, void *rbuf)
276 {
277 #ifdef GMX_MPI
278     int dum;
279
280     if (scount == 0)
281     {
282         /* MPI does not allow NULL pointers */
283         sbuf = &dum;
284     }
285     MPI_Gatherv(sbuf, scount, MPI_BYTE,
286                 rbuf, rcounts, disps, MPI_BYTE,
287                 DDMASTERRANK(dd), dd->mpi_comm_all);
288 #endif
289 }