2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
44 #include "types/commrec.h"
49 #include "gromacs/fileio/futil.h"
50 #include "gromacs/utility/basenetwork.h"
51 #include "gromacs/utility/cstringutil.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/gmxmpi.h"
54 #include "gromacs/utility/smalloc.h"
56 /* The source code in this file should be thread-safe.
57 Please keep it that way. */
59 void gmx_fill_commrec_from_mpi(t_commrec gmx_unused *cr)
62 gmx_call("gmx_fill_commrec_from_mpi");
64 if (!gmx_mpi_initialized())
66 gmx_comm("MPI has not been initialized properly");
69 cr->nnodes = gmx_node_num();
70 cr->nodeid = gmx_node_rank();
71 cr->sim_nodeid = cr->nodeid;
72 cr->mpi_comm_mysim = MPI_COMM_WORLD;
73 cr->mpi_comm_mygroup = MPI_COMM_WORLD;
78 t_commrec *init_commrec()
85 gmx_fill_commrec_from_mpi(cr);
87 cr->mpi_comm_mysim = NULL;
88 cr->mpi_comm_mygroup = NULL;
91 cr->nodeid = cr->sim_nodeid;
94 // TODO cr->duty should not be initialized here
95 cr->duty = (DUTY_PP | DUTY_PME);
97 #if defined GMX_MPI && !defined MPI_IN_PLACE_EXISTS
98 /* initialize the MPI_IN_PLACE replacement buffers */
100 cr->mpb->ibuf = NULL;
101 cr->mpb->libuf = NULL;
102 cr->mpb->fbuf = NULL;
103 cr->mpb->dbuf = NULL;
104 cr->mpb->ibuf_alloc = 0;
105 cr->mpb->libuf_alloc = 0;
106 cr->mpb->fbuf_alloc = 0;
107 cr->mpb->dbuf_alloc = 0;
113 t_commrec *reinitialize_commrec_for_this_thread(const t_commrec gmx_unused *cro)
115 #ifdef GMX_THREAD_MPI
118 /* make a thread-specific commrec */
120 /* now copy the whole thing, so settings like the number of PME nodes
124 /* and we start setting our own thread-specific values for things */
125 gmx_fill_commrec_from_mpi(cr);
127 // TODO cr->duty should not be initialized here
128 cr->duty = (DUTY_PP | DUTY_PME);
136 void gmx_setup_nodecomm(FILE gmx_unused *fplog, t_commrec *cr)
139 int n, rank, hostnum, ng, ni;
141 /* Many MPI implementations do not optimize MPI_Allreduce
142 * (and probably also other global communication calls)
143 * for multi-core nodes connected by a network.
144 * We can optimize such communication by using one MPI call
145 * within each node and one between the nodes.
146 * For MVAPICH2 and Intel MPI this reduces the time for
147 * the global_stat communication by 25%
148 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
149 * B. Hess, November 2007
155 #ifndef GMX_THREAD_MPI
157 MPI_Comm_size(cr->mpi_comm_mygroup, &n);
158 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
160 hostnum = gmx_hostname_num();
164 fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
168 /* The intra-node communicator, split on node number */
169 MPI_Comm_split(cr->mpi_comm_mygroup, hostnum, rank, &nc->comm_intra);
170 MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
173 fprintf(debug, "In gmx_setup_nodecomm: node rank %d rank_intra %d\n",
174 rank, nc->rank_intra);
176 /* The inter-node communicator, split on rank_intra.
177 * We actually only need the one for rank=0,
178 * but it is easier to create them all.
180 MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
181 /* Check if this really created two step communication */
182 MPI_Comm_size(nc->comm_inter, &ng);
183 MPI_Comm_size(nc->comm_intra, &ni);
186 fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
190 if (getenv("GMX_NO_NODECOMM") == NULL &&
191 ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
196 fprintf(fplog, "Using two step summing over %d groups of on average %.1f processes\n\n",
197 ng, (real)n/(real)ng);
199 if (nc->rank_intra > 0)
201 MPI_Comm_free(&nc->comm_inter);
206 /* One group or all processes in a separate group, use normal summing */
207 MPI_Comm_free(&nc->comm_inter);
208 MPI_Comm_free(&nc->comm_intra);
211 fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
216 /* tMPI runs only on a single node so just use the nodeid */
217 nc->rank_intra = cr->nodeid;
221 void gmx_init_intranode_counters(t_commrec *cr)
223 /* counters for PP+PME and PP-only processes on my physical node */
224 int nrank_intranode, rank_intranode;
225 int nrank_pp_intranode, rank_pp_intranode;
226 /* thread-MPI is not initialized when not running in parallel */
227 #if defined GMX_MPI && !defined GMX_THREAD_MPI
228 int nrank_world, rank_world;
229 int i, mynum, *num, *num_s, *num_pp, *num_pp_s;
231 MPI_Comm_size(MPI_COMM_WORLD, &nrank_world);
232 MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
234 /* Get the node number from the hostname to identify the nodes */
235 mynum = gmx_hostname_num();
237 /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
238 snew(num, nrank_world);
239 snew(num_s, nrank_world);
240 snew(num_pp, nrank_world);
241 snew(num_pp_s, nrank_world);
243 num_s[rank_world] = mynum;
244 num_pp_s[rank_world] = (cr->duty & DUTY_PP) ? mynum : -1;
246 MPI_Allreduce(num_s, num, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
247 MPI_Allreduce(num_pp_s, num_pp, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
251 nrank_pp_intranode = 0;
252 rank_pp_intranode = 0;
253 for (i = 0; i < nrank_world; i++)
263 if (num_pp[i] == mynum)
265 nrank_pp_intranode++;
266 if ((cr->duty & DUTY_PP) && i < rank_world)
277 /* Serial or thread-MPI code: we run within a single physical node */
278 nrank_intranode = cr->nnodes;
279 rank_intranode = cr->sim_nodeid;
280 nrank_pp_intranode = cr->nnodes - cr->npmenodes;
281 rank_pp_intranode = cr->nodeid;
287 if (cr->duty & DUTY_PP && cr->duty & DUTY_PME)
289 sprintf(sbuf, "PP+PME");
293 sprintf(sbuf, "%s", cr->duty & DUTY_PP ? "PP" : "PME");
295 fprintf(debug, "On %3s node %d: nrank_intranode=%d, rank_intranode=%d, "
296 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
297 sbuf, cr->sim_nodeid,
298 nrank_intranode, rank_intranode,
299 nrank_pp_intranode, rank_pp_intranode);
302 cr->nrank_intranode = nrank_intranode;
303 cr->rank_intranode = rank_intranode;
304 cr->nrank_pp_intranode = nrank_pp_intranode;
305 cr->rank_pp_intranode = rank_pp_intranode;
309 void gmx_barrier(const t_commrec gmx_unused *cr)
312 gmx_call("gmx_barrier");
314 MPI_Barrier(cr->mpi_comm_mygroup);
318 void gmx_bcast(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
321 gmx_call("gmx_bast");
323 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
327 void gmx_bcast_sim(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
330 gmx_call("gmx_bast");
332 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
336 void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unused *cr)
339 gmx_call("gmx_sumd");
341 #if defined(MPI_IN_PLACE_EXISTS)
344 if (cr->nc.rank_intra == 0)
346 /* Use two step summing. */
347 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
349 /* Sum the roots of the internal (intra) buffers. */
350 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
355 /* This is here because of the silly MPI specification
356 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
357 MPI_Reduce(r, NULL, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
359 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
363 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
364 cr->mpi_comm_mygroup);
369 if (nr > cr->mpb->dbuf_alloc)
371 cr->mpb->dbuf_alloc = nr;
372 srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
376 /* Use two step summing */
377 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
378 if (cr->nc.rank_intra == 0)
380 /* Sum with the buffers reversed */
381 MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
384 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
388 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
389 cr->mpi_comm_mygroup);
390 for (i = 0; i < nr; i++)
392 r[i] = cr->mpb->dbuf[i];
399 void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unused *cr)
402 gmx_call("gmx_sumf");
404 #if defined(MPI_IN_PLACE_EXISTS)
407 /* Use two step summing. */
408 if (cr->nc.rank_intra == 0)
410 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
412 /* Sum the roots of the internal (intra) buffers */
413 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
418 /* This is here because of the silly MPI specification
419 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
420 MPI_Reduce(r, NULL, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
422 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
426 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
431 if (nr > cr->mpb->fbuf_alloc)
433 cr->mpb->fbuf_alloc = nr;
434 srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
438 /* Use two step summing */
439 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
440 if (cr->nc.rank_intra == 0)
442 /* Sum with the buffers reversed */
443 MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
446 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
450 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
451 cr->mpi_comm_mygroup);
452 for (i = 0; i < nr; i++)
454 r[i] = cr->mpb->fbuf[i];
461 void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused *cr)
464 gmx_call("gmx_sumi");
466 #if defined(MPI_IN_PLACE_EXISTS)
469 /* Use two step summing */
470 if (cr->nc.rank_intra == 0)
472 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
473 /* Sum with the buffers reversed */
474 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
478 /* This is here because of the silly MPI specification
479 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
480 MPI_Reduce(r, NULL, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
482 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
486 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
491 if (nr > cr->mpb->ibuf_alloc)
493 cr->mpb->ibuf_alloc = nr;
494 srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
498 /* Use two step summing */
499 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
500 if (cr->nc.rank_intra == 0)
502 /* Sum with the buffers reversed */
503 MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
505 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
509 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
510 for (i = 0; i < nr; i++)
512 r[i] = cr->mpb->ibuf[i];
519 void gmx_sumli(int gmx_unused nr, gmx_int64_t gmx_unused r[], const t_commrec gmx_unused *cr)
522 gmx_call("gmx_sumli");
524 #if defined(MPI_IN_PLACE_EXISTS)
527 /* Use two step summing */
528 if (cr->nc.rank_intra == 0)
530 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, 0,
532 /* Sum with the buffers reversed */
533 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
538 /* This is here because of the silly MPI specification
539 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
540 MPI_Reduce(r, NULL, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
542 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
546 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
551 if (nr > cr->mpb->libuf_alloc)
553 cr->mpb->libuf_alloc = nr;
554 srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
558 /* Use two step summing */
559 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
561 if (cr->nc.rank_intra == 0)
563 /* Sum with the buffers reversed */
564 MPI_Allreduce(cr->mpb->libuf, r, nr, MPI_INT64_T, MPI_SUM,
567 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
571 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
572 cr->mpi_comm_mygroup);
573 for (i = 0; i < nr; i++)
575 r[i] = cr->mpb->libuf[i];
585 static void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
587 #if defined(MPI_IN_PLACE_EXISTS)
588 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
590 /* this function is only used in code that is not performance critical,
591 (during setup, when comm_rec is not the appropriate communication
592 structure), so this isn't as bad as it looks. */
597 MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
598 for (i = 0; i < nr; i++)
608 static void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
610 #if defined(MPI_IN_PLACE_EXISTS)
611 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
613 /* this function is only used in code that is not performance critical,
614 (during setup, when comm_rec is not the appropriate communication
615 structure), so this isn't as bad as it looks. */
620 MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
621 for (i = 0; i < nr; i++)
630 void gmx_sumd_sim(int gmx_unused nr, double gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
633 gmx_call("gmx_sumd_sim");
635 gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
639 void gmx_sumf_sim(int gmx_unused nr, float gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
642 gmx_call("gmx_sumf_sim");
644 gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
648 void gmx_sumi_sim(int gmx_unused nr, int gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
651 gmx_call("gmx_sumi_sim");
653 #if defined(MPI_IN_PLACE_EXISTS)
654 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
656 /* this is thread-unsafe, but it will do for now: */
659 if (nr > ms->mpb->ibuf_alloc)
661 ms->mpb->ibuf_alloc = nr;
662 srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
664 MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
665 for (i = 0; i < nr; i++)
667 r[i] = ms->mpb->ibuf[i];
673 void gmx_sumli_sim(int gmx_unused nr, gmx_int64_t gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
676 gmx_call("gmx_sumli_sim");
678 #if defined(MPI_IN_PLACE_EXISTS)
679 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
680 ms->mpi_comm_masters);
682 /* this is thread-unsafe, but it will do for now: */
685 if (nr > ms->mpb->libuf_alloc)
687 ms->mpb->libuf_alloc = nr;
688 srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
690 MPI_Allreduce(r, ms->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
691 ms->mpi_comm_masters);
692 for (i = 0; i < nr; i++)
694 r[i] = ms->mpb->libuf[i];
700 gmx_bool gmx_fexist_master(const char *fname, t_commrec *cr)
706 bExist = gmx_fexist(fname);
710 gmx_bcast(sizeof(bExist), &bExist, cr);