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
204 #include <spi/include/kernel/location.h>
206 static int bgq_nodenum(void)
209 Personality_t personality;
210 Kernel_GetPersonality(&personality, sizeof(personality));
211 /* Each MPI rank has a unique coordinate in a 6-dimensional space
212 (A,B,C,D,E,T), with dimensions A-E corresponding to different
213 physical nodes, and T within each node. Each node has sixteen
214 physical cores, each of which can have up to four hardware
215 threads, so 0 <= T <= 63 (but the maximum value of T depends on
216 the confituration of ranks and OpenMP threads per
217 node). However, T is irrelevant for computing a suitable return
218 value for gmx_hostname_num().
220 hostnum = personality.Network_Config.Acoord;
221 hostnum *= personality.Network_Config.Bnodes;
222 hostnum += personality.Network_Config.Bcoord;
223 hostnum *= personality.Network_Config.Cnodes;
224 hostnum += personality.Network_Config.Ccoord;
225 hostnum *= personality.Network_Config.Dnodes;
226 hostnum += personality.Network_Config.Dcoord;
227 hostnum *= personality.Network_Config.Enodes;
228 hostnum += personality.Network_Config.Ecoord;
233 "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",
234 personality.Network_Config.Acoord,
235 personality.Network_Config.Anodes,
236 personality.Network_Config.Bcoord,
237 personality.Network_Config.Bnodes,
238 personality.Network_Config.Ccoord,
239 personality.Network_Config.Cnodes,
240 personality.Network_Config.Dcoord,
241 personality.Network_Config.Dnodes,
242 personality.Network_Config.Ecoord,
243 personality.Network_Config.Enodes,
244 Kernel_ProcessorCoreID(),
246 Kernel_ProcessorID(),
248 Kernel_ProcessorThreadID(),
255 int gmx_physicalnode_id_hash(void)
262 #ifdef GMX_THREAD_MPI
263 /* thread-MPI currently puts the thread number in the process name,
264 * we might want to change this, as this is inconsistent with what
265 * most MPI implementations would do when running on a single node.
269 #ifdef GMX_TARGET_BGQ
270 hash = bgq_nodenum();
272 hash = mpi_hostname_hash();
279 fprintf(debug, "In gmx_physicalnode_id_hash: hash %d\n", hash);
285 void gmx_setup_nodecomm(FILE gmx_unused *fplog, t_commrec *cr)
288 int n, rank, nodehash, ng, ni;
290 /* Many MPI implementations do not optimize MPI_Allreduce
291 * (and probably also other global communication calls)
292 * for multi-core nodes connected by a network.
293 * We can optimize such communication by using one MPI call
294 * within each node and one between the nodes.
295 * For MVAPICH2 and Intel MPI this reduces the time for
296 * the global_stat communication by 25%
297 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
298 * B. Hess, November 2007
304 #ifndef GMX_THREAD_MPI
306 MPI_Comm_size(cr->mpi_comm_mygroup, &n);
307 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
309 nodehash = gmx_physicalnode_id_hash();
313 fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
317 /* The intra-node communicator, split on node number */
318 MPI_Comm_split(cr->mpi_comm_mygroup, nodehash, rank, &nc->comm_intra);
319 MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
322 fprintf(debug, "In gmx_setup_nodecomm: node rank %d rank_intra %d\n",
323 rank, nc->rank_intra);
325 /* The inter-node communicator, split on rank_intra.
326 * We actually only need the one for rank=0,
327 * but it is easier to create them all.
329 MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
330 /* Check if this really created two step communication */
331 MPI_Comm_size(nc->comm_inter, &ng);
332 MPI_Comm_size(nc->comm_intra, &ni);
335 fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
339 if (getenv("GMX_NO_NODECOMM") == NULL &&
340 ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
345 fprintf(fplog, "Using two step summing over %d groups of on average %.1f processes\n\n",
346 ng, (real)n/(real)ng);
348 if (nc->rank_intra > 0)
350 MPI_Comm_free(&nc->comm_inter);
355 /* One group or all processes in a separate group, use normal summing */
356 MPI_Comm_free(&nc->comm_inter);
357 MPI_Comm_free(&nc->comm_intra);
360 fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
365 /* tMPI runs only on a single node so just use the nodeid */
366 nc->rank_intra = cr->nodeid;
370 void gmx_init_intranode_counters(t_commrec *cr)
372 /* counters for PP+PME and PP-only processes on my physical node */
373 int nrank_intranode, rank_intranode;
374 int nrank_pp_intranode, rank_pp_intranode;
375 /* thread-MPI is not initialized when not running in parallel */
376 #if defined GMX_MPI && !defined GMX_THREAD_MPI
377 int nrank_world, rank_world;
378 int i, myhash, *hash, *hash_s, *hash_pp, *hash_pp_s;
380 MPI_Comm_size(MPI_COMM_WORLD, &nrank_world);
381 MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
383 /* Get a (hopefully unique) hash that identifies our physical node */
384 myhash = gmx_physicalnode_id_hash();
386 /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
387 snew(hash, nrank_world);
388 snew(hash_s, nrank_world);
389 snew(hash_pp, nrank_world);
390 snew(hash_pp_s, nrank_world);
392 hash_s[rank_world] = myhash;
393 hash_pp_s[rank_world] = (cr->duty & DUTY_PP) ? myhash : -1;
395 MPI_Allreduce(hash_s, hash, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
396 MPI_Allreduce(hash_pp_s, hash_pp, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
400 nrank_pp_intranode = 0;
401 rank_pp_intranode = 0;
402 for (i = 0; i < nrank_world; i++)
404 if (hash[i] == myhash)
412 if (hash_pp[i] == myhash)
414 nrank_pp_intranode++;
415 if ((cr->duty & DUTY_PP) && i < rank_world)
426 /* Serial or thread-MPI code: we run within a single physical node */
427 nrank_intranode = cr->nnodes;
428 rank_intranode = cr->sim_nodeid;
429 nrank_pp_intranode = cr->nnodes - cr->npmenodes;
430 rank_pp_intranode = cr->nodeid;
436 if (cr->duty & DUTY_PP && cr->duty & DUTY_PME)
438 sprintf(sbuf, "PP+PME");
442 sprintf(sbuf, "%s", cr->duty & DUTY_PP ? "PP" : "PME");
444 fprintf(debug, "On %3s node %d: nrank_intranode=%d, rank_intranode=%d, "
445 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
446 sbuf, cr->sim_nodeid,
447 nrank_intranode, rank_intranode,
448 nrank_pp_intranode, rank_pp_intranode);
451 cr->nrank_intranode = nrank_intranode;
452 cr->rank_intranode = rank_intranode;
453 cr->nrank_pp_intranode = nrank_pp_intranode;
454 cr->rank_pp_intranode = rank_pp_intranode;
458 void gmx_barrier(const t_commrec gmx_unused *cr)
461 gmx_call("gmx_barrier");
463 MPI_Barrier(cr->mpi_comm_mygroup);
467 void gmx_abort(int gmx_unused noderank, int gmx_unused nnodes, int gmx_unused errorno)
470 gmx_call("gmx_abort");
472 #ifdef GMX_THREAD_MPI
473 fprintf(stderr, "Halting program %s\n", ShortProgram());
479 fprintf(stderr, "Halting parallel program %s on CPU %d out of %d\n",
480 ShortProgram(), noderank, nnodes);
484 fprintf(stderr, "Halting program %s\n", ShortProgram());
488 MPI_Abort(MPI_COMM_WORLD, errorno);
494 void gmx_bcast(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
497 gmx_call("gmx_bast");
499 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
503 void gmx_bcast_sim(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
506 gmx_call("gmx_bast");
508 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
512 void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unused *cr)
515 gmx_call("gmx_sumd");
517 #if defined(MPI_IN_PLACE_EXISTS)
520 if (cr->nc.rank_intra == 0)
522 /* Use two step summing. */
523 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
525 /* Sum the roots of the internal (intra) buffers. */
526 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
531 /* This is here because of the silly MPI specification
532 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
533 MPI_Reduce(r, NULL, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
535 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
539 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
540 cr->mpi_comm_mygroup);
545 if (nr > cr->mpb->dbuf_alloc)
547 cr->mpb->dbuf_alloc = nr;
548 srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
552 /* Use two step summing */
553 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
554 if (cr->nc.rank_intra == 0)
556 /* Sum with the buffers reversed */
557 MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
560 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
564 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
565 cr->mpi_comm_mygroup);
566 for (i = 0; i < nr; i++)
568 r[i] = cr->mpb->dbuf[i];
575 void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unused *cr)
578 gmx_call("gmx_sumf");
580 #if defined(MPI_IN_PLACE_EXISTS)
583 /* Use two step summing. */
584 if (cr->nc.rank_intra == 0)
586 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
588 /* Sum the roots of the internal (intra) buffers */
589 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
594 /* This is here because of the silly MPI specification
595 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
596 MPI_Reduce(r, NULL, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
598 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
602 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
607 if (nr > cr->mpb->fbuf_alloc)
609 cr->mpb->fbuf_alloc = nr;
610 srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
614 /* Use two step summing */
615 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
616 if (cr->nc.rank_intra == 0)
618 /* Sum with the buffers reversed */
619 MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
622 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
626 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
627 cr->mpi_comm_mygroup);
628 for (i = 0; i < nr; i++)
630 r[i] = cr->mpb->fbuf[i];
637 void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused *cr)
640 gmx_call("gmx_sumi");
642 #if defined(MPI_IN_PLACE_EXISTS)
645 /* Use two step summing */
646 if (cr->nc.rank_intra == 0)
648 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
649 /* Sum with the buffers reversed */
650 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
654 /* This is here because of the silly MPI specification
655 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
656 MPI_Reduce(r, NULL, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
658 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
662 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
667 if (nr > cr->mpb->ibuf_alloc)
669 cr->mpb->ibuf_alloc = nr;
670 srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
674 /* Use two step summing */
675 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
676 if (cr->nc.rank_intra == 0)
678 /* Sum with the buffers reversed */
679 MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
681 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
685 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
686 for (i = 0; i < nr; i++)
688 r[i] = cr->mpb->ibuf[i];
695 void gmx_sumli(int gmx_unused nr, gmx_int64_t gmx_unused r[], const t_commrec gmx_unused *cr)
698 gmx_call("gmx_sumli");
700 #if defined(MPI_IN_PLACE_EXISTS)
703 /* Use two step summing */
704 if (cr->nc.rank_intra == 0)
706 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, 0,
708 /* Sum with the buffers reversed */
709 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
714 /* This is here because of the silly MPI specification
715 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
716 MPI_Reduce(r, NULL, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
718 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
722 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
727 if (nr > cr->mpb->libuf_alloc)
729 cr->mpb->libuf_alloc = nr;
730 srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
734 /* Use two step summing */
735 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
737 if (cr->nc.rank_intra == 0)
739 /* Sum with the buffers reversed */
740 MPI_Allreduce(cr->mpb->libuf, r, nr, MPI_INT64_T, MPI_SUM,
743 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
747 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
748 cr->mpi_comm_mygroup);
749 for (i = 0; i < nr; i++)
751 r[i] = cr->mpb->libuf[i];
761 static void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
763 #if defined(MPI_IN_PLACE_EXISTS)
764 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
766 /* this function is only used in code that is not performance critical,
767 (during setup, when comm_rec is not the appropriate communication
768 structure), so this isn't as bad as it looks. */
773 MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
774 for (i = 0; i < nr; i++)
784 static void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
786 #if defined(MPI_IN_PLACE_EXISTS)
787 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
789 /* this function is only used in code that is not performance critical,
790 (during setup, when comm_rec is not the appropriate communication
791 structure), so this isn't as bad as it looks. */
796 MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
797 for (i = 0; i < nr; i++)
806 void gmx_sumd_sim(int gmx_unused nr, double gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
809 gmx_call("gmx_sumd_sim");
811 gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
815 void gmx_sumf_sim(int gmx_unused nr, float gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
818 gmx_call("gmx_sumf_sim");
820 gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
824 void gmx_sumi_sim(int gmx_unused nr, int gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
827 gmx_call("gmx_sumi_sim");
829 #if defined(MPI_IN_PLACE_EXISTS)
830 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
832 /* this is thread-unsafe, but it will do for now: */
835 if (nr > ms->mpb->ibuf_alloc)
837 ms->mpb->ibuf_alloc = nr;
838 srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
840 MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
841 for (i = 0; i < nr; i++)
843 r[i] = ms->mpb->ibuf[i];
849 void gmx_sumli_sim(int gmx_unused nr, gmx_int64_t gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
852 gmx_call("gmx_sumli_sim");
854 #if defined(MPI_IN_PLACE_EXISTS)
855 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
856 ms->mpi_comm_masters);
858 /* this is thread-unsafe, but it will do for now: */
861 if (nr > ms->mpb->libuf_alloc)
863 ms->mpb->libuf_alloc = nr;
864 srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
866 MPI_Allreduce(r, ms->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
867 ms->mpi_comm_masters);
868 for (i = 0; i < nr; i++)
870 r[i] = ms->mpb->libuf[i];