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"
58 /* The source code in this file should be thread-safe.
59 Please keep it that way. */
61 gmx_bool gmx_mpi_initialized(void)
73 int gmx_setup(int *argc, char **argv, int *nnodes)
76 gmx_call("gmx_setup");
80 int resultlen; /* actual length of node name */
84 char mpi_hostname[MPI_MAX_PROCESSOR_NAME];
86 /* Call the MPI routines */
89 (void) fah_MPI_Init(argc, &argv);
91 (void) MPI_Init(argc, &argv);
94 (void) MPI_Comm_size( MPI_COMM_WORLD, &mpi_num_nodes );
95 (void) MPI_Comm_rank( MPI_COMM_WORLD, &mpi_my_rank );
96 (void) MPI_Get_processor_name( mpi_hostname, &resultlen );
101 fprintf(debug, "NNODES=%d, MYRANK=%d, HOSTNAME=%s\n",
102 mpi_num_nodes, mpi_my_rank, mpi_hostname);
106 *nnodes = mpi_num_nodes;
112 int gmx_node_num(void)
118 (void) MPI_Comm_size(MPI_COMM_WORLD, &i);
123 int gmx_node_rank(void)
129 (void) MPI_Comm_rank(MPI_COMM_WORLD, &i);
135 int gmx_hostname_num()
140 #ifdef GMX_THREAD_MPI
141 /* thread-MPI currently puts the thread number in the process name,
142 * we might want to change this, as this is inconsistent with what
143 * most MPI implementations would do when running on a single node.
147 int resultlen, hostnum, i, j;
148 char mpi_hostname[MPI_MAX_PROCESSOR_NAME], hostnum_str[MPI_MAX_PROCESSOR_NAME];
150 MPI_Get_processor_name(mpi_hostname, &resultlen);
151 /* This procedure can only differentiate nodes with host names
152 * that end on unique numbers.
156 /* Only parse the host name up to the first dot */
157 while (i < resultlen && mpi_hostname[i] != '.')
159 if (isdigit(mpi_hostname[i]))
161 hostnum_str[j++] = mpi_hostname[i];
165 hostnum_str[j] = '\0';
172 /* Use only the last 9 decimals, so we don't overflow an int */
173 hostnum = strtol(hostnum_str + max(0, j-9), NULL, 10);
178 fprintf(debug, "In gmx_setup_nodecomm: hostname '%s', hostnum %d\n",
179 mpi_hostname, hostnum);
186 void gmx_setup_nodecomm(FILE *fplog, t_commrec *cr)
189 int n, rank, hostnum, ng, ni;
191 /* Many MPI implementations do not optimize MPI_Allreduce
192 * (and probably also other global communication calls)
193 * for multi-core nodes connected by a network.
194 * We can optimize such communication by using one MPI call
195 * within each node and one between the nodes.
196 * For MVAPICH2 and Intel MPI this reduces the time for
197 * the global_stat communication by 25%
198 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
199 * B. Hess, November 2007
205 #ifndef GMX_THREAD_MPI
207 MPI_Comm_size(cr->mpi_comm_mygroup, &n);
208 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
210 hostnum = gmx_hostname_num();
214 fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
218 /* The intra-node communicator, split on node number */
219 MPI_Comm_split(cr->mpi_comm_mygroup, hostnum, rank, &nc->comm_intra);
220 MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
223 fprintf(debug, "In gmx_setup_nodecomm: node rank %d rank_intra %d\n",
224 rank, nc->rank_intra);
226 /* The inter-node communicator, split on rank_intra.
227 * We actually only need the one for rank=0,
228 * but it is easier to create them all.
230 MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
231 /* Check if this really created two step communication */
232 MPI_Comm_size(nc->comm_inter, &ng);
233 MPI_Comm_size(nc->comm_intra, &ni);
236 fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
240 if (getenv("GMX_NO_NODECOMM") == NULL &&
241 ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
246 fprintf(fplog, "Using two step summing over %d groups of on average %.1f processes\n\n",
247 ng, (real)n/(real)ng);
249 if (nc->rank_intra > 0)
251 MPI_Comm_free(&nc->comm_inter);
256 /* One group or all processes in a separate group, use normal summing */
257 MPI_Comm_free(&nc->comm_inter);
258 MPI_Comm_free(&nc->comm_intra);
261 fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
266 /* tMPI runs only on a single node so just use the nodeid */
267 nc->rank_intra = cr->nodeid;
271 void gmx_init_intranode_counters(t_commrec *cr)
273 /* counters for PP+PME and PP-only processes on my physical node */
274 int nrank_intranode, rank_intranode;
275 int nrank_pp_intranode, rank_pp_intranode;
276 /* thread-MPI is not initialized when not running in parallel */
277 #if defined GMX_MPI && !defined GMX_THREAD_MPI
278 int nrank_world, rank_world;
279 int i, mynum, *num, *num_s, *num_pp, *num_pp_s;
281 MPI_Comm_size(MPI_COMM_WORLD, &nrank_world);
282 MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
284 /* Get the node number from the hostname to identify the nodes */
285 mynum = gmx_hostname_num();
287 /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
288 snew(num, nrank_world);
289 snew(num_s, nrank_world);
290 snew(num_pp, nrank_world);
291 snew(num_pp_s, nrank_world);
293 num_s[rank_world] = mynum;
294 num_pp_s[rank_world] = (cr->duty & DUTY_PP) ? mynum : -1;
296 MPI_Allreduce(num_s, num, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
297 MPI_Allreduce(num_pp_s, num_pp, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
301 nrank_pp_intranode = 0;
302 rank_pp_intranode = 0;
303 for (i = 0; i < nrank_world; i++)
313 if ((cr->duty & DUTY_PP) && num_pp[i] == mynum)
315 nrank_pp_intranode++;
327 /* Serial or thread-MPI code: we run within a single physical node */
328 nrank_intranode = cr->nnodes;
329 rank_intranode = cr->sim_nodeid;
330 nrank_pp_intranode = cr->nnodes - cr->npmenodes;
331 rank_pp_intranode = cr->nodeid;
337 if (cr->duty & DUTY_PP && cr->duty & DUTY_PME)
339 sprintf(sbuf, "PP+PME");
343 sprintf(sbuf, "%s", cr->duty & DUTY_PP ? "PP" : "PME");
345 fprintf(debug, "On %3s node %d: nrank_intranode=%d, rank_intranode=%d, "
346 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
347 sbuf, cr->sim_nodeid,
348 nrank_intranode, rank_intranode,
349 nrank_pp_intranode, rank_pp_intranode);
352 cr->nrank_intranode = nrank_intranode;
353 cr->rank_intranode = rank_intranode;
354 cr->nrank_pp_intranode = nrank_pp_intranode;
355 cr->rank_pp_intranode = rank_pp_intranode;
359 void gmx_barrier(const t_commrec *cr)
362 gmx_call("gmx_barrier");
364 MPI_Barrier(cr->mpi_comm_mygroup);
368 void gmx_abort(int noderank, int nnodes, int errorno)
371 gmx_call("gmx_abort");
373 #ifdef GMX_THREAD_MPI
374 fprintf(stderr, "Halting program %s\n", ShortProgram());
380 fprintf(stderr, "Halting parallel program %s on CPU %d out of %d\n",
381 ShortProgram(), noderank, nnodes);
385 fprintf(stderr, "Halting program %s\n", ShortProgram());
389 MPI_Abort(MPI_COMM_WORLD, errorno);
395 void gmx_bcast(int nbytes, void *b, const t_commrec *cr)
398 gmx_call("gmx_bast");
400 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
404 void gmx_bcast_sim(int nbytes, void *b, const t_commrec *cr)
407 gmx_call("gmx_bast");
409 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
413 void gmx_sumd(int nr, double r[], const t_commrec *cr)
416 gmx_call("gmx_sumd");
418 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
421 if (cr->nc.rank_intra == 0)
423 /* Use two step summing. */
424 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
426 /* Sum the roots of the internal (intra) buffers. */
427 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
432 /* This is here because of the silly MPI specification
433 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
434 MPI_Reduce(r, NULL, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
436 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
440 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
441 cr->mpi_comm_mygroup);
446 if (nr > cr->mpb->dbuf_alloc)
448 cr->mpb->dbuf_alloc = nr;
449 srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
453 /* Use two step summing */
454 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
455 if (cr->nc.rank_intra == 0)
457 /* Sum with the buffers reversed */
458 MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
461 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
465 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
466 cr->mpi_comm_mygroup);
467 for (i = 0; i < nr; i++)
469 r[i] = cr->mpb->dbuf[i];
476 void gmx_sumf(int nr, float r[], const t_commrec *cr)
479 gmx_call("gmx_sumf");
481 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
484 /* Use two step summing. */
485 if (cr->nc.rank_intra == 0)
487 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
489 /* Sum the roots of the internal (intra) buffers */
490 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
495 /* This is here because of the silly MPI specification
496 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
497 MPI_Reduce(r, NULL, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
499 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
503 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
508 if (nr > cr->mpb->fbuf_alloc)
510 cr->mpb->fbuf_alloc = nr;
511 srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
515 /* Use two step summing */
516 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
517 if (cr->nc.rank_intra == 0)
519 /* Sum with the buffers reversed */
520 MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
523 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
527 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
528 cr->mpi_comm_mygroup);
529 for (i = 0; i < nr; i++)
531 r[i] = cr->mpb->fbuf[i];
538 void gmx_sumi(int nr, int r[], const t_commrec *cr)
541 gmx_call("gmx_sumi");
543 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
546 /* Use two step summing */
547 if (cr->nc.rank_intra == 0)
549 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
550 /* Sum with the buffers reversed */
551 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
555 /* This is here because of the silly MPI specification
556 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
557 MPI_Reduce(r, NULL, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
559 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
563 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
568 if (nr > cr->mpb->ibuf_alloc)
570 cr->mpb->ibuf_alloc = nr;
571 srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
575 /* Use two step summing */
576 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
577 if (cr->nc.rank_intra == 0)
579 /* Sum with the buffers reversed */
580 MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
582 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
586 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
587 for (i = 0; i < nr; i++)
589 r[i] = cr->mpb->ibuf[i];
596 void gmx_sumli(int nr, gmx_large_int_t r[], const t_commrec *cr)
599 gmx_call("gmx_sumli");
601 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
604 /* Use two step summing */
605 if (cr->nc.rank_intra == 0)
607 MPI_Reduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM, 0,
609 /* Sum with the buffers reversed */
610 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
615 /* This is here because of the silly MPI specification
616 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
617 MPI_Reduce(r, NULL, nr, GMX_MPI_LARGE_INT, MPI_SUM, 0, cr->nc.comm_intra);
619 MPI_Bcast(r, nr, GMX_MPI_LARGE_INT, 0, cr->nc.comm_intra);
623 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM, cr->mpi_comm_mygroup);
628 if (nr > cr->mpb->libuf_alloc)
630 cr->mpb->libuf_alloc = nr;
631 srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
635 /* Use two step summing */
636 MPI_Allreduce(r, cr->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
638 if (cr->nc.rank_intra == 0)
640 /* Sum with the buffers reversed */
641 MPI_Allreduce(cr->mpb->libuf, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
644 MPI_Bcast(r, nr, GMX_MPI_LARGE_INT, 0, cr->nc.comm_intra);
648 MPI_Allreduce(r, cr->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
649 cr->mpi_comm_mygroup);
650 for (i = 0; i < nr; i++)
652 r[i] = cr->mpb->libuf[i];
662 void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
664 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
665 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
667 /* this function is only used in code that is not performance critical,
668 (during setup, when comm_rec is not the appropriate communication
669 structure), so this isn't as bad as it looks. */
674 MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
675 for (i = 0; i < nr; i++)
685 void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
687 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
688 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
690 /* this function is only used in code that is not performance critical,
691 (during setup, when comm_rec is not the appropriate communication
692 structure), so this isn't as bad as it looks. */
697 MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
698 for (i = 0; i < nr; i++)
707 void gmx_sumd_sim(int nr, double r[], const gmx_multisim_t *ms)
710 gmx_call("gmx_sumd_sim");
712 gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
716 void gmx_sumf_sim(int nr, float r[], const gmx_multisim_t *ms)
719 gmx_call("gmx_sumf_sim");
721 gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
725 void gmx_sumi_sim(int nr, int r[], const gmx_multisim_t *ms)
728 gmx_call("gmx_sumi_sim");
730 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
731 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
733 /* this is thread-unsafe, but it will do for now: */
736 if (nr > ms->mpb->ibuf_alloc)
738 ms->mpb->ibuf_alloc = nr;
739 srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
741 MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
742 for (i = 0; i < nr; i++)
744 r[i] = ms->mpb->ibuf[i];
750 void gmx_sumli_sim(int nr, gmx_large_int_t r[], const gmx_multisim_t *ms)
753 gmx_call("gmx_sumli_sim");
755 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
756 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
757 ms->mpi_comm_masters);
759 /* this is thread-unsafe, but it will do for now: */
762 if (nr > ms->mpb->libuf_alloc)
764 ms->mpb->libuf_alloc = nr;
765 srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
767 MPI_Allreduce(r, ms->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
768 ms->mpi_comm_masters);
769 for (i = 0; i < nr; i++)
771 r[i] = ms->mpb->libuf[i];
778 void gmx_finalize_par(void)
781 /* Compiled without MPI, no MPI finalizing needed */
784 int initialized, finalized;
787 MPI_Initialized(&initialized);
792 /* just as a check; we don't want to finalize twice */
793 MPI_Finalized(&finalized);
799 /* We sync the processes here to try to avoid problems
800 * with buggy MPI implementations that could cause
801 * unfinished processes to terminate.
803 MPI_Barrier(MPI_COMM_WORLD);
806 if (DOMAINDECOMP(cr)) {
807 if (cr->npmenodes > 0 || cr->dd->bCartesian)
808 MPI_Comm_free(&cr->mpi_comm_mygroup);
809 if (cr->dd->bCartesian)
810 MPI_Comm_free(&cr->mpi_comm_mysim);
814 /* Apparently certain mpich implementations cause problems
815 * with MPI_Finalize. In that case comment out MPI_Finalize.
819 fprintf(debug, "Will call MPI_Finalize now\n");
822 ret = MPI_Finalize();
825 fprintf(debug, "Return code from MPI_Finalize = %d\n", ret);