2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2008
5 * Copyright (c) 2012,2013, by the GROMACS development team, led by
6 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
7 * others, as listed in the AUTHORS file in the top-level source
8 * directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
42 #include "domdec_network.h"
52 #define DDMASTERRANK(dd) (dd->masterrank)
55 void dd_sendrecv_int(const gmx_domdec_t *dd,
56 int ddimind,int direction,
64 rank_s = dd->neighbor[ddimind][direction==dddirForward ? 0 : 1];
65 rank_r = dd->neighbor[ddimind][direction==dddirForward ? 1 : 0];
69 MPI_Sendrecv(buf_s,n_s*sizeof(int),MPI_BYTE,rank_s,0,
70 buf_r,n_r*sizeof(int),MPI_BYTE,rank_r,0,
71 dd->mpi_comm_all,&stat);
75 MPI_Send( buf_s,n_s*sizeof(int),MPI_BYTE,rank_s,0,
80 MPI_Recv( buf_r,n_r*sizeof(int),MPI_BYTE,rank_r,0,
81 dd->mpi_comm_all,&stat);
87 void dd_sendrecv_real(const gmx_domdec_t *dd,
88 int ddimind,int direction,
96 rank_s = dd->neighbor[ddimind][direction==dddirForward ? 0 : 1];
97 rank_r = dd->neighbor[ddimind][direction==dddirForward ? 1 : 0];
101 MPI_Sendrecv(buf_s,n_s*sizeof(real),MPI_BYTE,rank_s,0,
102 buf_r,n_r*sizeof(real),MPI_BYTE,rank_r,0,
103 dd->mpi_comm_all,&stat);
107 MPI_Send( buf_s,n_s*sizeof(real),MPI_BYTE,rank_s,0,
112 MPI_Recv( buf_r,n_r*sizeof(real),MPI_BYTE,rank_r,0,
113 dd->mpi_comm_all,&stat);
119 void dd_sendrecv_rvec(const gmx_domdec_t *dd,
120 int ddimind,int direction,
128 rank_s = dd->neighbor[ddimind][direction==dddirForward ? 0 : 1];
129 rank_r = dd->neighbor[ddimind][direction==dddirForward ? 1 : 0];
133 MPI_Sendrecv(buf_s[0],n_s*sizeof(rvec),MPI_BYTE,rank_s,0,
134 buf_r[0],n_r*sizeof(rvec),MPI_BYTE,rank_r,0,
135 dd->mpi_comm_all,&stat);
139 MPI_Send( buf_s[0],n_s*sizeof(rvec),MPI_BYTE,rank_s,0,
143 MPI_Recv( buf_r[0],n_r*sizeof(rvec),MPI_BYTE,rank_r,0,
144 dd->mpi_comm_all,&stat);
150 void dd_sendrecv2_rvec(const gmx_domdec_t *dd,
152 rvec *buf_s_fw,int n_s_fw,
153 rvec *buf_r_fw,int n_r_fw,
154 rvec *buf_s_bw,int n_s_bw,
155 rvec *buf_r_bw,int n_r_bw)
158 int rank_fw,rank_bw,nreq;
162 rank_fw = dd->neighbor[ddimind][0];
163 rank_bw = dd->neighbor[ddimind][1];
167 /* Try to send and receive in two directions simultaneously.
168 * Should be faster, especially on machines
169 * with full 3D communication networks.
170 * However, it could be that communication libraries are
171 * optimized for MPI_Sendrecv and non-blocking MPI calls
173 * SendRecv2 can be turned on with the env.var. GMX_DD_SENDRECV2
178 MPI_Irecv(buf_r_fw[0],n_r_fw*sizeof(rvec),MPI_BYTE,
179 rank_bw,0,dd->mpi_comm_all,&req[nreq++]);
183 MPI_Irecv(buf_r_bw[0],n_r_bw*sizeof(rvec),MPI_BYTE,
184 rank_fw,1,dd->mpi_comm_all,&req[nreq++]);
188 MPI_Isend(buf_s_fw[0],n_s_fw*sizeof(rvec),MPI_BYTE,
189 rank_fw,0,dd->mpi_comm_all,&req[nreq++]);
193 MPI_Isend(buf_s_bw[0],n_s_bw*sizeof(rvec),MPI_BYTE,
194 rank_bw,1,dd->mpi_comm_all,&req[nreq++]);
198 MPI_Waitall(nreq,req,stat);
203 /* Communicate in two ordered phases.
204 * This is slower, even on a dual-core Opteron cluster
205 * with a single full-duplex network connection per machine.
208 MPI_Sendrecv(buf_s_fw[0],n_s_fw*sizeof(rvec),MPI_BYTE,rank_fw,0,
209 buf_r_fw[0],n_r_fw*sizeof(rvec),MPI_BYTE,rank_bw,0,
210 dd->mpi_comm_all,&stat[0]);
212 MPI_Sendrecv(buf_s_bw[0],n_s_bw*sizeof(rvec),MPI_BYTE,rank_bw,0,
213 buf_r_bw[0],n_r_bw*sizeof(rvec),MPI_BYTE,rank_fw,0,
214 dd->mpi_comm_all,&stat[0]);
219 /* IBM's BlueGene(/L) MPI_Bcast dereferences the data pointer
220 * even when 0 == nbytes, so we protect calls to it on BlueGene.
221 * Fortunately dd_bcast() and dd_bcastc() are only
222 * called during DD setup and partition.
225 void dd_bcast(gmx_domdec_t *dd,int nbytes,void *data)
232 MPI_Bcast(data,nbytes,MPI_BYTE,
233 DDMASTERRANK(dd),dd->mpi_comm_all);
240 void dd_bcastc(gmx_domdec_t *dd,int nbytes,void *src,void *dest)
244 memcpy(dest,src,nbytes);
251 MPI_Bcast(dest,nbytes,MPI_BYTE,
252 DDMASTERRANK(dd),dd->mpi_comm_all);
259 void dd_scatter(gmx_domdec_t *dd,int nbytes,void *src,void *dest)
262 MPI_Scatter(src,nbytes,MPI_BYTE,
263 dest,nbytes,MPI_BYTE,
264 DDMASTERRANK(dd),dd->mpi_comm_all);
268 void dd_gather(gmx_domdec_t *dd,int nbytes,void *src,void *dest)
271 MPI_Gather(src,nbytes,MPI_BYTE,
272 dest,nbytes,MPI_BYTE,
273 DDMASTERRANK(dd),dd->mpi_comm_all);
277 void dd_scatterv(gmx_domdec_t *dd,
278 int *scounts,int *disps,void *sbuf,
279 int rcount,void *rbuf)
286 /* MPI does not allow NULL pointers */
289 MPI_Scatterv(sbuf,scounts,disps,MPI_BYTE,
290 rbuf,rcount,MPI_BYTE,
291 DDMASTERRANK(dd),dd->mpi_comm_all);
295 void dd_gatherv(gmx_domdec_t *dd,
296 int scount,void *sbuf,
297 int *rcounts,int *disps,void *rbuf)
304 /* MPI does not allow NULL pointers */
307 MPI_Gatherv(sbuf,scount,MPI_BYTE,
308 rbuf,rcounts,disps,MPI_BYTE,
309 DDMASTERRANK(dd),dd->mpi_comm_all);