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)
130 if (nullptr != cr->dd)
133 // done_domdec(cr->dd);
135 done_mpi_in_place_buf(cr->mpb);
140 t_commrec *reinitialize_commrec_for_this_thread(const t_commrec *cro)
145 /* make a thread-specific commrec */
147 /* now copy the whole thing, so settings like the number of PME nodes
151 /* and we start setting our own thread-specific values for things */
152 gmx_fill_commrec_from_mpi(cr);
154 // TODO cr->duty should not be initialized here
155 cr->duty = (DUTY_PP | DUTY_PME);
159 GMX_UNUSED_VALUE(cro);
164 void gmx_setup_nodecomm(FILE gmx_unused *fplog, t_commrec *cr)
168 /* Many MPI implementations do not optimize MPI_Allreduce
169 * (and probably also other global communication calls)
170 * for multi-core nodes connected by a network.
171 * We can optimize such communication by using one MPI call
172 * within each node and one between the nodes.
173 * For MVAPICH2 and Intel MPI this reduces the time for
174 * the global_stat communication by 25%
175 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
176 * B. Hess, November 2007
186 // TODO PhysicalNodeCommunicator could be extended/used to handle
187 // the need for per-node per-group communicators.
188 MPI_Comm_size(cr->mpi_comm_mygroup, &n);
189 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
191 int nodehash = gmx_physicalnode_id_hash();
195 fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
199 /* The intra-node communicator, split on node number */
200 MPI_Comm_split(cr->mpi_comm_mygroup, nodehash, rank, &nc->comm_intra);
201 MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
204 fprintf(debug, "In gmx_setup_nodecomm: node ID %d rank within node %d\n",
205 rank, nc->rank_intra);
207 /* The inter-node communicator, split on rank_intra.
208 * We actually only need the one for rank=0,
209 * but it is easier to create them all.
211 MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
212 /* Check if this really created two step communication */
215 MPI_Comm_size(nc->comm_inter, &ng);
216 MPI_Comm_size(nc->comm_intra, &ni);
219 fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
223 if (getenv("GMX_NO_NODECOMM") == nullptr &&
224 ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
229 fprintf(fplog, "Using two step summing over %d groups of on average %.1f ranks\n\n",
230 ng, (real)n/(real)ng);
232 if (nc->rank_intra > 0)
234 MPI_Comm_free(&nc->comm_inter);
239 /* One group or all processes in a separate group, use normal summing */
240 MPI_Comm_free(&nc->comm_inter);
241 MPI_Comm_free(&nc->comm_intra);
244 fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
249 /* tMPI runs only on a single node so just use the nodeid */
250 nc->rank_intra = cr->nodeid;
254 void gmx_barrier(const t_commrec gmx_unused *cr)
257 gmx_call("gmx_barrier");
259 MPI_Barrier(cr->mpi_comm_mygroup);
263 void gmx_bcast(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
266 gmx_call("gmx_bast");
268 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
272 void gmx_bcast_sim(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
275 gmx_call("gmx_bast");
277 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
281 void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unused *cr)
284 gmx_call("gmx_sumd");
286 #if MPI_IN_PLACE_EXISTS
289 if (cr->nc.rank_intra == 0)
291 /* Use two step summing. */
292 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
294 /* Sum the roots of the internal (intra) buffers. */
295 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
300 /* This is here because of the silly MPI specification
301 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
302 MPI_Reduce(r, nullptr, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
304 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
308 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
309 cr->mpi_comm_mygroup);
314 if (nr > cr->mpb->dbuf_alloc)
316 cr->mpb->dbuf_alloc = nr;
317 srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
321 /* Use two step summing */
322 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
323 if (cr->nc.rank_intra == 0)
325 /* Sum with the buffers reversed */
326 MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
329 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
333 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
334 cr->mpi_comm_mygroup);
335 for (i = 0; i < nr; i++)
337 r[i] = cr->mpb->dbuf[i];
344 void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unused *cr)
347 gmx_call("gmx_sumf");
349 #if MPI_IN_PLACE_EXISTS
352 /* Use two step summing. */
353 if (cr->nc.rank_intra == 0)
355 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
357 /* Sum the roots of the internal (intra) buffers */
358 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
363 /* This is here because of the silly MPI specification
364 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
365 MPI_Reduce(r, nullptr, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
367 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
371 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
376 if (nr > cr->mpb->fbuf_alloc)
378 cr->mpb->fbuf_alloc = nr;
379 srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
383 /* Use two step summing */
384 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
385 if (cr->nc.rank_intra == 0)
387 /* Sum with the buffers reversed */
388 MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
391 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
395 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
396 cr->mpi_comm_mygroup);
397 for (i = 0; i < nr; i++)
399 r[i] = cr->mpb->fbuf[i];
406 void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused *cr)
409 gmx_call("gmx_sumi");
411 #if MPI_IN_PLACE_EXISTS
414 /* Use two step summing */
415 if (cr->nc.rank_intra == 0)
417 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
418 /* Sum with the buffers reversed */
419 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
423 /* This is here because of the silly MPI specification
424 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
425 MPI_Reduce(r, nullptr, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
427 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
431 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
436 if (nr > cr->mpb->ibuf_alloc)
438 cr->mpb->ibuf_alloc = nr;
439 srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
443 /* Use two step summing */
444 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
445 if (cr->nc.rank_intra == 0)
447 /* Sum with the buffers reversed */
448 MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
450 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
454 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
455 for (i = 0; i < nr; i++)
457 r[i] = cr->mpb->ibuf[i];
464 void gmx_sumli(int gmx_unused nr, int64_t gmx_unused r[], const t_commrec gmx_unused *cr)
467 gmx_call("gmx_sumli");
469 #if MPI_IN_PLACE_EXISTS
472 /* Use two step summing */
473 if (cr->nc.rank_intra == 0)
475 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, 0,
477 /* Sum with the buffers reversed */
478 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
483 /* This is here because of the silly MPI specification
484 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
485 MPI_Reduce(r, nullptr, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
487 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
491 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
496 if (nr > cr->mpb->libuf_alloc)
498 cr->mpb->libuf_alloc = nr;
499 srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
503 /* Use two step summing */
504 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
506 if (cr->nc.rank_intra == 0)
508 /* Sum with the buffers reversed */
509 MPI_Allreduce(cr->mpb->libuf, r, nr, MPI_INT64_T, MPI_SUM,
512 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
516 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
517 cr->mpi_comm_mygroup);
518 for (i = 0; i < nr; i++)
520 r[i] = cr->mpb->libuf[i];
530 static void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
532 #if MPI_IN_PLACE_EXISTS
533 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
535 /* this function is only used in code that is not performance critical,
536 (during setup, when comm_rec is not the appropriate communication
537 structure), so this isn't as bad as it looks. */
542 MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
543 for (i = 0; i < nr; i++)
553 static void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
555 #if MPI_IN_PLACE_EXISTS
556 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
558 /* this function is only used in code that is not performance critical,
559 (during setup, when comm_rec is not the appropriate communication
560 structure), so this isn't as bad as it looks. */
565 MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
566 for (i = 0; i < nr; i++)
575 void gmx_sumd_sim(int gmx_unused nr, double gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
578 gmx_call("gmx_sumd_sim");
580 gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
584 void gmx_sumf_sim(int gmx_unused nr, float gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
587 gmx_call("gmx_sumf_sim");
589 gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
593 void gmx_sumi_sim(int gmx_unused nr, int gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
596 gmx_call("gmx_sumi_sim");
598 #if MPI_IN_PLACE_EXISTS
599 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
601 /* this is thread-unsafe, but it will do for now: */
604 if (nr > ms->mpb->ibuf_alloc)
606 ms->mpb->ibuf_alloc = nr;
607 srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
609 MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
610 for (i = 0; i < nr; i++)
612 r[i] = ms->mpb->ibuf[i];
618 void gmx_sumli_sim(int gmx_unused nr, int64_t gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
621 gmx_call("gmx_sumli_sim");
623 #if MPI_IN_PLACE_EXISTS
624 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
625 ms->mpi_comm_masters);
627 /* this is thread-unsafe, but it will do for now: */
630 if (nr > ms->mpb->libuf_alloc)
632 ms->mpb->libuf_alloc = nr;
633 srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
635 MPI_Allreduce(r, ms->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
636 ms->mpi_comm_masters);
637 for (i = 0; i < nr; i++)
639 r[i] = ms->mpb->libuf[i];
645 void check_multi_int(FILE *log, const gmx_multisim_t *ms, int val,
650 gmx_bool bCompatible;
652 if (nullptr != log && !bQuiet)
654 fprintf(log, "Multi-checking %s ... ", name);
660 "check_multi_int called with a NULL communication pointer");
663 snew(ibuf, ms->nsim);
665 gmx_sumi_sim(ms->nsim, ibuf, ms);
668 for (p = 1; p < ms->nsim; p++)
670 bCompatible = bCompatible && (ibuf[p-1] == ibuf[p]);
675 if (nullptr != log && !bQuiet)
677 fprintf(log, "OK\n");
684 fprintf(log, "\n%s is not equal for all subsystems\n", name);
685 for (p = 0; p < ms->nsim; p++)
687 fprintf(log, " subsystem %d: %d\n", p, ibuf[p]);
690 gmx_fatal(FARGS, "The %d subsystems are not compatible\n", ms->nsim);
696 void check_multi_int64(FILE *log, const gmx_multisim_t *ms,
697 int64_t val, const char *name,
702 gmx_bool bCompatible;
704 if (nullptr != log && !bQuiet)
706 fprintf(log, "Multi-checking %s ... ", name);
712 "check_multi_int called with a NULL communication pointer");
715 snew(ibuf, ms->nsim);
717 gmx_sumli_sim(ms->nsim, ibuf, ms);
720 for (p = 1; p < ms->nsim; p++)
722 bCompatible = bCompatible && (ibuf[p-1] == ibuf[p]);
727 if (nullptr != log && !bQuiet)
729 fprintf(log, "OK\n");
734 // TODO Part of this error message would also be good to go to
735 // stderr (from one rank of one sim only)
738 fprintf(log, "\n%s is not equal for all subsystems\n", name);
739 for (p = 0; p < ms->nsim; p++)
742 /* first make the format string */
743 snprintf(strbuf, 255, " subsystem %%d: %s\n",
745 fprintf(log, strbuf, p, ibuf[p]);
748 gmx_fatal(FARGS, "The %d subsystems are not compatible\n", ms->nsim);
754 const char *opt2fn_master(const char *opt, int nfile, const t_filenm fnm[],
757 return SIMMASTER(cr) ? opt2fn(opt, nfile, fnm) : nullptr;
760 void gmx_fatal_collective(int f_errno, const char *file, int line,
761 MPI_Comm comm, gmx_bool bMaster,
762 gmx_fmtstr const char *fmt, ...)
768 /* Check if we are calling on all processes in MPI_COMM_WORLD */
769 MPI_Comm_compare(comm, MPI_COMM_WORLD, &result);
770 /* Any result except MPI_UNEQUAL allows us to call MPI_Finalize */
771 bFinalize = (result != MPI_UNEQUAL);
773 GMX_UNUSED_VALUE(comm);
778 gmx_fatal_mpi_va(f_errno, file, line, bMaster, bFinalize, fmt, ap);