1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * GROningen Mixture of Alchemy and Childrens' Stories
40 #include "gmx_fatal.h"
49 #include "gromacs/utility/gmxmpi.h"
52 /* The source code in this file should be thread-safe.
53 Please keep it that way. */
55 gmx_bool gmx_mpi_initialized(void)
67 void gmx_fill_commrec_from_mpi(t_commrec *cr)
70 gmx_call("gmx_fill_commrec_from_mpi");
72 if (!gmx_mpi_initialized())
74 gmx_comm("MPI has not been initialized properly");
77 cr->nnodes = gmx_node_num();
78 cr->nodeid = gmx_node_rank();
79 cr->sim_nodeid = cr->nodeid;
80 cr->mpi_comm_mysim = MPI_COMM_WORLD;
81 cr->mpi_comm_mygroup = MPI_COMM_WORLD;
86 t_commrec *init_commrec()
93 gmx_fill_commrec_from_mpi(cr);
95 cr->mpi_comm_mysim = NULL;
96 cr->mpi_comm_mygroup = NULL;
99 cr->nodeid = cr->sim_nodeid;
102 // TODO cr->duty should not be initialized here
103 cr->duty = (DUTY_PP | DUTY_PME);
105 #if defined GMX_LIB_MPI && !defined MPI_IN_PLACE_EXISTS
106 /* initialize the MPI_IN_PLACE replacement buffers */
108 cr->mpb->ibuf = NULL;
109 cr->mpb->libuf = NULL;
110 cr->mpb->fbuf = NULL;
111 cr->mpb->dbuf = NULL;
112 cr->mpb->ibuf_alloc = 0;
113 cr->mpb->libuf_alloc = 0;
114 cr->mpb->fbuf_alloc = 0;
115 cr->mpb->dbuf_alloc = 0;
121 t_commrec *reinitialize_commrec_for_this_thread(const t_commrec *cro)
123 #ifdef GMX_THREAD_MPI
126 /* make a thread-specific commrec */
128 /* now copy the whole thing, so settings like the number of PME nodes
132 /* and we start setting our own thread-specific values for things */
133 gmx_fill_commrec_from_mpi(cr);
135 // TODO cr->duty should not be initialized here
136 cr->duty = (DUTY_PP | DUTY_PME);
144 int gmx_node_num(void)
150 (void) MPI_Comm_size(MPI_COMM_WORLD, &i);
155 int gmx_node_rank(void)
161 (void) MPI_Comm_rank(MPI_COMM_WORLD, &i);
167 int gmx_hostname_num()
172 #ifdef GMX_THREAD_MPI
173 /* thread-MPI currently puts the thread number in the process name,
174 * we might want to change this, as this is inconsistent with what
175 * most MPI implementations would do when running on a single node.
179 int resultlen, hostnum, i, j;
180 char mpi_hostname[MPI_MAX_PROCESSOR_NAME], hostnum_str[MPI_MAX_PROCESSOR_NAME];
182 MPI_Get_processor_name(mpi_hostname, &resultlen);
183 /* This procedure can only differentiate nodes with host names
184 * that end on unique numbers.
188 /* Only parse the host name up to the first dot */
189 while (i < resultlen && mpi_hostname[i] != '.')
191 if (isdigit(mpi_hostname[i]))
193 hostnum_str[j++] = mpi_hostname[i];
197 hostnum_str[j] = '\0';
204 /* Use only the last 9 decimals, so we don't overflow an int */
205 hostnum = strtol(hostnum_str + max(0, j-9), NULL, 10);
210 fprintf(debug, "In gmx_setup_nodecomm: hostname '%s', hostnum %d\n",
211 mpi_hostname, hostnum);
218 void gmx_setup_nodecomm(FILE gmx_unused *fplog, t_commrec *cr)
221 int n, rank, hostnum, ng, ni;
223 /* Many MPI implementations do not optimize MPI_Allreduce
224 * (and probably also other global communication calls)
225 * for multi-core nodes connected by a network.
226 * We can optimize such communication by using one MPI call
227 * within each node and one between the nodes.
228 * For MVAPICH2 and Intel MPI this reduces the time for
229 * the global_stat communication by 25%
230 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
231 * B. Hess, November 2007
237 #ifndef GMX_THREAD_MPI
239 MPI_Comm_size(cr->mpi_comm_mygroup, &n);
240 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
242 hostnum = gmx_hostname_num();
246 fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
250 /* The intra-node communicator, split on node number */
251 MPI_Comm_split(cr->mpi_comm_mygroup, hostnum, rank, &nc->comm_intra);
252 MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
255 fprintf(debug, "In gmx_setup_nodecomm: node rank %d rank_intra %d\n",
256 rank, nc->rank_intra);
258 /* The inter-node communicator, split on rank_intra.
259 * We actually only need the one for rank=0,
260 * but it is easier to create them all.
262 MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
263 /* Check if this really created two step communication */
264 MPI_Comm_size(nc->comm_inter, &ng);
265 MPI_Comm_size(nc->comm_intra, &ni);
268 fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
272 if (getenv("GMX_NO_NODECOMM") == NULL &&
273 ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
278 fprintf(fplog, "Using two step summing over %d groups of on average %.1f processes\n\n",
279 ng, (real)n/(real)ng);
281 if (nc->rank_intra > 0)
283 MPI_Comm_free(&nc->comm_inter);
288 /* One group or all processes in a separate group, use normal summing */
289 MPI_Comm_free(&nc->comm_inter);
290 MPI_Comm_free(&nc->comm_intra);
293 fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
298 /* tMPI runs only on a single node so just use the nodeid */
299 nc->rank_intra = cr->nodeid;
303 void gmx_init_intranode_counters(t_commrec *cr)
305 /* counters for PP+PME and PP-only processes on my physical node */
306 int nrank_intranode, rank_intranode;
307 int nrank_pp_intranode, rank_pp_intranode;
308 /* thread-MPI is not initialized when not running in parallel */
309 #if defined GMX_MPI && !defined GMX_THREAD_MPI
310 int nrank_world, rank_world;
311 int i, mynum, *num, *num_s, *num_pp, *num_pp_s;
313 MPI_Comm_size(MPI_COMM_WORLD, &nrank_world);
314 MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
316 /* Get the node number from the hostname to identify the nodes */
317 mynum = gmx_hostname_num();
319 /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
320 snew(num, nrank_world);
321 snew(num_s, nrank_world);
322 snew(num_pp, nrank_world);
323 snew(num_pp_s, nrank_world);
325 num_s[rank_world] = mynum;
326 num_pp_s[rank_world] = (cr->duty & DUTY_PP) ? mynum : -1;
328 MPI_Allreduce(num_s, num, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
329 MPI_Allreduce(num_pp_s, num_pp, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
333 nrank_pp_intranode = 0;
334 rank_pp_intranode = 0;
335 for (i = 0; i < nrank_world; i++)
345 if ((cr->duty & DUTY_PP) && num_pp[i] == mynum)
347 nrank_pp_intranode++;
359 /* Serial or thread-MPI code: we run within a single physical node */
360 nrank_intranode = cr->nnodes;
361 rank_intranode = cr->sim_nodeid;
362 nrank_pp_intranode = cr->nnodes - cr->npmenodes;
363 rank_pp_intranode = cr->nodeid;
369 if (cr->duty & DUTY_PP && cr->duty & DUTY_PME)
371 sprintf(sbuf, "PP+PME");
375 sprintf(sbuf, "%s", cr->duty & DUTY_PP ? "PP" : "PME");
377 fprintf(debug, "On %3s node %d: nrank_intranode=%d, rank_intranode=%d, "
378 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
379 sbuf, cr->sim_nodeid,
380 nrank_intranode, rank_intranode,
381 nrank_pp_intranode, rank_pp_intranode);
384 cr->nrank_intranode = nrank_intranode;
385 cr->rank_intranode = rank_intranode;
386 cr->nrank_pp_intranode = nrank_pp_intranode;
387 cr->rank_pp_intranode = rank_pp_intranode;
391 void gmx_barrier(const t_commrec *cr)
394 gmx_call("gmx_barrier");
396 MPI_Barrier(cr->mpi_comm_mygroup);
400 void gmx_abort(int gmx_unused noderank, int gmx_unused nnodes, int gmx_unused errorno)
403 gmx_call("gmx_abort");
405 #ifdef GMX_THREAD_MPI
406 fprintf(stderr, "Halting program %s\n", ShortProgram());
412 fprintf(stderr, "Halting parallel program %s on CPU %d out of %d\n",
413 ShortProgram(), noderank, nnodes);
417 fprintf(stderr, "Halting program %s\n", ShortProgram());
421 MPI_Abort(MPI_COMM_WORLD, errorno);
427 void gmx_bcast(int nbytes, void *b, const t_commrec *cr)
430 gmx_call("gmx_bast");
432 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
436 void gmx_bcast_sim(int nbytes, void *b, const t_commrec *cr)
439 gmx_call("gmx_bast");
441 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
445 void gmx_sumd(int nr, double r[], const t_commrec *cr)
448 gmx_call("gmx_sumd");
450 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
453 if (cr->nc.rank_intra == 0)
455 /* Use two step summing. */
456 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
458 /* Sum the roots of the internal (intra) buffers. */
459 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
464 /* This is here because of the silly MPI specification
465 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
466 MPI_Reduce(r, NULL, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
468 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
472 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
473 cr->mpi_comm_mygroup);
478 if (nr > cr->mpb->dbuf_alloc)
480 cr->mpb->dbuf_alloc = nr;
481 srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
485 /* Use two step summing */
486 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
487 if (cr->nc.rank_intra == 0)
489 /* Sum with the buffers reversed */
490 MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
493 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
497 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
498 cr->mpi_comm_mygroup);
499 for (i = 0; i < nr; i++)
501 r[i] = cr->mpb->dbuf[i];
508 void gmx_sumf(int nr, float r[], const t_commrec *cr)
511 gmx_call("gmx_sumf");
513 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
516 /* Use two step summing. */
517 if (cr->nc.rank_intra == 0)
519 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
521 /* Sum the roots of the internal (intra) buffers */
522 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
527 /* This is here because of the silly MPI specification
528 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
529 MPI_Reduce(r, NULL, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
531 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
535 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
540 if (nr > cr->mpb->fbuf_alloc)
542 cr->mpb->fbuf_alloc = nr;
543 srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
547 /* Use two step summing */
548 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
549 if (cr->nc.rank_intra == 0)
551 /* Sum with the buffers reversed */
552 MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
555 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
559 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
560 cr->mpi_comm_mygroup);
561 for (i = 0; i < nr; i++)
563 r[i] = cr->mpb->fbuf[i];
570 void gmx_sumi(int nr, int r[], const t_commrec *cr)
573 gmx_call("gmx_sumi");
575 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
578 /* Use two step summing */
579 if (cr->nc.rank_intra == 0)
581 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
582 /* Sum with the buffers reversed */
583 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
587 /* This is here because of the silly MPI specification
588 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
589 MPI_Reduce(r, NULL, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
591 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
595 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
600 if (nr > cr->mpb->ibuf_alloc)
602 cr->mpb->ibuf_alloc = nr;
603 srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
607 /* Use two step summing */
608 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
609 if (cr->nc.rank_intra == 0)
611 /* Sum with the buffers reversed */
612 MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
614 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
618 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
619 for (i = 0; i < nr; i++)
621 r[i] = cr->mpb->ibuf[i];
628 void gmx_sumli(int nr, gmx_large_int_t r[], const t_commrec *cr)
631 gmx_call("gmx_sumli");
633 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
636 /* Use two step summing */
637 if (cr->nc.rank_intra == 0)
639 MPI_Reduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM, 0,
641 /* Sum with the buffers reversed */
642 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
647 /* This is here because of the silly MPI specification
648 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
649 MPI_Reduce(r, NULL, nr, GMX_MPI_LARGE_INT, MPI_SUM, 0, cr->nc.comm_intra);
651 MPI_Bcast(r, nr, GMX_MPI_LARGE_INT, 0, cr->nc.comm_intra);
655 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM, cr->mpi_comm_mygroup);
660 if (nr > cr->mpb->libuf_alloc)
662 cr->mpb->libuf_alloc = nr;
663 srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
667 /* Use two step summing */
668 MPI_Allreduce(r, cr->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
670 if (cr->nc.rank_intra == 0)
672 /* Sum with the buffers reversed */
673 MPI_Allreduce(cr->mpb->libuf, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
676 MPI_Bcast(r, nr, GMX_MPI_LARGE_INT, 0, cr->nc.comm_intra);
680 MPI_Allreduce(r, cr->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
681 cr->mpi_comm_mygroup);
682 for (i = 0; i < nr; i++)
684 r[i] = cr->mpb->libuf[i];
694 void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
696 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
697 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
699 /* this function is only used in code that is not performance critical,
700 (during setup, when comm_rec is not the appropriate communication
701 structure), so this isn't as bad as it looks. */
706 MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
707 for (i = 0; i < nr; i++)
717 void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
719 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
720 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
722 /* this function is only used in code that is not performance critical,
723 (during setup, when comm_rec is not the appropriate communication
724 structure), so this isn't as bad as it looks. */
729 MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
730 for (i = 0; i < nr; i++)
739 void gmx_sumd_sim(int nr, double r[], const gmx_multisim_t *ms)
742 gmx_call("gmx_sumd_sim");
744 gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
748 void gmx_sumf_sim(int nr, float r[], const gmx_multisim_t *ms)
751 gmx_call("gmx_sumf_sim");
753 gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
757 void gmx_sumi_sim(int nr, int r[], const gmx_multisim_t *ms)
760 gmx_call("gmx_sumi_sim");
762 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
763 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
765 /* this is thread-unsafe, but it will do for now: */
768 if (nr > ms->mpb->ibuf_alloc)
770 ms->mpb->ibuf_alloc = nr;
771 srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
773 MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
774 for (i = 0; i < nr; i++)
776 r[i] = ms->mpb->ibuf[i];
782 void gmx_sumli_sim(int nr, gmx_large_int_t r[], const gmx_multisim_t *ms)
785 gmx_call("gmx_sumli_sim");
787 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
788 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
789 ms->mpi_comm_masters);
791 /* this is thread-unsafe, but it will do for now: */
794 if (nr > ms->mpb->libuf_alloc)
796 ms->mpb->libuf_alloc = nr;
797 srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
799 MPI_Allreduce(r, ms->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
800 ms->mpi_comm_masters);
801 for (i = 0; i < nr; i++)
803 r[i] = ms->mpb->libuf[i];