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