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 cr->nnodes = gmx_node_num();
73 cr->nodeid = gmx_node_rank();
74 cr->sim_nodeid = cr->nodeid;
75 cr->mpi_comm_mysim = MPI_COMM_WORLD;
76 cr->mpi_comm_mygroup = MPI_COMM_WORLD;
81 t_commrec *init_commrec()
88 if (!gmx_mpi_initialized())
90 gmx_comm("MPI has not been initialized properly");
93 gmx_fill_commrec_from_mpi(cr);
95 /* These should never be accessed */
96 cr->mpi_comm_mysim = NULL;
97 cr->mpi_comm_mygroup = NULL;
100 cr->nodeid = cr->sim_nodeid;
103 // TODO this should be initialized elsewhere
104 cr->duty = (DUTY_PP | DUTY_PME);
106 #if defined GMX_LIB_MPI && !defined MPI_IN_PLACE_EXISTS
107 /* initialize the MPI_IN_PLACE replacement buffers */
109 cr->mpb->ibuf = NULL;
110 cr->mpb->libuf = NULL;
111 cr->mpb->fbuf = NULL;
112 cr->mpb->dbuf = NULL;
113 cr->mpb->ibuf_alloc = 0;
114 cr->mpb->libuf_alloc = 0;
115 cr->mpb->fbuf_alloc = 0;
116 cr->mpb->dbuf_alloc = 0;
122 t_commrec *init_par_threads(const t_commrec *cro)
124 #ifdef GMX_THREAD_MPI
128 /* make a thread-specific commrec */
130 /* now copy the whole thing, so settings like the number of PME nodes
134 /* and we start setting our own thread-specific values for things */
135 MPI_Initialized(&initialized);
138 gmx_comm("Initializing threads without comm");
141 /* No need to do MPI_Init() for thread-MPI. */
142 gmx_fill_commrec_from_mpi(cr);
144 // TODO cr->duty should not be initialized here
145 cr->duty = (DUTY_PP | DUTY_PME);
153 int gmx_node_num(void)
159 (void) MPI_Comm_size(MPI_COMM_WORLD, &i);
164 int gmx_node_rank(void)
170 (void) MPI_Comm_rank(MPI_COMM_WORLD, &i);
176 int gmx_hostname_num()
181 #ifdef GMX_THREAD_MPI
182 /* thread-MPI currently puts the thread number in the process name,
183 * we might want to change this, as this is inconsistent with what
184 * most MPI implementations would do when running on a single node.
188 int resultlen, hostnum, i, j;
189 char mpi_hostname[MPI_MAX_PROCESSOR_NAME], hostnum_str[MPI_MAX_PROCESSOR_NAME];
191 MPI_Get_processor_name(mpi_hostname, &resultlen);
192 /* This procedure can only differentiate nodes with host names
193 * that end on unique numbers.
197 /* Only parse the host name up to the first dot */
198 while (i < resultlen && mpi_hostname[i] != '.')
200 if (isdigit(mpi_hostname[i]))
202 hostnum_str[j++] = mpi_hostname[i];
206 hostnum_str[j] = '\0';
213 /* Use only the last 9 decimals, so we don't overflow an int */
214 hostnum = strtol(hostnum_str + max(0, j-9), NULL, 10);
219 fprintf(debug, "In gmx_setup_nodecomm: hostname '%s', hostnum %d\n",
220 mpi_hostname, hostnum);
227 void gmx_setup_nodecomm(FILE gmx_unused *fplog, t_commrec *cr)
230 int n, rank, hostnum, ng, ni;
232 /* Many MPI implementations do not optimize MPI_Allreduce
233 * (and probably also other global communication calls)
234 * for multi-core nodes connected by a network.
235 * We can optimize such communication by using one MPI call
236 * within each node and one between the nodes.
237 * For MVAPICH2 and Intel MPI this reduces the time for
238 * the global_stat communication by 25%
239 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
240 * B. Hess, November 2007
246 #ifndef GMX_THREAD_MPI
248 MPI_Comm_size(cr->mpi_comm_mygroup, &n);
249 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
251 hostnum = gmx_hostname_num();
255 fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
259 /* The intra-node communicator, split on node number */
260 MPI_Comm_split(cr->mpi_comm_mygroup, hostnum, rank, &nc->comm_intra);
261 MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
264 fprintf(debug, "In gmx_setup_nodecomm: node rank %d rank_intra %d\n",
265 rank, nc->rank_intra);
267 /* The inter-node communicator, split on rank_intra.
268 * We actually only need the one for rank=0,
269 * but it is easier to create them all.
271 MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
272 /* Check if this really created two step communication */
273 MPI_Comm_size(nc->comm_inter, &ng);
274 MPI_Comm_size(nc->comm_intra, &ni);
277 fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
281 if (getenv("GMX_NO_NODECOMM") == NULL &&
282 ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
287 fprintf(fplog, "Using two step summing over %d groups of on average %.1f processes\n\n",
288 ng, (real)n/(real)ng);
290 if (nc->rank_intra > 0)
292 MPI_Comm_free(&nc->comm_inter);
297 /* One group or all processes in a separate group, use normal summing */
298 MPI_Comm_free(&nc->comm_inter);
299 MPI_Comm_free(&nc->comm_intra);
302 fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
307 /* tMPI runs only on a single node so just use the nodeid */
308 nc->rank_intra = cr->nodeid;
312 void gmx_init_intranode_counters(t_commrec *cr)
314 /* counters for PP+PME and PP-only processes on my physical node */
315 int nrank_intranode, rank_intranode;
316 int nrank_pp_intranode, rank_pp_intranode;
317 /* thread-MPI is not initialized when not running in parallel */
318 #if defined GMX_MPI && !defined GMX_THREAD_MPI
319 int nrank_world, rank_world;
320 int i, mynum, *num, *num_s, *num_pp, *num_pp_s;
322 MPI_Comm_size(MPI_COMM_WORLD, &nrank_world);
323 MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
325 /* Get the node number from the hostname to identify the nodes */
326 mynum = gmx_hostname_num();
328 /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
329 snew(num, nrank_world);
330 snew(num_s, nrank_world);
331 snew(num_pp, nrank_world);
332 snew(num_pp_s, nrank_world);
334 num_s[rank_world] = mynum;
335 num_pp_s[rank_world] = (cr->duty & DUTY_PP) ? mynum : -1;
337 MPI_Allreduce(num_s, num, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
338 MPI_Allreduce(num_pp_s, num_pp, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
342 nrank_pp_intranode = 0;
343 rank_pp_intranode = 0;
344 for (i = 0; i < nrank_world; i++)
354 if ((cr->duty & DUTY_PP) && num_pp[i] == mynum)
356 nrank_pp_intranode++;
368 /* Serial or thread-MPI code: we run within a single physical node */
369 nrank_intranode = cr->nnodes;
370 rank_intranode = cr->sim_nodeid;
371 nrank_pp_intranode = cr->nnodes - cr->npmenodes;
372 rank_pp_intranode = cr->nodeid;
378 if (cr->duty & DUTY_PP && cr->duty & DUTY_PME)
380 sprintf(sbuf, "PP+PME");
384 sprintf(sbuf, "%s", cr->duty & DUTY_PP ? "PP" : "PME");
386 fprintf(debug, "On %3s node %d: nrank_intranode=%d, rank_intranode=%d, "
387 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
388 sbuf, cr->sim_nodeid,
389 nrank_intranode, rank_intranode,
390 nrank_pp_intranode, rank_pp_intranode);
393 cr->nrank_intranode = nrank_intranode;
394 cr->rank_intranode = rank_intranode;
395 cr->nrank_pp_intranode = nrank_pp_intranode;
396 cr->rank_pp_intranode = rank_pp_intranode;
400 void gmx_barrier(const t_commrec *cr)
403 gmx_call("gmx_barrier");
405 MPI_Barrier(cr->mpi_comm_mygroup);
409 void gmx_abort(int gmx_unused noderank, int gmx_unused nnodes, int gmx_unused errorno)
412 gmx_call("gmx_abort");
414 #ifdef GMX_THREAD_MPI
415 fprintf(stderr, "Halting program %s\n", ShortProgram());
421 fprintf(stderr, "Halting parallel program %s on CPU %d out of %d\n",
422 ShortProgram(), noderank, nnodes);
426 fprintf(stderr, "Halting program %s\n", ShortProgram());
430 MPI_Abort(MPI_COMM_WORLD, errorno);
436 void gmx_bcast(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_mygroup);
445 void gmx_bcast_sim(int nbytes, void *b, const t_commrec *cr)
448 gmx_call("gmx_bast");
450 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
454 void gmx_sumd(int nr, double r[], const t_commrec *cr)
457 gmx_call("gmx_sumd");
459 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
462 if (cr->nc.rank_intra == 0)
464 /* Use two step summing. */
465 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
467 /* Sum the roots of the internal (intra) buffers. */
468 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
473 /* This is here because of the silly MPI specification
474 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
475 MPI_Reduce(r, NULL, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
477 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
481 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
482 cr->mpi_comm_mygroup);
487 if (nr > cr->mpb->dbuf_alloc)
489 cr->mpb->dbuf_alloc = nr;
490 srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
494 /* Use two step summing */
495 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
496 if (cr->nc.rank_intra == 0)
498 /* Sum with the buffers reversed */
499 MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
502 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
506 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
507 cr->mpi_comm_mygroup);
508 for (i = 0; i < nr; i++)
510 r[i] = cr->mpb->dbuf[i];
517 void gmx_sumf(int nr, float r[], const t_commrec *cr)
520 gmx_call("gmx_sumf");
522 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
525 /* Use two step summing. */
526 if (cr->nc.rank_intra == 0)
528 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
530 /* Sum the roots of the internal (intra) buffers */
531 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
536 /* This is here because of the silly MPI specification
537 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
538 MPI_Reduce(r, NULL, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
540 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
544 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
549 if (nr > cr->mpb->fbuf_alloc)
551 cr->mpb->fbuf_alloc = nr;
552 srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
556 /* Use two step summing */
557 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
558 if (cr->nc.rank_intra == 0)
560 /* Sum with the buffers reversed */
561 MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
564 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
568 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
569 cr->mpi_comm_mygroup);
570 for (i = 0; i < nr; i++)
572 r[i] = cr->mpb->fbuf[i];
579 void gmx_sumi(int nr, int r[], const t_commrec *cr)
582 gmx_call("gmx_sumi");
584 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
587 /* Use two step summing */
588 if (cr->nc.rank_intra == 0)
590 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
591 /* Sum with the buffers reversed */
592 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
596 /* This is here because of the silly MPI specification
597 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
598 MPI_Reduce(r, NULL, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
600 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
604 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
609 if (nr > cr->mpb->ibuf_alloc)
611 cr->mpb->ibuf_alloc = nr;
612 srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
616 /* Use two step summing */
617 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
618 if (cr->nc.rank_intra == 0)
620 /* Sum with the buffers reversed */
621 MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
623 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
627 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
628 for (i = 0; i < nr; i++)
630 r[i] = cr->mpb->ibuf[i];
637 void gmx_sumli(int nr, gmx_large_int_t r[], const t_commrec *cr)
640 gmx_call("gmx_sumli");
642 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
645 /* Use two step summing */
646 if (cr->nc.rank_intra == 0)
648 MPI_Reduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM, 0,
650 /* Sum with the buffers reversed */
651 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
656 /* This is here because of the silly MPI specification
657 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
658 MPI_Reduce(r, NULL, nr, GMX_MPI_LARGE_INT, MPI_SUM, 0, cr->nc.comm_intra);
660 MPI_Bcast(r, nr, GMX_MPI_LARGE_INT, 0, cr->nc.comm_intra);
664 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM, cr->mpi_comm_mygroup);
669 if (nr > cr->mpb->libuf_alloc)
671 cr->mpb->libuf_alloc = nr;
672 srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
676 /* Use two step summing */
677 MPI_Allreduce(r, cr->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
679 if (cr->nc.rank_intra == 0)
681 /* Sum with the buffers reversed */
682 MPI_Allreduce(cr->mpb->libuf, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
685 MPI_Bcast(r, nr, GMX_MPI_LARGE_INT, 0, cr->nc.comm_intra);
689 MPI_Allreduce(r, cr->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
690 cr->mpi_comm_mygroup);
691 for (i = 0; i < nr; i++)
693 r[i] = cr->mpb->libuf[i];
703 void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
705 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
706 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
708 /* this function is only used in code that is not performance critical,
709 (during setup, when comm_rec is not the appropriate communication
710 structure), so this isn't as bad as it looks. */
715 MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
716 for (i = 0; i < nr; i++)
726 void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
728 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
729 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
731 /* this function is only used in code that is not performance critical,
732 (during setup, when comm_rec is not the appropriate communication
733 structure), so this isn't as bad as it looks. */
738 MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
739 for (i = 0; i < nr; i++)
748 void gmx_sumd_sim(int nr, double r[], const gmx_multisim_t *ms)
751 gmx_call("gmx_sumd_sim");
753 gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
757 void gmx_sumf_sim(int nr, float r[], const gmx_multisim_t *ms)
760 gmx_call("gmx_sumf_sim");
762 gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
766 void gmx_sumi_sim(int nr, int r[], const gmx_multisim_t *ms)
769 gmx_call("gmx_sumi_sim");
771 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
772 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
774 /* this is thread-unsafe, but it will do for now: */
777 if (nr > ms->mpb->ibuf_alloc)
779 ms->mpb->ibuf_alloc = nr;
780 srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
782 MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
783 for (i = 0; i < nr; i++)
785 r[i] = ms->mpb->ibuf[i];
791 void gmx_sumli_sim(int nr, gmx_large_int_t r[], const gmx_multisim_t *ms)
794 gmx_call("gmx_sumli_sim");
796 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
797 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
798 ms->mpi_comm_masters);
800 /* this is thread-unsafe, but it will do for now: */
803 if (nr > ms->mpb->libuf_alloc)
805 ms->mpb->libuf_alloc = nr;
806 srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
808 MPI_Allreduce(r, ms->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
809 ms->mpi_comm_masters);
810 for (i = 0; i < nr; i++)
812 r[i] = ms->mpb->libuf[i];