2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
42 #include "gmx_fatal.h"
44 #include "gromacs/utility/smalloc.h"
45 #include "types/commrec.h"
52 #include "gromacs/utility/gmxmpi.h"
55 /* The source code in this file should be thread-safe.
56 Please keep it that way. */
58 gmx_bool gmx_mpi_initialized(void)
70 void gmx_fill_commrec_from_mpi(t_commrec gmx_unused *cr)
73 gmx_call("gmx_fill_commrec_from_mpi");
75 if (!gmx_mpi_initialized())
77 gmx_comm("MPI has not been initialized properly");
80 cr->nnodes = gmx_node_num();
81 cr->nodeid = gmx_node_rank();
82 cr->sim_nodeid = cr->nodeid;
83 cr->mpi_comm_mysim = MPI_COMM_WORLD;
84 cr->mpi_comm_mygroup = MPI_COMM_WORLD;
89 t_commrec *init_commrec()
96 gmx_fill_commrec_from_mpi(cr);
98 cr->mpi_comm_mysim = NULL;
99 cr->mpi_comm_mygroup = NULL;
102 cr->nodeid = cr->sim_nodeid;
105 // TODO cr->duty should not be initialized here
106 cr->duty = (DUTY_PP | DUTY_PME);
108 #if defined GMX_MPI && !defined MPI_IN_PLACE_EXISTS
109 /* initialize the MPI_IN_PLACE replacement buffers */
111 cr->mpb->ibuf = NULL;
112 cr->mpb->libuf = NULL;
113 cr->mpb->fbuf = NULL;
114 cr->mpb->dbuf = NULL;
115 cr->mpb->ibuf_alloc = 0;
116 cr->mpb->libuf_alloc = 0;
117 cr->mpb->fbuf_alloc = 0;
118 cr->mpb->dbuf_alloc = 0;
124 t_commrec *reinitialize_commrec_for_this_thread(const t_commrec gmx_unused *cro)
126 #ifdef GMX_THREAD_MPI
129 /* make a thread-specific commrec */
131 /* now copy the whole thing, so settings like the number of PME nodes
135 /* and we start setting our own thread-specific values for things */
136 gmx_fill_commrec_from_mpi(cr);
138 // TODO cr->duty should not be initialized here
139 cr->duty = (DUTY_PP | DUTY_PME);
147 int gmx_node_num(void)
153 (void) MPI_Comm_size(MPI_COMM_WORLD, &i);
158 int gmx_node_rank(void)
164 (void) MPI_Comm_rank(MPI_COMM_WORLD, &i);
169 #if defined GMX_LIB_MPI && defined GMX_TARGET_BGQ
170 #include <spi/include/kernel/location.h>
173 int gmx_physicalnode_id_hash(void)
178 /* We have a single physical node */
182 char mpi_hostname[MPI_MAX_PROCESSOR_NAME];
184 /* This procedure can only differentiate nodes with different names.
185 * Architectures where different physical nodes have identical names,
186 * such as IBM Blue Gene, should use an architecture specific solution.
188 MPI_Get_processor_name(mpi_hostname, &resultlen);
190 /* The string hash function returns an unsigned int. We cast to an int.
191 * Negative numbers are converted to positive by setting the sign bit to 0.
192 * This makes the hash one bit smaller.
193 * A 63-bit hash (with 64-bit int) should be enough for unique node hashes,
194 * even on a million node machine. 31 bits might not be enough though!
197 (int)gmx_string_fullhash_func(mpi_hostname, gmx_string_hash_init);
207 /* TODO: this function should be fully replaced by gmx_physicalnode_id_hash */
208 int gmx_hostname_num()
213 #ifdef GMX_THREAD_MPI
214 /* thread-MPI currently puts the thread number in the process name,
215 * we might want to change this, as this is inconsistent with what
216 * most MPI implementations would do when running on a single node.
220 int resultlen, hostnum, i, j;
221 char mpi_hostname[MPI_MAX_PROCESSOR_NAME], hostnum_str[MPI_MAX_PROCESSOR_NAME];
223 MPI_Get_processor_name(mpi_hostname, &resultlen);
224 #ifdef GMX_TARGET_BGQ
225 Personality_t personality;
226 Kernel_GetPersonality(&personality, sizeof(personality));
227 /* Each MPI rank has a unique coordinate in a 6-dimensional space
228 (A,B,C,D,E,T), with dimensions A-E corresponding to different
229 physical nodes, and T within each node. Each node has sixteen
230 physical cores, each of which can have up to four hardware
231 threads, so 0 <= T <= 63 (but the maximum value of T depends on
232 the confituration of ranks and OpenMP threads per
233 node). However, T is irrelevant for computing a suitable return
234 value for gmx_hostname_num().
236 hostnum = personality.Network_Config.Acoord;
237 hostnum *= personality.Network_Config.Bnodes;
238 hostnum += personality.Network_Config.Bcoord;
239 hostnum *= personality.Network_Config.Cnodes;
240 hostnum += personality.Network_Config.Ccoord;
241 hostnum *= personality.Network_Config.Dnodes;
242 hostnum += personality.Network_Config.Dcoord;
243 hostnum *= personality.Network_Config.Enodes;
244 hostnum += personality.Network_Config.Ecoord;
246 /* This procedure can only differentiate nodes with host names
247 * that end on unique numbers.
251 /* Only parse the host name up to the first dot */
252 while (i < resultlen && mpi_hostname[i] != '.')
254 if (isdigit(mpi_hostname[i]))
256 hostnum_str[j++] = mpi_hostname[i];
260 hostnum_str[j] = '\0';
267 /* Use only the last 9 decimals, so we don't overflow an int */
268 hostnum = strtol(hostnum_str + max(0, j-9), NULL, 10);
274 fprintf(debug, "In gmx_hostname_num: hostname '%s', hostnum %d\n",
275 mpi_hostname, hostnum);
276 #ifdef GMX_TARGET_BGQ
278 "Torus ID A: %d / %d B: %d / %d C: %d / %d D: %d / %d E: %d / %d\nNode ID T: %d / %d core: %d / %d hardware thread: %d / %d\n",
279 personality.Network_Config.Acoord,
280 personality.Network_Config.Anodes,
281 personality.Network_Config.Bcoord,
282 personality.Network_Config.Bnodes,
283 personality.Network_Config.Ccoord,
284 personality.Network_Config.Cnodes,
285 personality.Network_Config.Dcoord,
286 personality.Network_Config.Dnodes,
287 personality.Network_Config.Ecoord,
288 personality.Network_Config.Enodes,
289 Kernel_ProcessorCoreID(),
291 Kernel_ProcessorID(),
293 Kernel_ProcessorThreadID(),
302 void gmx_setup_nodecomm(FILE gmx_unused *fplog, t_commrec *cr)
305 int n, rank, hostnum, ng, ni;
307 /* Many MPI implementations do not optimize MPI_Allreduce
308 * (and probably also other global communication calls)
309 * for multi-core nodes connected by a network.
310 * We can optimize such communication by using one MPI call
311 * within each node and one between the nodes.
312 * For MVAPICH2 and Intel MPI this reduces the time for
313 * the global_stat communication by 25%
314 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
315 * B. Hess, November 2007
321 #ifndef GMX_THREAD_MPI
323 MPI_Comm_size(cr->mpi_comm_mygroup, &n);
324 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
326 hostnum = gmx_hostname_num();
330 fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
334 /* The intra-node communicator, split on node number */
335 MPI_Comm_split(cr->mpi_comm_mygroup, hostnum, rank, &nc->comm_intra);
336 MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
339 fprintf(debug, "In gmx_setup_nodecomm: node rank %d rank_intra %d\n",
340 rank, nc->rank_intra);
342 /* The inter-node communicator, split on rank_intra.
343 * We actually only need the one for rank=0,
344 * but it is easier to create them all.
346 MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
347 /* Check if this really created two step communication */
348 MPI_Comm_size(nc->comm_inter, &ng);
349 MPI_Comm_size(nc->comm_intra, &ni);
352 fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
356 if (getenv("GMX_NO_NODECOMM") == NULL &&
357 ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
362 fprintf(fplog, "Using two step summing over %d groups of on average %.1f processes\n\n",
363 ng, (real)n/(real)ng);
365 if (nc->rank_intra > 0)
367 MPI_Comm_free(&nc->comm_inter);
372 /* One group or all processes in a separate group, use normal summing */
373 MPI_Comm_free(&nc->comm_inter);
374 MPI_Comm_free(&nc->comm_intra);
377 fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
382 /* tMPI runs only on a single node so just use the nodeid */
383 nc->rank_intra = cr->nodeid;
387 void gmx_init_intranode_counters(t_commrec *cr)
389 /* counters for PP+PME and PP-only processes on my physical node */
390 int nrank_intranode, rank_intranode;
391 int nrank_pp_intranode, rank_pp_intranode;
392 /* thread-MPI is not initialized when not running in parallel */
393 #if defined GMX_MPI && !defined GMX_THREAD_MPI
394 int nrank_world, rank_world;
395 int i, mynum, *num, *num_s, *num_pp, *num_pp_s;
397 MPI_Comm_size(MPI_COMM_WORLD, &nrank_world);
398 MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
400 /* Get the node number from the hostname to identify the nodes */
401 mynum = gmx_hostname_num();
403 /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
404 snew(num, nrank_world);
405 snew(num_s, nrank_world);
406 snew(num_pp, nrank_world);
407 snew(num_pp_s, nrank_world);
409 num_s[rank_world] = mynum;
410 num_pp_s[rank_world] = (cr->duty & DUTY_PP) ? mynum : -1;
412 MPI_Allreduce(num_s, num, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
413 MPI_Allreduce(num_pp_s, num_pp, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
417 nrank_pp_intranode = 0;
418 rank_pp_intranode = 0;
419 for (i = 0; i < nrank_world; i++)
429 if (num_pp[i] == mynum)
431 nrank_pp_intranode++;
432 if ((cr->duty & DUTY_PP) && i < rank_world)
443 /* Serial or thread-MPI code: we run within a single physical node */
444 nrank_intranode = cr->nnodes;
445 rank_intranode = cr->sim_nodeid;
446 nrank_pp_intranode = cr->nnodes - cr->npmenodes;
447 rank_pp_intranode = cr->nodeid;
453 if (cr->duty & DUTY_PP && cr->duty & DUTY_PME)
455 sprintf(sbuf, "PP+PME");
459 sprintf(sbuf, "%s", cr->duty & DUTY_PP ? "PP" : "PME");
461 fprintf(debug, "On %3s node %d: nrank_intranode=%d, rank_intranode=%d, "
462 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
463 sbuf, cr->sim_nodeid,
464 nrank_intranode, rank_intranode,
465 nrank_pp_intranode, rank_pp_intranode);
468 cr->nrank_intranode = nrank_intranode;
469 cr->rank_intranode = rank_intranode;
470 cr->nrank_pp_intranode = nrank_pp_intranode;
471 cr->rank_pp_intranode = rank_pp_intranode;
475 void gmx_barrier(const t_commrec gmx_unused *cr)
478 gmx_call("gmx_barrier");
480 MPI_Barrier(cr->mpi_comm_mygroup);
484 void gmx_abort(int gmx_unused noderank, int gmx_unused nnodes, int gmx_unused errorno)
487 gmx_call("gmx_abort");
489 #ifdef GMX_THREAD_MPI
490 fprintf(stderr, "Halting program %s\n", ShortProgram());
496 fprintf(stderr, "Halting parallel program %s on CPU %d out of %d\n",
497 ShortProgram(), noderank, nnodes);
501 fprintf(stderr, "Halting program %s\n", ShortProgram());
505 MPI_Abort(MPI_COMM_WORLD, errorno);
511 void gmx_bcast(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
514 gmx_call("gmx_bast");
516 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
520 void gmx_bcast_sim(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
523 gmx_call("gmx_bast");
525 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
529 void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unused *cr)
532 gmx_call("gmx_sumd");
534 #if defined(MPI_IN_PLACE_EXISTS)
537 if (cr->nc.rank_intra == 0)
539 /* Use two step summing. */
540 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
542 /* Sum the roots of the internal (intra) buffers. */
543 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
548 /* This is here because of the silly MPI specification
549 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
550 MPI_Reduce(r, NULL, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
552 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
556 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
557 cr->mpi_comm_mygroup);
562 if (nr > cr->mpb->dbuf_alloc)
564 cr->mpb->dbuf_alloc = nr;
565 srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
569 /* Use two step summing */
570 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
571 if (cr->nc.rank_intra == 0)
573 /* Sum with the buffers reversed */
574 MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
577 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
581 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
582 cr->mpi_comm_mygroup);
583 for (i = 0; i < nr; i++)
585 r[i] = cr->mpb->dbuf[i];
592 void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unused *cr)
595 gmx_call("gmx_sumf");
597 #if defined(MPI_IN_PLACE_EXISTS)
600 /* Use two step summing. */
601 if (cr->nc.rank_intra == 0)
603 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
605 /* Sum the roots of the internal (intra) buffers */
606 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
611 /* This is here because of the silly MPI specification
612 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
613 MPI_Reduce(r, NULL, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
615 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
619 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
624 if (nr > cr->mpb->fbuf_alloc)
626 cr->mpb->fbuf_alloc = nr;
627 srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
631 /* Use two step summing */
632 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
633 if (cr->nc.rank_intra == 0)
635 /* Sum with the buffers reversed */
636 MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
639 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
643 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
644 cr->mpi_comm_mygroup);
645 for (i = 0; i < nr; i++)
647 r[i] = cr->mpb->fbuf[i];
654 void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused *cr)
657 gmx_call("gmx_sumi");
659 #if defined(MPI_IN_PLACE_EXISTS)
662 /* Use two step summing */
663 if (cr->nc.rank_intra == 0)
665 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
666 /* Sum with the buffers reversed */
667 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
671 /* This is here because of the silly MPI specification
672 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
673 MPI_Reduce(r, NULL, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
675 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
679 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
684 if (nr > cr->mpb->ibuf_alloc)
686 cr->mpb->ibuf_alloc = nr;
687 srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
691 /* Use two step summing */
692 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
693 if (cr->nc.rank_intra == 0)
695 /* Sum with the buffers reversed */
696 MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
698 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
702 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
703 for (i = 0; i < nr; i++)
705 r[i] = cr->mpb->ibuf[i];
712 void gmx_sumli(int gmx_unused nr, gmx_int64_t gmx_unused r[], const t_commrec gmx_unused *cr)
715 gmx_call("gmx_sumli");
717 #if defined(MPI_IN_PLACE_EXISTS)
720 /* Use two step summing */
721 if (cr->nc.rank_intra == 0)
723 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, 0,
725 /* Sum with the buffers reversed */
726 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
731 /* This is here because of the silly MPI specification
732 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
733 MPI_Reduce(r, NULL, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
735 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
739 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
744 if (nr > cr->mpb->libuf_alloc)
746 cr->mpb->libuf_alloc = nr;
747 srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
751 /* Use two step summing */
752 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
754 if (cr->nc.rank_intra == 0)
756 /* Sum with the buffers reversed */
757 MPI_Allreduce(cr->mpb->libuf, r, nr, MPI_INT64_T, MPI_SUM,
760 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
764 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
765 cr->mpi_comm_mygroup);
766 for (i = 0; i < nr; i++)
768 r[i] = cr->mpb->libuf[i];
778 static void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
780 #if defined(MPI_IN_PLACE_EXISTS)
781 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
783 /* this function is only used in code that is not performance critical,
784 (during setup, when comm_rec is not the appropriate communication
785 structure), so this isn't as bad as it looks. */
790 MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
791 for (i = 0; i < nr; i++)
801 static void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
803 #if defined(MPI_IN_PLACE_EXISTS)
804 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
806 /* this function is only used in code that is not performance critical,
807 (during setup, when comm_rec is not the appropriate communication
808 structure), so this isn't as bad as it looks. */
813 MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
814 for (i = 0; i < nr; i++)
823 void gmx_sumd_sim(int gmx_unused nr, double gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
826 gmx_call("gmx_sumd_sim");
828 gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
832 void gmx_sumf_sim(int gmx_unused nr, float gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
835 gmx_call("gmx_sumf_sim");
837 gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
841 void gmx_sumi_sim(int gmx_unused nr, int gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
844 gmx_call("gmx_sumi_sim");
846 #if defined(MPI_IN_PLACE_EXISTS)
847 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
849 /* this is thread-unsafe, but it will do for now: */
852 if (nr > ms->mpb->ibuf_alloc)
854 ms->mpb->ibuf_alloc = nr;
855 srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
857 MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
858 for (i = 0; i < nr; i++)
860 r[i] = ms->mpb->ibuf[i];
866 void gmx_sumli_sim(int gmx_unused nr, gmx_int64_t gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
869 gmx_call("gmx_sumli_sim");
871 #if defined(MPI_IN_PLACE_EXISTS)
872 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
873 ms->mpi_comm_masters);
875 /* this is thread-unsafe, but it will do for now: */
878 if (nr > ms->mpb->libuf_alloc)
880 ms->mpb->libuf_alloc = nr;
881 srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
883 MPI_Allreduce(r, ms->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
884 ms->mpi_comm_masters);
885 for (i = 0; i < nr; i++)
887 r[i] = ms->mpb->libuf[i];