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.
46 #include "types/commrec.h"
50 #include "gromacs/utility/basenetwork.h"
51 #include "gromacs/utility/cstringutil.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/futil.h"
54 #include "gromacs/utility/gmxmpi.h"
55 #include "gromacs/utility/smalloc.h"
57 /* The source code in this file should be thread-safe.
58 Please keep it that way. */
60 void gmx_fill_commrec_from_mpi(t_commrec gmx_unused *cr)
63 gmx_call("gmx_fill_commrec_from_mpi");
65 if (!gmx_mpi_initialized())
67 gmx_comm("MPI has not been initialized properly");
70 cr->nnodes = gmx_node_num();
71 cr->nodeid = gmx_node_rank();
72 cr->sim_nodeid = cr->nodeid;
73 cr->mpi_comm_mysim = MPI_COMM_WORLD;
74 cr->mpi_comm_mygroup = MPI_COMM_WORLD;
79 t_commrec *init_commrec()
86 gmx_fill_commrec_from_mpi(cr);
88 cr->mpi_comm_mysim = NULL;
89 cr->mpi_comm_mygroup = NULL;
92 cr->nodeid = cr->sim_nodeid;
95 // TODO cr->duty should not be initialized here
96 cr->duty = (DUTY_PP | DUTY_PME);
98 #if defined GMX_MPI && !defined MPI_IN_PLACE_EXISTS
99 /* initialize the MPI_IN_PLACE replacement buffers */
101 cr->mpb->ibuf = NULL;
102 cr->mpb->libuf = NULL;
103 cr->mpb->fbuf = NULL;
104 cr->mpb->dbuf = NULL;
105 cr->mpb->ibuf_alloc = 0;
106 cr->mpb->libuf_alloc = 0;
107 cr->mpb->fbuf_alloc = 0;
108 cr->mpb->dbuf_alloc = 0;
114 t_commrec *reinitialize_commrec_for_this_thread(const t_commrec gmx_unused *cro)
116 #ifdef GMX_THREAD_MPI
119 /* make a thread-specific commrec */
121 /* now copy the whole thing, so settings like the number of PME nodes
125 /* and we start setting our own thread-specific values for things */
126 gmx_fill_commrec_from_mpi(cr);
128 // TODO cr->duty should not be initialized here
129 cr->duty = (DUTY_PP | DUTY_PME);
137 void gmx_setup_nodecomm(FILE gmx_unused *fplog, t_commrec *cr)
140 int n, rank, nodehash, ng, ni;
142 /* Many MPI implementations do not optimize MPI_Allreduce
143 * (and probably also other global communication calls)
144 * for multi-core nodes connected by a network.
145 * We can optimize such communication by using one MPI call
146 * within each node and one between the nodes.
147 * For MVAPICH2 and Intel MPI this reduces the time for
148 * the global_stat communication by 25%
149 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
150 * B. Hess, November 2007
156 #ifndef GMX_THREAD_MPI
158 MPI_Comm_size(cr->mpi_comm_mygroup, &n);
159 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
161 nodehash = gmx_physicalnode_id_hash();
165 fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
169 /* The intra-node communicator, split on node number */
170 MPI_Comm_split(cr->mpi_comm_mygroup, nodehash, rank, &nc->comm_intra);
171 MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
174 fprintf(debug, "In gmx_setup_nodecomm: node ID %d rank within node %d\n",
175 rank, nc->rank_intra);
177 /* The inter-node communicator, split on rank_intra.
178 * We actually only need the one for rank=0,
179 * but it is easier to create them all.
181 MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
182 /* Check if this really created two step communication */
183 MPI_Comm_size(nc->comm_inter, &ng);
184 MPI_Comm_size(nc->comm_intra, &ni);
187 fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
191 if (getenv("GMX_NO_NODECOMM") == NULL &&
192 ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
197 fprintf(fplog, "Using two step summing over %d groups of on average %.1f ranks\n\n",
198 ng, (real)n/(real)ng);
200 if (nc->rank_intra > 0)
202 MPI_Comm_free(&nc->comm_inter);
207 /* One group or all processes in a separate group, use normal summing */
208 MPI_Comm_free(&nc->comm_inter);
209 MPI_Comm_free(&nc->comm_intra);
212 fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
217 /* tMPI runs only on a single node so just use the nodeid */
218 nc->rank_intra = cr->nodeid;
222 void gmx_init_intranode_counters(t_commrec *cr)
224 /* counters for PP+PME and PP-only processes on my physical node */
225 int nrank_intranode, rank_intranode;
226 int nrank_pp_intranode, rank_pp_intranode;
227 /* thread-MPI is not initialized when not running in parallel */
228 #if defined GMX_MPI && !defined GMX_THREAD_MPI
229 int nrank_world, rank_world;
230 int i, myhash, *hash, *hash_s, *hash_pp, *hash_pp_s;
232 MPI_Comm_size(MPI_COMM_WORLD, &nrank_world);
233 MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
235 /* Get a (hopefully unique) hash that identifies our physical node */
236 myhash = gmx_physicalnode_id_hash();
238 /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
239 snew(hash, nrank_world);
240 snew(hash_s, nrank_world);
241 snew(hash_pp, nrank_world);
242 snew(hash_pp_s, nrank_world);
244 hash_s[rank_world] = myhash;
245 hash_pp_s[rank_world] = (cr->duty & DUTY_PP) ? myhash : -1;
247 MPI_Allreduce(hash_s, hash, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
248 MPI_Allreduce(hash_pp_s, hash_pp, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
252 nrank_pp_intranode = 0;
253 rank_pp_intranode = 0;
254 for (i = 0; i < nrank_world; i++)
256 if (hash[i] == myhash)
264 if (hash_pp[i] == myhash)
266 nrank_pp_intranode++;
267 if ((cr->duty & DUTY_PP) && i < rank_world)
278 /* Serial or thread-MPI code: we run within a single physical node */
279 nrank_intranode = cr->nnodes;
280 rank_intranode = cr->sim_nodeid;
281 nrank_pp_intranode = cr->nnodes - cr->npmenodes;
282 rank_pp_intranode = cr->nodeid;
288 if (cr->duty & DUTY_PP && cr->duty & DUTY_PME)
290 sprintf(sbuf, "PP+PME");
294 sprintf(sbuf, "%s", cr->duty & DUTY_PP ? "PP" : "PME");
296 fprintf(debug, "On %3s rank %d: nrank_intranode=%d, rank_intranode=%d, "
297 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
298 sbuf, cr->sim_nodeid,
299 nrank_intranode, rank_intranode,
300 nrank_pp_intranode, rank_pp_intranode);
303 cr->nrank_intranode = nrank_intranode;
304 cr->rank_intranode = rank_intranode;
305 cr->nrank_pp_intranode = nrank_pp_intranode;
306 cr->rank_pp_intranode = rank_pp_intranode;
310 void gmx_barrier(const t_commrec gmx_unused *cr)
313 gmx_call("gmx_barrier");
315 MPI_Barrier(cr->mpi_comm_mygroup);
319 void gmx_bcast(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
322 gmx_call("gmx_bast");
324 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
328 void gmx_bcast_sim(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
331 gmx_call("gmx_bast");
333 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
337 void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unused *cr)
340 gmx_call("gmx_sumd");
342 #if defined(MPI_IN_PLACE_EXISTS)
345 if (cr->nc.rank_intra == 0)
347 /* Use two step summing. */
348 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
350 /* Sum the roots of the internal (intra) buffers. */
351 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
356 /* This is here because of the silly MPI specification
357 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
358 MPI_Reduce(r, NULL, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
360 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
364 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
365 cr->mpi_comm_mygroup);
370 if (nr > cr->mpb->dbuf_alloc)
372 cr->mpb->dbuf_alloc = nr;
373 srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
377 /* Use two step summing */
378 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
379 if (cr->nc.rank_intra == 0)
381 /* Sum with the buffers reversed */
382 MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
385 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
389 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
390 cr->mpi_comm_mygroup);
391 for (i = 0; i < nr; i++)
393 r[i] = cr->mpb->dbuf[i];
400 void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unused *cr)
403 gmx_call("gmx_sumf");
405 #if defined(MPI_IN_PLACE_EXISTS)
408 /* Use two step summing. */
409 if (cr->nc.rank_intra == 0)
411 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
413 /* Sum the roots of the internal (intra) buffers */
414 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
419 /* This is here because of the silly MPI specification
420 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
421 MPI_Reduce(r, NULL, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
423 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
427 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
432 if (nr > cr->mpb->fbuf_alloc)
434 cr->mpb->fbuf_alloc = nr;
435 srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
439 /* Use two step summing */
440 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
441 if (cr->nc.rank_intra == 0)
443 /* Sum with the buffers reversed */
444 MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
447 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
451 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
452 cr->mpi_comm_mygroup);
453 for (i = 0; i < nr; i++)
455 r[i] = cr->mpb->fbuf[i];
462 void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused *cr)
465 gmx_call("gmx_sumi");
467 #if defined(MPI_IN_PLACE_EXISTS)
470 /* Use two step summing */
471 if (cr->nc.rank_intra == 0)
473 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
474 /* Sum with the buffers reversed */
475 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
479 /* This is here because of the silly MPI specification
480 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
481 MPI_Reduce(r, NULL, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
483 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
487 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
492 if (nr > cr->mpb->ibuf_alloc)
494 cr->mpb->ibuf_alloc = nr;
495 srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
499 /* Use two step summing */
500 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
501 if (cr->nc.rank_intra == 0)
503 /* Sum with the buffers reversed */
504 MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
506 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
510 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
511 for (i = 0; i < nr; i++)
513 r[i] = cr->mpb->ibuf[i];
520 void gmx_sumli(int gmx_unused nr, gmx_int64_t gmx_unused r[], const t_commrec gmx_unused *cr)
523 gmx_call("gmx_sumli");
525 #if defined(MPI_IN_PLACE_EXISTS)
528 /* Use two step summing */
529 if (cr->nc.rank_intra == 0)
531 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, 0,
533 /* Sum with the buffers reversed */
534 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
539 /* This is here because of the silly MPI specification
540 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
541 MPI_Reduce(r, NULL, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
543 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
547 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
552 if (nr > cr->mpb->libuf_alloc)
554 cr->mpb->libuf_alloc = nr;
555 srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
559 /* Use two step summing */
560 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
562 if (cr->nc.rank_intra == 0)
564 /* Sum with the buffers reversed */
565 MPI_Allreduce(cr->mpb->libuf, r, nr, MPI_INT64_T, MPI_SUM,
568 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
572 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
573 cr->mpi_comm_mygroup);
574 for (i = 0; i < nr; i++)
576 r[i] = cr->mpb->libuf[i];
586 static void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
588 #if defined(MPI_IN_PLACE_EXISTS)
589 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
591 /* this function is only used in code that is not performance critical,
592 (during setup, when comm_rec is not the appropriate communication
593 structure), so this isn't as bad as it looks. */
598 MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
599 for (i = 0; i < nr; i++)
609 static void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
611 #if defined(MPI_IN_PLACE_EXISTS)
612 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
614 /* this function is only used in code that is not performance critical,
615 (during setup, when comm_rec is not the appropriate communication
616 structure), so this isn't as bad as it looks. */
621 MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
622 for (i = 0; i < nr; i++)
631 void gmx_sumd_sim(int gmx_unused nr, double gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
634 gmx_call("gmx_sumd_sim");
636 gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
640 void gmx_sumf_sim(int gmx_unused nr, float gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
643 gmx_call("gmx_sumf_sim");
645 gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
649 void gmx_sumi_sim(int gmx_unused nr, int gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
652 gmx_call("gmx_sumi_sim");
654 #if defined(MPI_IN_PLACE_EXISTS)
655 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
657 /* this is thread-unsafe, but it will do for now: */
660 if (nr > ms->mpb->ibuf_alloc)
662 ms->mpb->ibuf_alloc = nr;
663 srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
665 MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
666 for (i = 0; i < nr; i++)
668 r[i] = ms->mpb->ibuf[i];
674 void gmx_sumli_sim(int gmx_unused nr, gmx_int64_t gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
677 gmx_call("gmx_sumli_sim");
679 #if defined(MPI_IN_PLACE_EXISTS)
680 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
681 ms->mpi_comm_masters);
683 /* this is thread-unsafe, but it will do for now: */
686 if (nr > ms->mpb->libuf_alloc)
688 ms->mpb->libuf_alloc = nr;
689 srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
691 MPI_Allreduce(r, ms->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
692 ms->mpi_comm_masters);
693 for (i = 0; i < nr; i++)
695 r[i] = ms->mpb->libuf[i];
701 gmx_bool gmx_fexist_master(const char *fname, t_commrec *cr)
707 bExist = gmx_fexist(fname);
711 gmx_bcast(sizeof(bExist), &bExist, cr);
716 void gmx_fatal_collective(int f_errno, const char *file, int line,
717 const t_commrec *cr, gmx_domdec_t *dd,
718 const char *fmt, ...)
721 gmx_bool bMaster, bFinalize;
724 /* Check if we are calling on all processes in MPI_COMM_WORLD */
727 MPI_Comm_compare(cr->mpi_comm_mysim, MPI_COMM_WORLD, &result);
731 MPI_Comm_compare(dd->mpi_comm_all, MPI_COMM_WORLD, &result);
733 /* Any result except MPI_UNEQUAL allows us to call MPI_Finalize */
734 bFinalize = (result != MPI_UNEQUAL);
738 bMaster = (cr != NULL && MASTER(cr)) || (dd != NULL && DDMASTER(dd));
741 gmx_fatal_mpi_va(f_errno, file, line, bMaster, bFinalize, fmt, ap);