2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
38 #include "gromacs/legacyheaders/domdec_network.h"
44 #include "gromacs/legacyheaders/types/commrec.h"
45 #include "gromacs/utility/gmxmpi.h"
48 #define DDMASTERRANK(dd) (dd->masterrank)
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)
60 rank_s = dd->neighbor[ddimind][direction == dddirForward ? 0 : 1];
61 rank_r = dd->neighbor[ddimind][direction == dddirForward ? 1 : 0];
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);
71 MPI_Send( buf_s, n_s*sizeof(int), MPI_BYTE, rank_s, 0,
76 MPI_Recv( buf_r, n_r*sizeof(int), MPI_BYTE, rank_r, 0,
77 dd->mpi_comm_all, &stat);
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)
92 rank_s = dd->neighbor[ddimind][direction == dddirForward ? 0 : 1];
93 rank_r = dd->neighbor[ddimind][direction == dddirForward ? 1 : 0];
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);
103 MPI_Send( buf_s, n_s*sizeof(real), MPI_BYTE, rank_s, 0,
108 MPI_Recv( buf_r, n_r*sizeof(real), MPI_BYTE, rank_r, 0,
109 dd->mpi_comm_all, &stat);
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)
124 rank_s = dd->neighbor[ddimind][direction == dddirForward ? 0 : 1];
125 rank_r = dd->neighbor[ddimind][direction == dddirForward ? 1 : 0];
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);
135 MPI_Send( buf_s[0], n_s*sizeof(rvec), MPI_BYTE, rank_s, 0,
140 MPI_Recv( buf_r[0], n_r*sizeof(rvec), MPI_BYTE, rank_r, 0,
141 dd->mpi_comm_all, &stat);
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)
155 int rank_fw, rank_bw, nreq;
159 rank_fw = dd->neighbor[ddimind][0];
160 rank_bw = dd->neighbor[ddimind][1];
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
170 * SendRecv2 can be turned on with the env.var. GMX_DD_SENDRECV2
175 MPI_Irecv(buf_r_fw[0], n_r_fw*sizeof(rvec), MPI_BYTE,
176 rank_bw, 0, dd->mpi_comm_all, &req[nreq++]);
180 MPI_Irecv(buf_r_bw[0], n_r_bw*sizeof(rvec), MPI_BYTE,
181 rank_fw, 1, dd->mpi_comm_all, &req[nreq++]);
185 MPI_Isend(buf_s_fw[0], n_s_fw*sizeof(rvec), MPI_BYTE,
186 rank_fw, 0, dd->mpi_comm_all, &req[nreq++]);
190 MPI_Isend(buf_s_bw[0], n_s_bw*sizeof(rvec), MPI_BYTE,
191 rank_bw, 1, dd->mpi_comm_all, &req[nreq++]);
195 MPI_Waitall(nreq, req, stat);
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.
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]);
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]);
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.
222 void dd_bcast(gmx_domdec_t gmx_unused *dd, int gmx_unused nbytes, void gmx_unused *data)
229 MPI_Bcast(data, nbytes, MPI_BYTE,
230 DDMASTERRANK(dd), dd->mpi_comm_all);
237 void dd_bcastc(gmx_domdec_t *dd, int nbytes, void *src, void *dest)
241 memcpy(dest, src, nbytes);
248 MPI_Bcast(dest, nbytes, MPI_BYTE,
249 DDMASTERRANK(dd), dd->mpi_comm_all);
256 void dd_scatter(gmx_domdec_t gmx_unused *dd, int gmx_unused nbytes, void gmx_unused *src, void gmx_unused *dest)
259 MPI_Scatter(src, nbytes, MPI_BYTE,
260 dest, nbytes, MPI_BYTE,
261 DDMASTERRANK(dd), dd->mpi_comm_all);
265 void dd_gather(gmx_domdec_t gmx_unused *dd, int gmx_unused nbytes, void gmx_unused *src, void gmx_unused *dest)
268 MPI_Gather(src, nbytes, MPI_BYTE,
269 dest, nbytes, MPI_BYTE,
270 DDMASTERRANK(dd), dd->mpi_comm_all);
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)
283 /* MPI does not allow NULL pointers */
286 MPI_Scatterv(sbuf, scounts, disps, MPI_BYTE,
287 rbuf, rcount, MPI_BYTE,
288 DDMASTERRANK(dd), dd->mpi_comm_all);
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)
301 /* MPI does not allow NULL pointers */
304 MPI_Gatherv(sbuf, scount, MPI_BYTE,
305 rbuf, rcounts, disps, MPI_BYTE,
306 DDMASTERRANK(dd), dd->mpi_comm_all);