1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This file is part of Gromacs Copyright (c) 1991-2008
5 * David van der Spoel, Erik Lindahl, Berk Hess, University of Groningen.
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.
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
16 * Gnomes, ROck Monsters And Chili Sauce
24 #include "domdec_network.h"
26 #include "gromacs/utility/gmxmpi.h"
29 #define DDMASTERRANK(dd) (dd->masterrank)
32 void dd_sendrecv_int(const gmx_domdec_t *dd,
33 int ddimind, int direction,
41 rank_s = dd->neighbor[ddimind][direction == dddirForward ? 0 : 1];
42 rank_r = dd->neighbor[ddimind][direction == dddirForward ? 1 : 0];
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);
52 MPI_Send( buf_s, n_s*sizeof(int), MPI_BYTE, rank_s, 0,
57 MPI_Recv( buf_r, n_r*sizeof(int), MPI_BYTE, rank_r, 0,
58 dd->mpi_comm_all, &stat);
64 void dd_sendrecv_real(const gmx_domdec_t *dd,
65 int ddimind, int direction,
73 rank_s = dd->neighbor[ddimind][direction == dddirForward ? 0 : 1];
74 rank_r = dd->neighbor[ddimind][direction == dddirForward ? 1 : 0];
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);
84 MPI_Send( buf_s, n_s*sizeof(real), MPI_BYTE, rank_s, 0,
89 MPI_Recv( buf_r, n_r*sizeof(real), MPI_BYTE, rank_r, 0,
90 dd->mpi_comm_all, &stat);
96 void dd_sendrecv_rvec(const gmx_domdec_t *dd,
97 int ddimind, int direction,
105 rank_s = dd->neighbor[ddimind][direction == dddirForward ? 0 : 1];
106 rank_r = dd->neighbor[ddimind][direction == dddirForward ? 1 : 0];
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);
116 MPI_Send( buf_s[0], n_s*sizeof(rvec), MPI_BYTE, rank_s, 0,
121 MPI_Recv( buf_r[0], n_r*sizeof(rvec), MPI_BYTE, rank_r, 0,
122 dd->mpi_comm_all, &stat);
128 void dd_sendrecv2_rvec(const gmx_domdec_t *dd,
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)
136 int rank_fw, rank_bw, nreq;
140 rank_fw = dd->neighbor[ddimind][0];
141 rank_bw = dd->neighbor[ddimind][1];
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
151 * SendRecv2 can be turned on with the env.var. GMX_DD_SENDRECV2
156 MPI_Irecv(buf_r_fw[0], n_r_fw*sizeof(rvec), MPI_BYTE,
157 rank_bw, 0, dd->mpi_comm_all, &req[nreq++]);
161 MPI_Irecv(buf_r_bw[0], n_r_bw*sizeof(rvec), MPI_BYTE,
162 rank_fw, 1, dd->mpi_comm_all, &req[nreq++]);
166 MPI_Isend(buf_s_fw[0], n_s_fw*sizeof(rvec), MPI_BYTE,
167 rank_fw, 0, dd->mpi_comm_all, &req[nreq++]);
171 MPI_Isend(buf_s_bw[0], n_s_bw*sizeof(rvec), MPI_BYTE,
172 rank_bw, 1, dd->mpi_comm_all, &req[nreq++]);
176 MPI_Waitall(nreq, req, stat);
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.
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]);
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]);
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.
203 void dd_bcast(gmx_domdec_t *dd, int nbytes, void *data)
210 MPI_Bcast(data, nbytes, MPI_BYTE,
211 DDMASTERRANK(dd), dd->mpi_comm_all);
218 void dd_bcastc(gmx_domdec_t *dd, int nbytes, void *src, void *dest)
222 memcpy(dest, src, nbytes);
229 MPI_Bcast(dest, nbytes, MPI_BYTE,
230 DDMASTERRANK(dd), dd->mpi_comm_all);
237 void dd_scatter(gmx_domdec_t *dd, int nbytes, void *src, void *dest)
240 MPI_Scatter(src, nbytes, MPI_BYTE,
241 dest, nbytes, MPI_BYTE,
242 DDMASTERRANK(dd), dd->mpi_comm_all);
246 void dd_gather(gmx_domdec_t *dd, int nbytes, void *src, void *dest)
249 MPI_Gather(src, nbytes, MPI_BYTE,
250 dest, nbytes, MPI_BYTE,
251 DDMASTERRANK(dd), dd->mpi_comm_all);
255 void dd_scatterv(gmx_domdec_t *dd,
256 int *scounts, int *disps, void *sbuf,
257 int rcount, void *rbuf)
264 /* MPI does not allow NULL pointers */
267 MPI_Scatterv(sbuf, scounts, disps, MPI_BYTE,
268 rbuf, rcount, MPI_BYTE,
269 DDMASTERRANK(dd), dd->mpi_comm_all);
273 void dd_gatherv(gmx_domdec_t *dd,
274 int scount, void *sbuf,
275 int *rcounts, int *disps, void *rbuf)
282 /* MPI does not allow NULL pointers */
285 MPI_Gatherv(sbuf, scount, MPI_BYTE,
286 rbuf, rcounts, disps, MPI_BYTE,
287 DDMASTERRANK(dd), dd->mpi_comm_all);