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_do_mpi_init(int gmx_unused *argc, char gmx_unused ***argv)
70 gmx_call("gmx_do_mpi_init");
72 if (!gmx_mpi_initialized())
76 (void) fah_MPI_Init(argc, argv);
78 (void) MPI_Init(argc, argv);
85 void gmx_fill_commrec_from_mpi(t_commrec *cr)
88 gmx_call("gmx_fill_commrec_from_mpi");
91 int resultlen; /* actual length of node name */
95 char mpi_hostname[MPI_MAX_PROCESSOR_NAME];
97 mpi_num_nodes = gmx_node_num();
98 mpi_my_rank = gmx_node_rank();
99 (void) MPI_Get_processor_name( mpi_hostname, &resultlen );
104 fprintf(debug, "NNODES=%d, MYRANK=%d, HOSTNAME=%s\n",
105 mpi_num_nodes, mpi_my_rank, mpi_hostname);
109 cr->nnodes = mpi_num_nodes;
110 cr->nodeid = mpi_my_rank;
111 cr->sim_nodeid = mpi_my_rank;
112 cr->mpi_comm_mysim = MPI_COMM_WORLD;
113 cr->mpi_comm_mygroup = MPI_COMM_WORLD;
117 int gmx_node_num(void)
123 (void) MPI_Comm_size(MPI_COMM_WORLD, &i);
128 int gmx_node_rank(void)
134 (void) MPI_Comm_rank(MPI_COMM_WORLD, &i);
140 int gmx_hostname_num()
145 #ifdef GMX_THREAD_MPI
146 /* thread-MPI currently puts the thread number in the process name,
147 * we might want to change this, as this is inconsistent with what
148 * most MPI implementations would do when running on a single node.
152 int resultlen, hostnum, i, j;
153 char mpi_hostname[MPI_MAX_PROCESSOR_NAME], hostnum_str[MPI_MAX_PROCESSOR_NAME];
155 MPI_Get_processor_name(mpi_hostname, &resultlen);
156 /* This procedure can only differentiate nodes with host names
157 * that end on unique numbers.
161 /* Only parse the host name up to the first dot */
162 while (i < resultlen && mpi_hostname[i] != '.')
164 if (isdigit(mpi_hostname[i]))
166 hostnum_str[j++] = mpi_hostname[i];
170 hostnum_str[j] = '\0';
177 /* Use only the last 9 decimals, so we don't overflow an int */
178 hostnum = strtol(hostnum_str + max(0, j-9), NULL, 10);
183 fprintf(debug, "In gmx_setup_nodecomm: hostname '%s', hostnum %d\n",
184 mpi_hostname, hostnum);
191 void gmx_setup_nodecomm(FILE gmx_unused *fplog, t_commrec *cr)
194 int n, rank, hostnum, ng, ni;
196 /* Many MPI implementations do not optimize MPI_Allreduce
197 * (and probably also other global communication calls)
198 * for multi-core nodes connected by a network.
199 * We can optimize such communication by using one MPI call
200 * within each node and one between the nodes.
201 * For MVAPICH2 and Intel MPI this reduces the time for
202 * the global_stat communication by 25%
203 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
204 * B. Hess, November 2007
210 #ifndef GMX_THREAD_MPI
212 MPI_Comm_size(cr->mpi_comm_mygroup, &n);
213 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
215 hostnum = gmx_hostname_num();
219 fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
223 /* The intra-node communicator, split on node number */
224 MPI_Comm_split(cr->mpi_comm_mygroup, hostnum, rank, &nc->comm_intra);
225 MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
228 fprintf(debug, "In gmx_setup_nodecomm: node rank %d rank_intra %d\n",
229 rank, nc->rank_intra);
231 /* The inter-node communicator, split on rank_intra.
232 * We actually only need the one for rank=0,
233 * but it is easier to create them all.
235 MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
236 /* Check if this really created two step communication */
237 MPI_Comm_size(nc->comm_inter, &ng);
238 MPI_Comm_size(nc->comm_intra, &ni);
241 fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
245 if (getenv("GMX_NO_NODECOMM") == NULL &&
246 ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
251 fprintf(fplog, "Using two step summing over %d groups of on average %.1f processes\n\n",
252 ng, (real)n/(real)ng);
254 if (nc->rank_intra > 0)
256 MPI_Comm_free(&nc->comm_inter);
261 /* One group or all processes in a separate group, use normal summing */
262 MPI_Comm_free(&nc->comm_inter);
263 MPI_Comm_free(&nc->comm_intra);
266 fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
271 /* tMPI runs only on a single node so just use the nodeid */
272 nc->rank_intra = cr->nodeid;
276 void gmx_init_intranode_counters(t_commrec *cr)
278 /* counters for PP+PME and PP-only processes on my physical node */
279 int nrank_intranode, rank_intranode;
280 int nrank_pp_intranode, rank_pp_intranode;
281 /* thread-MPI is not initialized when not running in parallel */
282 #if defined GMX_MPI && !defined GMX_THREAD_MPI
283 int nrank_world, rank_world;
284 int i, mynum, *num, *num_s, *num_pp, *num_pp_s;
286 MPI_Comm_size(MPI_COMM_WORLD, &nrank_world);
287 MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
289 /* Get the node number from the hostname to identify the nodes */
290 mynum = gmx_hostname_num();
292 /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
293 snew(num, nrank_world);
294 snew(num_s, nrank_world);
295 snew(num_pp, nrank_world);
296 snew(num_pp_s, nrank_world);
298 num_s[rank_world] = mynum;
299 num_pp_s[rank_world] = (cr->duty & DUTY_PP) ? mynum : -1;
301 MPI_Allreduce(num_s, num, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
302 MPI_Allreduce(num_pp_s, num_pp, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
306 nrank_pp_intranode = 0;
307 rank_pp_intranode = 0;
308 for (i = 0; i < nrank_world; i++)
318 if ((cr->duty & DUTY_PP) && num_pp[i] == mynum)
320 nrank_pp_intranode++;
332 /* Serial or thread-MPI code: we run within a single physical node */
333 nrank_intranode = cr->nnodes;
334 rank_intranode = cr->sim_nodeid;
335 nrank_pp_intranode = cr->nnodes - cr->npmenodes;
336 rank_pp_intranode = cr->nodeid;
342 if (cr->duty & DUTY_PP && cr->duty & DUTY_PME)
344 sprintf(sbuf, "PP+PME");
348 sprintf(sbuf, "%s", cr->duty & DUTY_PP ? "PP" : "PME");
350 fprintf(debug, "On %3s node %d: nrank_intranode=%d, rank_intranode=%d, "
351 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
352 sbuf, cr->sim_nodeid,
353 nrank_intranode, rank_intranode,
354 nrank_pp_intranode, rank_pp_intranode);
357 cr->nrank_intranode = nrank_intranode;
358 cr->rank_intranode = rank_intranode;
359 cr->nrank_pp_intranode = nrank_pp_intranode;
360 cr->rank_pp_intranode = rank_pp_intranode;
364 void gmx_barrier(const t_commrec *cr)
367 gmx_call("gmx_barrier");
369 MPI_Barrier(cr->mpi_comm_mygroup);
373 void gmx_abort(int gmx_unused noderank, int gmx_unused nnodes, int gmx_unused errorno)
376 gmx_call("gmx_abort");
378 #ifdef GMX_THREAD_MPI
379 fprintf(stderr, "Halting program %s\n", ShortProgram());
385 fprintf(stderr, "Halting parallel program %s on CPU %d out of %d\n",
386 ShortProgram(), noderank, nnodes);
390 fprintf(stderr, "Halting program %s\n", ShortProgram());
394 MPI_Abort(MPI_COMM_WORLD, errorno);
400 void gmx_bcast(int nbytes, void *b, const t_commrec *cr)
403 gmx_call("gmx_bast");
405 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
409 void gmx_bcast_sim(int nbytes, void *b, const t_commrec *cr)
412 gmx_call("gmx_bast");
414 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
418 void gmx_sumd(int nr, double r[], const t_commrec *cr)
421 gmx_call("gmx_sumd");
423 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
426 if (cr->nc.rank_intra == 0)
428 /* Use two step summing. */
429 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
431 /* Sum the roots of the internal (intra) buffers. */
432 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
437 /* This is here because of the silly MPI specification
438 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
439 MPI_Reduce(r, NULL, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
441 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
445 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
446 cr->mpi_comm_mygroup);
451 if (nr > cr->mpb->dbuf_alloc)
453 cr->mpb->dbuf_alloc = nr;
454 srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
458 /* Use two step summing */
459 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
460 if (cr->nc.rank_intra == 0)
462 /* Sum with the buffers reversed */
463 MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
466 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
470 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
471 cr->mpi_comm_mygroup);
472 for (i = 0; i < nr; i++)
474 r[i] = cr->mpb->dbuf[i];
481 void gmx_sumf(int nr, float r[], const t_commrec *cr)
484 gmx_call("gmx_sumf");
486 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
489 /* Use two step summing. */
490 if (cr->nc.rank_intra == 0)
492 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
494 /* Sum the roots of the internal (intra) buffers */
495 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
500 /* This is here because of the silly MPI specification
501 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
502 MPI_Reduce(r, NULL, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
504 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
508 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
513 if (nr > cr->mpb->fbuf_alloc)
515 cr->mpb->fbuf_alloc = nr;
516 srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
520 /* Use two step summing */
521 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
522 if (cr->nc.rank_intra == 0)
524 /* Sum with the buffers reversed */
525 MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
528 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
532 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
533 cr->mpi_comm_mygroup);
534 for (i = 0; i < nr; i++)
536 r[i] = cr->mpb->fbuf[i];
543 void gmx_sumi(int nr, int r[], const t_commrec *cr)
546 gmx_call("gmx_sumi");
548 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
551 /* Use two step summing */
552 if (cr->nc.rank_intra == 0)
554 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
555 /* Sum with the buffers reversed */
556 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
560 /* This is here because of the silly MPI specification
561 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
562 MPI_Reduce(r, NULL, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
564 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
568 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
573 if (nr > cr->mpb->ibuf_alloc)
575 cr->mpb->ibuf_alloc = nr;
576 srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
580 /* Use two step summing */
581 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
582 if (cr->nc.rank_intra == 0)
584 /* Sum with the buffers reversed */
585 MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
587 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
591 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
592 for (i = 0; i < nr; i++)
594 r[i] = cr->mpb->ibuf[i];
601 void gmx_sumli(int nr, gmx_large_int_t r[], const t_commrec *cr)
604 gmx_call("gmx_sumli");
606 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
609 /* Use two step summing */
610 if (cr->nc.rank_intra == 0)
612 MPI_Reduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM, 0,
614 /* Sum with the buffers reversed */
615 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
620 /* This is here because of the silly MPI specification
621 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
622 MPI_Reduce(r, NULL, nr, GMX_MPI_LARGE_INT, MPI_SUM, 0, cr->nc.comm_intra);
624 MPI_Bcast(r, nr, GMX_MPI_LARGE_INT, 0, cr->nc.comm_intra);
628 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM, cr->mpi_comm_mygroup);
633 if (nr > cr->mpb->libuf_alloc)
635 cr->mpb->libuf_alloc = nr;
636 srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
640 /* Use two step summing */
641 MPI_Allreduce(r, cr->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
643 if (cr->nc.rank_intra == 0)
645 /* Sum with the buffers reversed */
646 MPI_Allreduce(cr->mpb->libuf, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
649 MPI_Bcast(r, nr, GMX_MPI_LARGE_INT, 0, cr->nc.comm_intra);
653 MPI_Allreduce(r, cr->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
654 cr->mpi_comm_mygroup);
655 for (i = 0; i < nr; i++)
657 r[i] = cr->mpb->libuf[i];
667 void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
669 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
670 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
672 /* this function is only used in code that is not performance critical,
673 (during setup, when comm_rec is not the appropriate communication
674 structure), so this isn't as bad as it looks. */
679 MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
680 for (i = 0; i < nr; i++)
690 void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
692 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
693 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
695 /* this function is only used in code that is not performance critical,
696 (during setup, when comm_rec is not the appropriate communication
697 structure), so this isn't as bad as it looks. */
702 MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
703 for (i = 0; i < nr; i++)
712 void gmx_sumd_sim(int nr, double r[], const gmx_multisim_t *ms)
715 gmx_call("gmx_sumd_sim");
717 gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
721 void gmx_sumf_sim(int nr, float r[], const gmx_multisim_t *ms)
724 gmx_call("gmx_sumf_sim");
726 gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
730 void gmx_sumi_sim(int nr, int r[], const gmx_multisim_t *ms)
733 gmx_call("gmx_sumi_sim");
735 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
736 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
738 /* this is thread-unsafe, but it will do for now: */
741 if (nr > ms->mpb->ibuf_alloc)
743 ms->mpb->ibuf_alloc = nr;
744 srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
746 MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
747 for (i = 0; i < nr; i++)
749 r[i] = ms->mpb->ibuf[i];
755 void gmx_sumli_sim(int nr, gmx_large_int_t r[], const gmx_multisim_t *ms)
758 gmx_call("gmx_sumli_sim");
760 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
761 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
762 ms->mpi_comm_masters);
764 /* this is thread-unsafe, but it will do for now: */
767 if (nr > ms->mpb->libuf_alloc)
769 ms->mpb->libuf_alloc = nr;
770 srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
772 MPI_Allreduce(r, ms->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
773 ms->mpi_comm_masters);
774 for (i = 0; i < nr; i++)
776 r[i] = ms->mpb->libuf[i];
783 void gmx_finalize_mpi(void)
786 /* Compiled without MPI, no MPI finalizing needed */
792 if (!gmx_mpi_initialized())
796 /* just as a check; we don't want to finalize twice */
797 MPI_Finalized(&finalized);
803 /* We sync the processes here to try to avoid problems
804 * with buggy MPI implementations that could cause
805 * unfinished processes to terminate.
807 MPI_Barrier(MPI_COMM_WORLD);
810 if (DOMAINDECOMP(cr)) {
811 if (cr->npmenodes > 0 || cr->dd->bCartesian)
812 MPI_Comm_free(&cr->mpi_comm_mygroup);
813 if (cr->dd->bCartesian)
814 MPI_Comm_free(&cr->mpi_comm_mysim);
818 /* Apparently certain mpich implementations cause problems
819 * with MPI_Finalize. In that case comment out MPI_Finalize.
823 fprintf(debug, "Will call MPI_Finalize now\n");
826 ret = MPI_Finalize();
829 fprintf(debug, "Return code from MPI_Finalize = %d\n", ret);