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,2015,2016,2017,2018, 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.
48 #include "gromacs/commandline/filenm.h"
49 #include "gromacs/mdtypes/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 *cr)
63 gmx_call("gmx_fill_commrec_from_mpi");
66 if (!gmx_mpi_initialized())
68 gmx_comm("MPI has not been initialized properly");
71 cr->nnodes = gmx_node_num();
72 cr->nodeid = gmx_node_rank();
73 cr->sim_nodeid = cr->nodeid;
74 cr->mpi_comm_mysim = MPI_COMM_WORLD;
75 cr->mpi_comm_mygroup = MPI_COMM_WORLD;
79 t_commrec *init_commrec()
86 gmx_fill_commrec_from_mpi(cr);
88 cr->mpi_comm_mysim = MPI_COMM_NULL;
89 cr->mpi_comm_mygroup = MPI_COMM_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 GMX_MPI && !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 void done_mpi_in_place_buf(mpi_in_place_buf_t *buf)
126 void done_commrec(t_commrec *cr)
128 if (nullptr != cr->dd)
131 // done_domdec(cr->dd);
133 done_mpi_in_place_buf(cr->mpb);
137 t_commrec *reinitialize_commrec_for_this_thread(const t_commrec *cro)
142 /* make a thread-specific commrec */
144 /* now copy the whole thing, so settings like the number of PME nodes
148 /* and we start setting our own thread-specific values for things */
149 gmx_fill_commrec_from_mpi(cr);
151 // TODO cr->duty should not be initialized here
152 cr->duty = (DUTY_PP | DUTY_PME);
156 GMX_UNUSED_VALUE(cro);
161 void gmx_setup_nodecomm(FILE gmx_unused *fplog, t_commrec *cr)
165 /* Many MPI implementations do not optimize MPI_Allreduce
166 * (and probably also other global communication calls)
167 * for multi-core nodes connected by a network.
168 * We can optimize such communication by using one MPI call
169 * within each node and one between the nodes.
170 * For MVAPICH2 and Intel MPI this reduces the time for
171 * the global_stat communication by 25%
172 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
173 * B. Hess, November 2007
183 // TODO PhysicalNodeCommunicator could be extended/used to handle
184 // the need for per-node per-group communicators.
185 MPI_Comm_size(cr->mpi_comm_mygroup, &n);
186 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
188 int nodehash = gmx_physicalnode_id_hash();
192 fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
196 /* The intra-node communicator, split on node number */
197 MPI_Comm_split(cr->mpi_comm_mygroup, nodehash, rank, &nc->comm_intra);
198 MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
201 fprintf(debug, "In gmx_setup_nodecomm: node ID %d rank within node %d\n",
202 rank, nc->rank_intra);
204 /* The inter-node communicator, split on rank_intra.
205 * We actually only need the one for rank=0,
206 * but it is easier to create them all.
208 MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
209 /* Check if this really created two step communication */
212 MPI_Comm_size(nc->comm_inter, &ng);
213 MPI_Comm_size(nc->comm_intra, &ni);
216 fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
220 if (getenv("GMX_NO_NODECOMM") == nullptr &&
221 ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
226 fprintf(fplog, "Using two step summing over %d groups of on average %.1f ranks\n\n",
227 ng, (real)n/(real)ng);
229 if (nc->rank_intra > 0)
231 MPI_Comm_free(&nc->comm_inter);
236 /* One group or all processes in a separate group, use normal summing */
237 MPI_Comm_free(&nc->comm_inter);
238 MPI_Comm_free(&nc->comm_intra);
241 fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
246 /* tMPI runs only on a single node so just use the nodeid */
247 nc->rank_intra = cr->nodeid;
251 void gmx_barrier(const t_commrec gmx_unused *cr)
254 gmx_call("gmx_barrier");
256 MPI_Barrier(cr->mpi_comm_mygroup);
260 void gmx_bcast(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
263 gmx_call("gmx_bast");
265 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
269 void gmx_bcast_sim(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
272 gmx_call("gmx_bast");
274 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
278 void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unused *cr)
281 gmx_call("gmx_sumd");
283 #if MPI_IN_PLACE_EXISTS
286 if (cr->nc.rank_intra == 0)
288 /* Use two step summing. */
289 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
291 /* Sum the roots of the internal (intra) buffers. */
292 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
297 /* This is here because of the silly MPI specification
298 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
299 MPI_Reduce(r, nullptr, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
301 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
305 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
306 cr->mpi_comm_mygroup);
311 if (nr > cr->mpb->dbuf_alloc)
313 cr->mpb->dbuf_alloc = nr;
314 srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
318 /* Use two step summing */
319 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
320 if (cr->nc.rank_intra == 0)
322 /* Sum with the buffers reversed */
323 MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
326 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
330 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
331 cr->mpi_comm_mygroup);
332 for (i = 0; i < nr; i++)
334 r[i] = cr->mpb->dbuf[i];
341 void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unused *cr)
344 gmx_call("gmx_sumf");
346 #if MPI_IN_PLACE_EXISTS
349 /* Use two step summing. */
350 if (cr->nc.rank_intra == 0)
352 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
354 /* Sum the roots of the internal (intra) buffers */
355 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
360 /* This is here because of the silly MPI specification
361 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
362 MPI_Reduce(r, nullptr, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
364 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
368 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
373 if (nr > cr->mpb->fbuf_alloc)
375 cr->mpb->fbuf_alloc = nr;
376 srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
380 /* Use two step summing */
381 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
382 if (cr->nc.rank_intra == 0)
384 /* Sum with the buffers reversed */
385 MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
388 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
392 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
393 cr->mpi_comm_mygroup);
394 for (i = 0; i < nr; i++)
396 r[i] = cr->mpb->fbuf[i];
403 void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused *cr)
406 gmx_call("gmx_sumi");
408 #if MPI_IN_PLACE_EXISTS
411 /* Use two step summing */
412 if (cr->nc.rank_intra == 0)
414 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
415 /* Sum with the buffers reversed */
416 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
420 /* This is here because of the silly MPI specification
421 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
422 MPI_Reduce(r, nullptr, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
424 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
428 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
433 if (nr > cr->mpb->ibuf_alloc)
435 cr->mpb->ibuf_alloc = nr;
436 srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
440 /* Use two step summing */
441 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
442 if (cr->nc.rank_intra == 0)
444 /* Sum with the buffers reversed */
445 MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
447 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
451 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
452 for (i = 0; i < nr; i++)
454 r[i] = cr->mpb->ibuf[i];
461 void gmx_sumli(int gmx_unused nr, int64_t gmx_unused r[], const t_commrec gmx_unused *cr)
464 gmx_call("gmx_sumli");
466 #if 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_INT64_T, MPI_SUM, 0,
474 /* Sum with the buffers reversed */
475 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
480 /* This is here because of the silly MPI specification
481 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
482 MPI_Reduce(r, nullptr, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
484 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
488 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
493 if (nr > cr->mpb->libuf_alloc)
495 cr->mpb->libuf_alloc = nr;
496 srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
500 /* Use two step summing */
501 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
503 if (cr->nc.rank_intra == 0)
505 /* Sum with the buffers reversed */
506 MPI_Allreduce(cr->mpb->libuf, r, nr, MPI_INT64_T, MPI_SUM,
509 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
513 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
514 cr->mpi_comm_mygroup);
515 for (i = 0; i < nr; i++)
517 r[i] = cr->mpb->libuf[i];
527 static void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
529 #if MPI_IN_PLACE_EXISTS
530 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
532 /* this function is only used in code that is not performance critical,
533 (during setup, when comm_rec is not the appropriate communication
534 structure), so this isn't as bad as it looks. */
539 MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
540 for (i = 0; i < nr; i++)
550 static void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
552 #if MPI_IN_PLACE_EXISTS
553 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
555 /* this function is only used in code that is not performance critical,
556 (during setup, when comm_rec is not the appropriate communication
557 structure), so this isn't as bad as it looks. */
562 MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
563 for (i = 0; i < nr; i++)
572 void gmx_sumd_sim(int gmx_unused nr, double gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
575 gmx_call("gmx_sumd_sim");
577 gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
581 void gmx_sumf_sim(int gmx_unused nr, float gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
584 gmx_call("gmx_sumf_sim");
586 gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
590 void gmx_sumi_sim(int gmx_unused nr, int gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
593 gmx_call("gmx_sumi_sim");
595 #if MPI_IN_PLACE_EXISTS
596 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
598 /* this is thread-unsafe, but it will do for now: */
601 if (nr > ms->mpb->ibuf_alloc)
603 ms->mpb->ibuf_alloc = nr;
604 srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
606 MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
607 for (i = 0; i < nr; i++)
609 r[i] = ms->mpb->ibuf[i];
615 void gmx_sumli_sim(int gmx_unused nr, int64_t gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
618 gmx_call("gmx_sumli_sim");
620 #if MPI_IN_PLACE_EXISTS
621 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
622 ms->mpi_comm_masters);
624 /* this is thread-unsafe, but it will do for now: */
627 if (nr > ms->mpb->libuf_alloc)
629 ms->mpb->libuf_alloc = nr;
630 srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
632 MPI_Allreduce(r, ms->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
633 ms->mpi_comm_masters);
634 for (i = 0; i < nr; i++)
636 r[i] = ms->mpb->libuf[i];
642 void check_multi_int(FILE *log, const gmx_multisim_t *ms, int val,
647 gmx_bool bCompatible;
649 if (nullptr != log && !bQuiet)
651 fprintf(log, "Multi-checking %s ... ", name);
657 "check_multi_int called with a NULL communication pointer");
660 snew(ibuf, ms->nsim);
662 gmx_sumi_sim(ms->nsim, ibuf, ms);
665 for (p = 1; p < ms->nsim; p++)
667 bCompatible = bCompatible && (ibuf[p-1] == ibuf[p]);
672 if (nullptr != log && !bQuiet)
674 fprintf(log, "OK\n");
681 fprintf(log, "\n%s is not equal for all subsystems\n", name);
682 for (p = 0; p < ms->nsim; p++)
684 fprintf(log, " subsystem %d: %d\n", p, ibuf[p]);
687 gmx_fatal(FARGS, "The %d subsystems are not compatible\n", ms->nsim);
693 void check_multi_int64(FILE *log, const gmx_multisim_t *ms,
694 int64_t val, const char *name,
699 gmx_bool bCompatible;
701 if (nullptr != log && !bQuiet)
703 fprintf(log, "Multi-checking %s ... ", name);
709 "check_multi_int called with a NULL communication pointer");
712 snew(ibuf, ms->nsim);
714 gmx_sumli_sim(ms->nsim, ibuf, ms);
717 for (p = 1; p < ms->nsim; p++)
719 bCompatible = bCompatible && (ibuf[p-1] == ibuf[p]);
724 if (nullptr != log && !bQuiet)
726 fprintf(log, "OK\n");
731 // TODO Part of this error message would also be good to go to
732 // stderr (from one rank of one sim only)
735 fprintf(log, "\n%s is not equal for all subsystems\n", name);
736 for (p = 0; p < ms->nsim; p++)
739 /* first make the format string */
740 snprintf(strbuf, 255, " subsystem %%d: %s\n",
742 fprintf(log, strbuf, p, ibuf[p]);
745 gmx_fatal(FARGS, "The %d subsystems are not compatible\n", ms->nsim);
751 const char *opt2fn_master(const char *opt, int nfile, const t_filenm fnm[],
754 return SIMMASTER(cr) ? opt2fn(opt, nfile, fnm) : nullptr;
757 void gmx_fatal_collective(int f_errno, const char *file, int line,
758 MPI_Comm comm, gmx_bool bMaster,
759 gmx_fmtstr const char *fmt, ...)
765 /* Check if we are calling on all processes in MPI_COMM_WORLD */
766 MPI_Comm_compare(comm, MPI_COMM_WORLD, &result);
767 /* Any result except MPI_UNEQUAL allows us to call MPI_Finalize */
768 bFinalize = (result != MPI_UNEQUAL);
770 GMX_UNUSED_VALUE(comm);
775 gmx_fatal_mpi_va(f_errno, file, line, bMaster, bFinalize, fmt, ap);