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 int gmx_setup(int gmx_unused *argc, char gmx_unused **argv, int *nnodes)
70 gmx_call("gmx_setup");
74 int resultlen; /* actual length of node name */
78 char mpi_hostname[MPI_MAX_PROCESSOR_NAME];
80 /* Call the MPI routines */
83 (void) fah_MPI_Init(argc, &argv);
85 (void) MPI_Init(argc, &argv);
88 (void) MPI_Comm_size( MPI_COMM_WORLD, &mpi_num_nodes );
89 (void) MPI_Comm_rank( MPI_COMM_WORLD, &mpi_my_rank );
90 (void) MPI_Get_processor_name( mpi_hostname, &resultlen );
95 fprintf(debug, "NNODES=%d, MYRANK=%d, HOSTNAME=%s\n",
96 mpi_num_nodes, mpi_my_rank, mpi_hostname);
100 *nnodes = mpi_num_nodes;
106 int gmx_node_num(void)
112 (void) MPI_Comm_size(MPI_COMM_WORLD, &i);
117 int gmx_node_rank(void)
123 (void) MPI_Comm_rank(MPI_COMM_WORLD, &i);
129 int gmx_hostname_num()
134 #ifdef GMX_THREAD_MPI
135 /* thread-MPI currently puts the thread number in the process name,
136 * we might want to change this, as this is inconsistent with what
137 * most MPI implementations would do when running on a single node.
141 int resultlen, hostnum, i, j;
142 char mpi_hostname[MPI_MAX_PROCESSOR_NAME], hostnum_str[MPI_MAX_PROCESSOR_NAME];
144 MPI_Get_processor_name(mpi_hostname, &resultlen);
145 /* This procedure can only differentiate nodes with host names
146 * that end on unique numbers.
150 /* Only parse the host name up to the first dot */
151 while (i < resultlen && mpi_hostname[i] != '.')
153 if (isdigit(mpi_hostname[i]))
155 hostnum_str[j++] = mpi_hostname[i];
159 hostnum_str[j] = '\0';
166 /* Use only the last 9 decimals, so we don't overflow an int */
167 hostnum = strtol(hostnum_str + max(0, j-9), NULL, 10);
172 fprintf(debug, "In gmx_setup_nodecomm: hostname '%s', hostnum %d\n",
173 mpi_hostname, hostnum);
180 void gmx_setup_nodecomm(FILE gmx_unused *fplog, t_commrec *cr)
183 int n, rank, hostnum, ng, ni;
185 /* Many MPI implementations do not optimize MPI_Allreduce
186 * (and probably also other global communication calls)
187 * for multi-core nodes connected by a network.
188 * We can optimize such communication by using one MPI call
189 * within each node and one between the nodes.
190 * For MVAPICH2 and Intel MPI this reduces the time for
191 * the global_stat communication by 25%
192 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
193 * B. Hess, November 2007
199 #ifndef GMX_THREAD_MPI
201 MPI_Comm_size(cr->mpi_comm_mygroup, &n);
202 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
204 hostnum = gmx_hostname_num();
208 fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
212 /* The intra-node communicator, split on node number */
213 MPI_Comm_split(cr->mpi_comm_mygroup, hostnum, rank, &nc->comm_intra);
214 MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
217 fprintf(debug, "In gmx_setup_nodecomm: node rank %d rank_intra %d\n",
218 rank, nc->rank_intra);
220 /* The inter-node communicator, split on rank_intra.
221 * We actually only need the one for rank=0,
222 * but it is easier to create them all.
224 MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
225 /* Check if this really created two step communication */
226 MPI_Comm_size(nc->comm_inter, &ng);
227 MPI_Comm_size(nc->comm_intra, &ni);
230 fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
234 if (getenv("GMX_NO_NODECOMM") == NULL &&
235 ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
240 fprintf(fplog, "Using two step summing over %d groups of on average %.1f processes\n\n",
241 ng, (real)n/(real)ng);
243 if (nc->rank_intra > 0)
245 MPI_Comm_free(&nc->comm_inter);
250 /* One group or all processes in a separate group, use normal summing */
251 MPI_Comm_free(&nc->comm_inter);
252 MPI_Comm_free(&nc->comm_intra);
255 fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
260 /* tMPI runs only on a single node so just use the nodeid */
261 nc->rank_intra = cr->nodeid;
265 void gmx_init_intranode_counters(t_commrec *cr)
267 /* counters for PP+PME and PP-only processes on my physical node */
268 int nrank_intranode, rank_intranode;
269 int nrank_pp_intranode, rank_pp_intranode;
270 /* thread-MPI is not initialized when not running in parallel */
271 #if defined GMX_MPI && !defined GMX_THREAD_MPI
272 int nrank_world, rank_world;
273 int i, mynum, *num, *num_s, *num_pp, *num_pp_s;
275 MPI_Comm_size(MPI_COMM_WORLD, &nrank_world);
276 MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
278 /* Get the node number from the hostname to identify the nodes */
279 mynum = gmx_hostname_num();
281 /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
282 snew(num, nrank_world);
283 snew(num_s, nrank_world);
284 snew(num_pp, nrank_world);
285 snew(num_pp_s, nrank_world);
287 num_s[rank_world] = mynum;
288 num_pp_s[rank_world] = (cr->duty & DUTY_PP) ? mynum : -1;
290 MPI_Allreduce(num_s, num, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
291 MPI_Allreduce(num_pp_s, num_pp, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
295 nrank_pp_intranode = 0;
296 rank_pp_intranode = 0;
297 for (i = 0; i < nrank_world; i++)
307 if ((cr->duty & DUTY_PP) && num_pp[i] == mynum)
309 nrank_pp_intranode++;
321 /* Serial or thread-MPI code: we run within a single physical node */
322 nrank_intranode = cr->nnodes;
323 rank_intranode = cr->sim_nodeid;
324 nrank_pp_intranode = cr->nnodes - cr->npmenodes;
325 rank_pp_intranode = cr->nodeid;
331 if (cr->duty & DUTY_PP && cr->duty & DUTY_PME)
333 sprintf(sbuf, "PP+PME");
337 sprintf(sbuf, "%s", cr->duty & DUTY_PP ? "PP" : "PME");
339 fprintf(debug, "On %3s node %d: nrank_intranode=%d, rank_intranode=%d, "
340 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
341 sbuf, cr->sim_nodeid,
342 nrank_intranode, rank_intranode,
343 nrank_pp_intranode, rank_pp_intranode);
346 cr->nrank_intranode = nrank_intranode;
347 cr->rank_intranode = rank_intranode;
348 cr->nrank_pp_intranode = nrank_pp_intranode;
349 cr->rank_pp_intranode = rank_pp_intranode;
353 void gmx_barrier(const t_commrec *cr)
356 gmx_call("gmx_barrier");
358 MPI_Barrier(cr->mpi_comm_mygroup);
362 void gmx_abort(int gmx_unused noderank, int gmx_unused nnodes, int gmx_unused errorno)
365 gmx_call("gmx_abort");
367 #ifdef GMX_THREAD_MPI
368 fprintf(stderr, "Halting program %s\n", ShortProgram());
374 fprintf(stderr, "Halting parallel program %s on CPU %d out of %d\n",
375 ShortProgram(), noderank, nnodes);
379 fprintf(stderr, "Halting program %s\n", ShortProgram());
383 MPI_Abort(MPI_COMM_WORLD, errorno);
389 void gmx_bcast(int nbytes, void *b, const t_commrec *cr)
392 gmx_call("gmx_bast");
394 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
398 void gmx_bcast_sim(int nbytes, void *b, const t_commrec *cr)
401 gmx_call("gmx_bast");
403 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
407 void gmx_sumd(int nr, double r[], const t_commrec *cr)
410 gmx_call("gmx_sumd");
412 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
415 if (cr->nc.rank_intra == 0)
417 /* Use two step summing. */
418 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
420 /* Sum the roots of the internal (intra) buffers. */
421 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
426 /* This is here because of the silly MPI specification
427 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
428 MPI_Reduce(r, NULL, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
430 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
434 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
435 cr->mpi_comm_mygroup);
440 if (nr > cr->mpb->dbuf_alloc)
442 cr->mpb->dbuf_alloc = nr;
443 srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
447 /* Use two step summing */
448 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
449 if (cr->nc.rank_intra == 0)
451 /* Sum with the buffers reversed */
452 MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
455 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
459 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
460 cr->mpi_comm_mygroup);
461 for (i = 0; i < nr; i++)
463 r[i] = cr->mpb->dbuf[i];
470 void gmx_sumf(int nr, float r[], const t_commrec *cr)
473 gmx_call("gmx_sumf");
475 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
478 /* Use two step summing. */
479 if (cr->nc.rank_intra == 0)
481 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
483 /* Sum the roots of the internal (intra) buffers */
484 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
489 /* This is here because of the silly MPI specification
490 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
491 MPI_Reduce(r, NULL, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
493 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
497 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
502 if (nr > cr->mpb->fbuf_alloc)
504 cr->mpb->fbuf_alloc = nr;
505 srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
509 /* Use two step summing */
510 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
511 if (cr->nc.rank_intra == 0)
513 /* Sum with the buffers reversed */
514 MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
517 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
521 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
522 cr->mpi_comm_mygroup);
523 for (i = 0; i < nr; i++)
525 r[i] = cr->mpb->fbuf[i];
532 void gmx_sumi(int nr, int r[], const t_commrec *cr)
535 gmx_call("gmx_sumi");
537 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
540 /* Use two step summing */
541 if (cr->nc.rank_intra == 0)
543 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
544 /* Sum with the buffers reversed */
545 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
549 /* This is here because of the silly MPI specification
550 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
551 MPI_Reduce(r, NULL, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
553 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
557 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
562 if (nr > cr->mpb->ibuf_alloc)
564 cr->mpb->ibuf_alloc = nr;
565 srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
569 /* Use two step summing */
570 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
571 if (cr->nc.rank_intra == 0)
573 /* Sum with the buffers reversed */
574 MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
576 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
580 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
581 for (i = 0; i < nr; i++)
583 r[i] = cr->mpb->ibuf[i];
590 void gmx_sumli(int nr, gmx_large_int_t r[], const t_commrec *cr)
593 gmx_call("gmx_sumli");
595 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
598 /* Use two step summing */
599 if (cr->nc.rank_intra == 0)
601 MPI_Reduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM, 0,
603 /* Sum with the buffers reversed */
604 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
609 /* This is here because of the silly MPI specification
610 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
611 MPI_Reduce(r, NULL, nr, GMX_MPI_LARGE_INT, MPI_SUM, 0, cr->nc.comm_intra);
613 MPI_Bcast(r, nr, GMX_MPI_LARGE_INT, 0, cr->nc.comm_intra);
617 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM, cr->mpi_comm_mygroup);
622 if (nr > cr->mpb->libuf_alloc)
624 cr->mpb->libuf_alloc = nr;
625 srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
629 /* Use two step summing */
630 MPI_Allreduce(r, cr->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
632 if (cr->nc.rank_intra == 0)
634 /* Sum with the buffers reversed */
635 MPI_Allreduce(cr->mpb->libuf, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
638 MPI_Bcast(r, nr, GMX_MPI_LARGE_INT, 0, cr->nc.comm_intra);
642 MPI_Allreduce(r, cr->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
643 cr->mpi_comm_mygroup);
644 for (i = 0; i < nr; i++)
646 r[i] = cr->mpb->libuf[i];
656 void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
658 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
659 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
661 /* this function is only used in code that is not performance critical,
662 (during setup, when comm_rec is not the appropriate communication
663 structure), so this isn't as bad as it looks. */
668 MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
669 for (i = 0; i < nr; i++)
679 void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
681 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
682 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
684 /* this function is only used in code that is not performance critical,
685 (during setup, when comm_rec is not the appropriate communication
686 structure), so this isn't as bad as it looks. */
691 MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
692 for (i = 0; i < nr; i++)
701 void gmx_sumd_sim(int nr, double r[], const gmx_multisim_t *ms)
704 gmx_call("gmx_sumd_sim");
706 gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
710 void gmx_sumf_sim(int nr, float r[], const gmx_multisim_t *ms)
713 gmx_call("gmx_sumf_sim");
715 gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
719 void gmx_sumi_sim(int nr, int r[], const gmx_multisim_t *ms)
722 gmx_call("gmx_sumi_sim");
724 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
725 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
727 /* this is thread-unsafe, but it will do for now: */
730 if (nr > ms->mpb->ibuf_alloc)
732 ms->mpb->ibuf_alloc = nr;
733 srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
735 MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
736 for (i = 0; i < nr; i++)
738 r[i] = ms->mpb->ibuf[i];
744 void gmx_sumli_sim(int nr, gmx_large_int_t r[], const gmx_multisim_t *ms)
747 gmx_call("gmx_sumli_sim");
749 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
750 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
751 ms->mpi_comm_masters);
753 /* this is thread-unsafe, but it will do for now: */
756 if (nr > ms->mpb->libuf_alloc)
758 ms->mpb->libuf_alloc = nr;
759 srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
761 MPI_Allreduce(r, ms->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
762 ms->mpi_comm_masters);
763 for (i = 0; i < nr; i++)
765 r[i] = ms->mpb->libuf[i];
772 void gmx_finalize_par(void)
775 /* Compiled without MPI, no MPI finalizing needed */
778 int initialized, finalized;
781 MPI_Initialized(&initialized);
786 /* just as a check; we don't want to finalize twice */
787 MPI_Finalized(&finalized);
793 /* We sync the processes here to try to avoid problems
794 * with buggy MPI implementations that could cause
795 * unfinished processes to terminate.
797 MPI_Barrier(MPI_COMM_WORLD);
800 if (DOMAINDECOMP(cr)) {
801 if (cr->npmenodes > 0 || cr->dd->bCartesian)
802 MPI_Comm_free(&cr->mpi_comm_mygroup);
803 if (cr->dd->bCartesian)
804 MPI_Comm_free(&cr->mpi_comm_mysim);
808 /* Apparently certain mpich implementations cause problems
809 * with MPI_Finalize. In that case comment out MPI_Finalize.
813 fprintf(debug, "Will call MPI_Finalize now\n");
816 ret = MPI_Finalize();
819 fprintf(debug, "Return code from MPI_Finalize = %d\n", ret);