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, 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.
39 #include "gromacs/legacyheaders/network.h"
48 #include "gromacs/legacyheaders/copyrite.h"
49 #include "gromacs/legacyheaders/macros.h"
50 #include "gromacs/legacyheaders/types/commrec.h"
51 #include "gromacs/utility/basenetwork.h"
52 #include "gromacs/utility/cstringutil.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/futil.h"
55 #include "gromacs/utility/gmxmpi.h"
56 #include "gromacs/utility/smalloc.h"
58 /* The source code in this file should be thread-safe.
59 Please keep it that way. */
61 void gmx_fill_commrec_from_mpi(t_commrec gmx_unused *cr)
64 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;
80 t_commrec *init_commrec()
87 gmx_fill_commrec_from_mpi(cr);
89 cr->mpi_comm_mysim = NULL;
90 cr->mpi_comm_mygroup = NULL;
93 cr->nodeid = cr->sim_nodeid;
96 // TODO cr->duty should not be initialized here
97 cr->duty = (DUTY_PP | DUTY_PME);
99 #if defined GMX_MPI && !defined MPI_IN_PLACE_EXISTS
100 /* initialize the MPI_IN_PLACE replacement buffers */
102 cr->mpb->ibuf = NULL;
103 cr->mpb->libuf = NULL;
104 cr->mpb->fbuf = NULL;
105 cr->mpb->dbuf = NULL;
106 cr->mpb->ibuf_alloc = 0;
107 cr->mpb->libuf_alloc = 0;
108 cr->mpb->fbuf_alloc = 0;
109 cr->mpb->dbuf_alloc = 0;
115 t_commrec *reinitialize_commrec_for_this_thread(const t_commrec gmx_unused *cro)
117 #ifdef GMX_THREAD_MPI
120 /* make a thread-specific commrec */
122 /* now copy the whole thing, so settings like the number of PME nodes
126 /* and we start setting our own thread-specific values for things */
127 gmx_fill_commrec_from_mpi(cr);
129 // TODO cr->duty should not be initialized here
130 cr->duty = (DUTY_PP | DUTY_PME);
138 void gmx_setup_nodecomm(FILE gmx_unused *fplog, t_commrec *cr)
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
160 MPI_Comm_size(cr->mpi_comm_mygroup, &n);
161 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
163 int nodehash = gmx_physicalnode_id_hash();
167 fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
171 /* The intra-node communicator, split on node number */
172 MPI_Comm_split(cr->mpi_comm_mygroup, nodehash, rank, &nc->comm_intra);
173 MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
176 fprintf(debug, "In gmx_setup_nodecomm: node ID %d rank within node %d\n",
177 rank, nc->rank_intra);
179 /* The inter-node communicator, split on rank_intra.
180 * We actually only need the one for rank=0,
181 * but it is easier to create them all.
183 MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
184 /* Check if this really created two step communication */
187 MPI_Comm_size(nc->comm_inter, &ng);
188 MPI_Comm_size(nc->comm_intra, &ni);
191 fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
195 if (getenv("GMX_NO_NODECOMM") == NULL &&
196 ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
201 fprintf(fplog, "Using two step summing over %d groups of on average %.1f ranks\n\n",
202 ng, (real)n/(real)ng);
204 if (nc->rank_intra > 0)
206 MPI_Comm_free(&nc->comm_inter);
211 /* One group or all processes in a separate group, use normal summing */
212 MPI_Comm_free(&nc->comm_inter);
213 MPI_Comm_free(&nc->comm_intra);
216 fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
221 /* tMPI runs only on a single node so just use the nodeid */
222 nc->rank_intra = cr->nodeid;
226 void gmx_init_intranode_counters(t_commrec *cr)
228 /* counters for PP+PME and PP-only processes on my physical node */
229 int nrank_intranode, rank_intranode;
230 int nrank_pp_intranode, rank_pp_intranode;
231 /* thread-MPI is not initialized when not running in parallel */
232 #if defined GMX_MPI && !defined GMX_THREAD_MPI
233 int nrank_world, rank_world;
234 int i, myhash, *hash, *hash_s, *hash_pp, *hash_pp_s;
236 MPI_Comm_size(MPI_COMM_WORLD, &nrank_world);
237 MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
239 /* Get a (hopefully unique) hash that identifies our physical node */
240 myhash = gmx_physicalnode_id_hash();
242 /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
243 snew(hash, nrank_world);
244 snew(hash_s, nrank_world);
245 snew(hash_pp, nrank_world);
246 snew(hash_pp_s, nrank_world);
248 hash_s[rank_world] = myhash;
249 hash_pp_s[rank_world] = (cr->duty & DUTY_PP) ? myhash : -1;
251 MPI_Allreduce(hash_s, hash, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
252 MPI_Allreduce(hash_pp_s, hash_pp, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
256 nrank_pp_intranode = 0;
257 rank_pp_intranode = 0;
258 for (i = 0; i < nrank_world; i++)
260 if (hash[i] == myhash)
268 if (hash_pp[i] == myhash)
270 nrank_pp_intranode++;
271 if ((cr->duty & DUTY_PP) && i < rank_world)
282 /* Serial or thread-MPI code: we run within a single physical node */
283 nrank_intranode = cr->nnodes;
284 rank_intranode = cr->sim_nodeid;
285 nrank_pp_intranode = cr->nnodes - cr->npmenodes;
286 rank_pp_intranode = cr->nodeid;
292 if (cr->duty & DUTY_PP && cr->duty & DUTY_PME)
294 sprintf(sbuf, "PP+PME");
298 sprintf(sbuf, "%s", cr->duty & DUTY_PP ? "PP" : "PME");
300 fprintf(debug, "On %3s rank %d: nrank_intranode=%d, rank_intranode=%d, "
301 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
302 sbuf, cr->sim_nodeid,
303 nrank_intranode, rank_intranode,
304 nrank_pp_intranode, rank_pp_intranode);
307 cr->nrank_intranode = nrank_intranode;
308 cr->rank_intranode = rank_intranode;
309 cr->nrank_pp_intranode = nrank_pp_intranode;
310 cr->rank_pp_intranode = rank_pp_intranode;
314 void gmx_barrier(const t_commrec gmx_unused *cr)
317 gmx_call("gmx_barrier");
319 MPI_Barrier(cr->mpi_comm_mygroup);
323 void gmx_bcast(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
326 gmx_call("gmx_bast");
328 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
332 void gmx_bcast_sim(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
335 gmx_call("gmx_bast");
337 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
341 void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unused *cr)
344 gmx_call("gmx_sumd");
346 #if defined(MPI_IN_PLACE_EXISTS)
349 if (cr->nc.rank_intra == 0)
351 /* Use two step summing. */
352 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
354 /* Sum the roots of the internal (intra) buffers. */
355 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, 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, NULL, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
364 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
368 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
369 cr->mpi_comm_mygroup);
374 if (nr > cr->mpb->dbuf_alloc)
376 cr->mpb->dbuf_alloc = nr;
377 srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
381 /* Use two step summing */
382 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
383 if (cr->nc.rank_intra == 0)
385 /* Sum with the buffers reversed */
386 MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
389 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
393 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
394 cr->mpi_comm_mygroup);
395 for (i = 0; i < nr; i++)
397 r[i] = cr->mpb->dbuf[i];
404 void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unused *cr)
407 gmx_call("gmx_sumf");
409 #if defined(MPI_IN_PLACE_EXISTS)
412 /* Use two step summing. */
413 if (cr->nc.rank_intra == 0)
415 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
417 /* Sum the roots of the internal (intra) buffers */
418 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
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, NULL, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
427 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
431 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
436 if (nr > cr->mpb->fbuf_alloc)
438 cr->mpb->fbuf_alloc = nr;
439 srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
443 /* Use two step summing */
444 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
445 if (cr->nc.rank_intra == 0)
447 /* Sum with the buffers reversed */
448 MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
451 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
455 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
456 cr->mpi_comm_mygroup);
457 for (i = 0; i < nr; i++)
459 r[i] = cr->mpb->fbuf[i];
466 void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused *cr)
469 gmx_call("gmx_sumi");
471 #if defined(MPI_IN_PLACE_EXISTS)
474 /* Use two step summing */
475 if (cr->nc.rank_intra == 0)
477 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
478 /* Sum with the buffers reversed */
479 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
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, NULL, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
487 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
491 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
496 if (nr > cr->mpb->ibuf_alloc)
498 cr->mpb->ibuf_alloc = nr;
499 srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
503 /* Use two step summing */
504 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
505 if (cr->nc.rank_intra == 0)
507 /* Sum with the buffers reversed */
508 MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
510 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
514 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
515 for (i = 0; i < nr; i++)
517 r[i] = cr->mpb->ibuf[i];
524 void gmx_sumli(int gmx_unused nr, gmx_int64_t gmx_unused r[], const t_commrec gmx_unused *cr)
527 gmx_call("gmx_sumli");
529 #if defined(MPI_IN_PLACE_EXISTS)
532 /* Use two step summing */
533 if (cr->nc.rank_intra == 0)
535 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, 0,
537 /* Sum with the buffers reversed */
538 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
543 /* This is here because of the silly MPI specification
544 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
545 MPI_Reduce(r, NULL, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
547 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
551 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
556 if (nr > cr->mpb->libuf_alloc)
558 cr->mpb->libuf_alloc = nr;
559 srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
563 /* Use two step summing */
564 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
566 if (cr->nc.rank_intra == 0)
568 /* Sum with the buffers reversed */
569 MPI_Allreduce(cr->mpb->libuf, r, nr, MPI_INT64_T, MPI_SUM,
572 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
576 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
577 cr->mpi_comm_mygroup);
578 for (i = 0; i < nr; i++)
580 r[i] = cr->mpb->libuf[i];
590 static void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
592 #if defined(MPI_IN_PLACE_EXISTS)
593 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
595 /* this function is only used in code that is not performance critical,
596 (during setup, when comm_rec is not the appropriate communication
597 structure), so this isn't as bad as it looks. */
602 MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
603 for (i = 0; i < nr; i++)
613 static void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
615 #if defined(MPI_IN_PLACE_EXISTS)
616 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
618 /* this function is only used in code that is not performance critical,
619 (during setup, when comm_rec is not the appropriate communication
620 structure), so this isn't as bad as it looks. */
625 MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
626 for (i = 0; i < nr; i++)
635 void gmx_sumd_sim(int gmx_unused nr, double gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
638 gmx_call("gmx_sumd_sim");
640 gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
644 void gmx_sumf_sim(int gmx_unused nr, float gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
647 gmx_call("gmx_sumf_sim");
649 gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
653 void gmx_sumi_sim(int gmx_unused nr, int gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
656 gmx_call("gmx_sumi_sim");
658 #if defined(MPI_IN_PLACE_EXISTS)
659 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
661 /* this is thread-unsafe, but it will do for now: */
664 if (nr > ms->mpb->ibuf_alloc)
666 ms->mpb->ibuf_alloc = nr;
667 srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
669 MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
670 for (i = 0; i < nr; i++)
672 r[i] = ms->mpb->ibuf[i];
678 void gmx_sumli_sim(int gmx_unused nr, gmx_int64_t gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
681 gmx_call("gmx_sumli_sim");
683 #if defined(MPI_IN_PLACE_EXISTS)
684 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
685 ms->mpi_comm_masters);
687 /* this is thread-unsafe, but it will do for now: */
690 if (nr > ms->mpb->libuf_alloc)
692 ms->mpb->libuf_alloc = nr;
693 srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
695 MPI_Allreduce(r, ms->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
696 ms->mpi_comm_masters);
697 for (i = 0; i < nr; i++)
699 r[i] = ms->mpb->libuf[i];
705 gmx_bool gmx_fexist_master(const char *fname, t_commrec *cr)
711 bExist = gmx_fexist(fname);
715 gmx_bcast(sizeof(bExist), &bExist, cr);
720 void gmx_fatal_collective(int f_errno, const char *file, int line,
721 const t_commrec *cr, gmx_domdec_t *dd,
722 const char *fmt, ...)
725 gmx_bool bMaster, bFinalize;
728 /* Check if we are calling on all processes in MPI_COMM_WORLD */
731 MPI_Comm_compare(cr->mpi_comm_mysim, MPI_COMM_WORLD, &result);
735 MPI_Comm_compare(dd->mpi_comm_all, MPI_COMM_WORLD, &result);
737 /* Any result except MPI_UNEQUAL allows us to call MPI_Finalize */
738 bFinalize = (result != MPI_UNEQUAL);
742 bMaster = (cr != NULL && MASTER(cr)) || (dd != NULL && DDMASTER(dd));
745 gmx_fatal_mpi_va(f_errno, file, line, bMaster, bFinalize, fmt, ap);