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"
50 #include "gromacs/utility/cstringutil.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 static int mpi_hostname_hash(void)
174 /* We have a single physical node */
178 char mpi_hostname[MPI_MAX_PROCESSOR_NAME];
180 /* This procedure can only differentiate nodes with different names.
181 * Architectures where different physical nodes have identical names,
182 * such as IBM Blue Gene, should use an architecture specific solution.
184 MPI_Get_processor_name(mpi_hostname, &resultlen);
186 /* The string hash function returns an unsigned int. We cast to an int.
187 * Negative numbers are converted to positive by setting the sign bit to 0.
188 * This makes the hash one bit smaller.
189 * A 63-bit hash (with 64-bit int) should be enough for unique node hashes,
190 * even on a million node machine. 31 bits might not be enough though!
193 (int)gmx_string_fullhash_func(mpi_hostname, gmx_string_hash_init);
203 #if defined GMX_LIB_MPI && defined GMX_TARGET_BGQ
206 /* IBM's declaration of this function in
207 * /bgsys/drivers/V1R2M2/ppc64/spi/include/kernel/process.h
208 * erroneously fails to specify __INLINE__, despite
209 * /bgsys/drivers/V1R2M2/ppc64/spi/include/kernel/cnk/process_impl.h
210 * specifiying __INLINE__, so bgclang thinks they are different enough
211 * to complain about. */
212 static uint64_t Kernel_GetJobID();
215 #include <spi/include/kernel/location.h>
217 static int bgq_nodenum(void)
220 Personality_t personality;
221 Kernel_GetPersonality(&personality, sizeof(personality));
222 /* Each MPI rank has a unique coordinate in a 6-dimensional space
223 (A,B,C,D,E,T), with dimensions A-E corresponding to different
224 physical nodes, and T within each node. Each node has sixteen
225 physical cores, each of which can have up to four hardware
226 threads, so 0 <= T <= 63 (but the maximum value of T depends on
227 the confituration of ranks and OpenMP threads per
228 node). However, T is irrelevant for computing a suitable return
229 value for gmx_hostname_num().
231 hostnum = personality.Network_Config.Acoord;
232 hostnum *= personality.Network_Config.Bnodes;
233 hostnum += personality.Network_Config.Bcoord;
234 hostnum *= personality.Network_Config.Cnodes;
235 hostnum += personality.Network_Config.Ccoord;
236 hostnum *= personality.Network_Config.Dnodes;
237 hostnum += personality.Network_Config.Dcoord;
238 hostnum *= personality.Network_Config.Enodes;
239 hostnum += personality.Network_Config.Ecoord;
244 "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",
245 personality.Network_Config.Acoord,
246 personality.Network_Config.Anodes,
247 personality.Network_Config.Bcoord,
248 personality.Network_Config.Bnodes,
249 personality.Network_Config.Ccoord,
250 personality.Network_Config.Cnodes,
251 personality.Network_Config.Dcoord,
252 personality.Network_Config.Dnodes,
253 personality.Network_Config.Ecoord,
254 personality.Network_Config.Enodes,
255 Kernel_ProcessorCoreID(),
257 Kernel_ProcessorID(),
259 Kernel_ProcessorThreadID(),
266 int gmx_physicalnode_id_hash(void)
273 #ifdef GMX_THREAD_MPI
274 /* thread-MPI currently puts the thread number in the process name,
275 * we might want to change this, as this is inconsistent with what
276 * most MPI implementations would do when running on a single node.
280 #ifdef GMX_TARGET_BGQ
281 hash = bgq_nodenum();
283 hash = mpi_hostname_hash();
290 fprintf(debug, "In gmx_physicalnode_id_hash: hash %d\n", hash);
296 void gmx_setup_nodecomm(FILE gmx_unused *fplog, t_commrec *cr)
299 int n, rank, nodehash, ng, ni;
301 /* Many MPI implementations do not optimize MPI_Allreduce
302 * (and probably also other global communication calls)
303 * for multi-core nodes connected by a network.
304 * We can optimize such communication by using one MPI call
305 * within each node and one between the nodes.
306 * For MVAPICH2 and Intel MPI this reduces the time for
307 * the global_stat communication by 25%
308 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
309 * B. Hess, November 2007
315 #ifndef GMX_THREAD_MPI
317 MPI_Comm_size(cr->mpi_comm_mygroup, &n);
318 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
320 nodehash = gmx_physicalnode_id_hash();
324 fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
328 /* The intra-node communicator, split on node number */
329 MPI_Comm_split(cr->mpi_comm_mygroup, nodehash, rank, &nc->comm_intra);
330 MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
333 fprintf(debug, "In gmx_setup_nodecomm: node ID %d rank within node %d\n",
334 rank, nc->rank_intra);
336 /* The inter-node communicator, split on rank_intra.
337 * We actually only need the one for rank=0,
338 * but it is easier to create them all.
340 MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
341 /* Check if this really created two step communication */
342 MPI_Comm_size(nc->comm_inter, &ng);
343 MPI_Comm_size(nc->comm_intra, &ni);
346 fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
350 if (getenv("GMX_NO_NODECOMM") == NULL &&
351 ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
356 fprintf(fplog, "Using two step summing over %d groups of on average %.1f ranks\n\n",
357 ng, (real)n/(real)ng);
359 if (nc->rank_intra > 0)
361 MPI_Comm_free(&nc->comm_inter);
366 /* One group or all processes in a separate group, use normal summing */
367 MPI_Comm_free(&nc->comm_inter);
368 MPI_Comm_free(&nc->comm_intra);
371 fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
376 /* tMPI runs only on a single node so just use the nodeid */
377 nc->rank_intra = cr->nodeid;
381 void gmx_init_intranode_counters(t_commrec *cr)
383 /* counters for PP+PME and PP-only processes on my physical node */
384 int nrank_intranode, rank_intranode;
385 int nrank_pp_intranode, rank_pp_intranode;
386 /* thread-MPI is not initialized when not running in parallel */
387 #if defined GMX_MPI && !defined GMX_THREAD_MPI
388 int nrank_world, rank_world;
389 int i, myhash, *hash, *hash_s, *hash_pp, *hash_pp_s;
391 MPI_Comm_size(MPI_COMM_WORLD, &nrank_world);
392 MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
394 /* Get a (hopefully unique) hash that identifies our physical node */
395 myhash = gmx_physicalnode_id_hash();
397 /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
398 snew(hash, nrank_world);
399 snew(hash_s, nrank_world);
400 snew(hash_pp, nrank_world);
401 snew(hash_pp_s, nrank_world);
403 hash_s[rank_world] = myhash;
404 hash_pp_s[rank_world] = (cr->duty & DUTY_PP) ? myhash : -1;
406 MPI_Allreduce(hash_s, hash, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
407 MPI_Allreduce(hash_pp_s, hash_pp, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
411 nrank_pp_intranode = 0;
412 rank_pp_intranode = 0;
413 for (i = 0; i < nrank_world; i++)
415 if (hash[i] == myhash)
423 if (hash_pp[i] == myhash)
425 nrank_pp_intranode++;
426 if ((cr->duty & DUTY_PP) && i < rank_world)
437 /* Serial or thread-MPI code: we run within a single physical node */
438 nrank_intranode = cr->nnodes;
439 rank_intranode = cr->sim_nodeid;
440 nrank_pp_intranode = cr->nnodes - cr->npmenodes;
441 rank_pp_intranode = cr->nodeid;
447 if (cr->duty & DUTY_PP && cr->duty & DUTY_PME)
449 sprintf(sbuf, "PP+PME");
453 sprintf(sbuf, "%s", cr->duty & DUTY_PP ? "PP" : "PME");
455 fprintf(debug, "On %3s rank %d: nrank_intranode=%d, rank_intranode=%d, "
456 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
457 sbuf, cr->sim_nodeid,
458 nrank_intranode, rank_intranode,
459 nrank_pp_intranode, rank_pp_intranode);
462 cr->nrank_intranode = nrank_intranode;
463 cr->rank_intranode = rank_intranode;
464 cr->nrank_pp_intranode = nrank_pp_intranode;
465 cr->rank_pp_intranode = rank_pp_intranode;
469 void gmx_barrier(const t_commrec gmx_unused *cr)
472 gmx_call("gmx_barrier");
474 MPI_Barrier(cr->mpi_comm_mygroup);
478 void gmx_abort(int gmx_unused noderank, int gmx_unused nnodes, int gmx_unused errorno)
481 gmx_call("gmx_abort");
483 #ifdef GMX_THREAD_MPI
484 fprintf(stderr, "Halting program %s\n", ShortProgram());
490 fprintf(stderr, "Halting parallel program %s on CPU %d out of %d\n",
491 ShortProgram(), noderank, nnodes);
495 fprintf(stderr, "Halting program %s\n", ShortProgram());
499 MPI_Abort(MPI_COMM_WORLD, errorno);
505 void gmx_bcast(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
508 gmx_call("gmx_bast");
510 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
514 void gmx_bcast_sim(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
517 gmx_call("gmx_bast");
519 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
523 void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unused *cr)
526 gmx_call("gmx_sumd");
528 #if defined(MPI_IN_PLACE_EXISTS)
531 if (cr->nc.rank_intra == 0)
533 /* Use two step summing. */
534 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
536 /* Sum the roots of the internal (intra) buffers. */
537 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
542 /* This is here because of the silly MPI specification
543 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
544 MPI_Reduce(r, NULL, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
546 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
550 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
551 cr->mpi_comm_mygroup);
556 if (nr > cr->mpb->dbuf_alloc)
558 cr->mpb->dbuf_alloc = nr;
559 srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
563 /* Use two step summing */
564 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
565 if (cr->nc.rank_intra == 0)
567 /* Sum with the buffers reversed */
568 MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
571 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
575 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
576 cr->mpi_comm_mygroup);
577 for (i = 0; i < nr; i++)
579 r[i] = cr->mpb->dbuf[i];
586 void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unused *cr)
589 gmx_call("gmx_sumf");
591 #if defined(MPI_IN_PLACE_EXISTS)
594 /* Use two step summing. */
595 if (cr->nc.rank_intra == 0)
597 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
599 /* Sum the roots of the internal (intra) buffers */
600 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
605 /* This is here because of the silly MPI specification
606 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
607 MPI_Reduce(r, NULL, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
609 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
613 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
618 if (nr > cr->mpb->fbuf_alloc)
620 cr->mpb->fbuf_alloc = nr;
621 srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
625 /* Use two step summing */
626 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
627 if (cr->nc.rank_intra == 0)
629 /* Sum with the buffers reversed */
630 MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
633 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
637 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
638 cr->mpi_comm_mygroup);
639 for (i = 0; i < nr; i++)
641 r[i] = cr->mpb->fbuf[i];
648 void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused *cr)
651 gmx_call("gmx_sumi");
653 #if defined(MPI_IN_PLACE_EXISTS)
656 /* Use two step summing */
657 if (cr->nc.rank_intra == 0)
659 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
660 /* Sum with the buffers reversed */
661 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
665 /* This is here because of the silly MPI specification
666 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
667 MPI_Reduce(r, NULL, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
669 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
673 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
678 if (nr > cr->mpb->ibuf_alloc)
680 cr->mpb->ibuf_alloc = nr;
681 srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
685 /* Use two step summing */
686 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
687 if (cr->nc.rank_intra == 0)
689 /* Sum with the buffers reversed */
690 MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
692 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
696 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
697 for (i = 0; i < nr; i++)
699 r[i] = cr->mpb->ibuf[i];
706 void gmx_sumli(int gmx_unused nr, gmx_int64_t gmx_unused r[], const t_commrec gmx_unused *cr)
709 gmx_call("gmx_sumli");
711 #if defined(MPI_IN_PLACE_EXISTS)
714 /* Use two step summing */
715 if (cr->nc.rank_intra == 0)
717 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, 0,
719 /* Sum with the buffers reversed */
720 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
725 /* This is here because of the silly MPI specification
726 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
727 MPI_Reduce(r, NULL, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
729 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
733 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
738 if (nr > cr->mpb->libuf_alloc)
740 cr->mpb->libuf_alloc = nr;
741 srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
745 /* Use two step summing */
746 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
748 if (cr->nc.rank_intra == 0)
750 /* Sum with the buffers reversed */
751 MPI_Allreduce(cr->mpb->libuf, r, nr, MPI_INT64_T, MPI_SUM,
754 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
758 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
759 cr->mpi_comm_mygroup);
760 for (i = 0; i < nr; i++)
762 r[i] = cr->mpb->libuf[i];
772 static void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
774 #if defined(MPI_IN_PLACE_EXISTS)
775 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
777 /* this function is only used in code that is not performance critical,
778 (during setup, when comm_rec is not the appropriate communication
779 structure), so this isn't as bad as it looks. */
784 MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
785 for (i = 0; i < nr; i++)
795 static void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
797 #if defined(MPI_IN_PLACE_EXISTS)
798 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
800 /* this function is only used in code that is not performance critical,
801 (during setup, when comm_rec is not the appropriate communication
802 structure), so this isn't as bad as it looks. */
807 MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
808 for (i = 0; i < nr; i++)
817 void gmx_sumd_sim(int gmx_unused nr, double gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
820 gmx_call("gmx_sumd_sim");
822 gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
826 void gmx_sumf_sim(int gmx_unused nr, float gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
829 gmx_call("gmx_sumf_sim");
831 gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
835 void gmx_sumi_sim(int gmx_unused nr, int gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
838 gmx_call("gmx_sumi_sim");
840 #if defined(MPI_IN_PLACE_EXISTS)
841 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
843 /* this is thread-unsafe, but it will do for now: */
846 if (nr > ms->mpb->ibuf_alloc)
848 ms->mpb->ibuf_alloc = nr;
849 srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
851 MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
852 for (i = 0; i < nr; i++)
854 r[i] = ms->mpb->ibuf[i];
860 void gmx_sumli_sim(int gmx_unused nr, gmx_int64_t gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
863 gmx_call("gmx_sumli_sim");
865 #if defined(MPI_IN_PLACE_EXISTS)
866 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
867 ms->mpi_comm_masters);
869 /* this is thread-unsafe, but it will do for now: */
872 if (nr > ms->mpb->libuf_alloc)
874 ms->mpb->libuf_alloc = nr;
875 srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
877 MPI_Allreduce(r, ms->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
878 ms->mpi_comm_masters);
879 for (i = 0; i < nr; i++)
881 r[i] = ms->mpb->libuf[i];