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 void gmx_fill_commrec_from_mpi(t_commrec gmx_unused *cr)
70 gmx_call("gmx_fill_commrec_from_mpi");
72 if (!gmx_mpi_initialized())
74 gmx_comm("MPI has not been initialized properly");
77 cr->nnodes = gmx_node_num();
78 cr->nodeid = gmx_node_rank();
79 cr->sim_nodeid = cr->nodeid;
80 cr->mpi_comm_mysim = MPI_COMM_WORLD;
81 cr->mpi_comm_mygroup = MPI_COMM_WORLD;
86 t_commrec *init_commrec()
93 gmx_fill_commrec_from_mpi(cr);
95 cr->mpi_comm_mysim = NULL;
96 cr->mpi_comm_mygroup = NULL;
99 cr->nodeid = cr->sim_nodeid;
102 // TODO cr->duty should not be initialized here
103 cr->duty = (DUTY_PP | DUTY_PME);
105 #if defined GMX_MPI && !defined MPI_IN_PLACE_EXISTS
106 /* initialize the MPI_IN_PLACE replacement buffers */
108 cr->mpb->ibuf = NULL;
109 cr->mpb->libuf = NULL;
110 cr->mpb->fbuf = NULL;
111 cr->mpb->dbuf = NULL;
112 cr->mpb->ibuf_alloc = 0;
113 cr->mpb->libuf_alloc = 0;
114 cr->mpb->fbuf_alloc = 0;
115 cr->mpb->dbuf_alloc = 0;
121 t_commrec *reinitialize_commrec_for_this_thread(const t_commrec gmx_unused *cro)
123 #ifdef GMX_THREAD_MPI
126 /* make a thread-specific commrec */
128 /* now copy the whole thing, so settings like the number of PME nodes
132 /* and we start setting our own thread-specific values for things */
133 gmx_fill_commrec_from_mpi(cr);
135 // TODO cr->duty should not be initialized here
136 cr->duty = (DUTY_PP | DUTY_PME);
144 int gmx_node_num(void)
150 (void) MPI_Comm_size(MPI_COMM_WORLD, &i);
155 int gmx_node_rank(void)
161 (void) MPI_Comm_rank(MPI_COMM_WORLD, &i);
166 #if defined GMX_LIB_MPI && defined GMX_TARGET_BGQ
167 #include <spi/include/kernel/location.h>
170 int gmx_physicalnode_id_hash(void)
175 /* We have a single physical node */
179 char mpi_hostname[MPI_MAX_PROCESSOR_NAME];
181 /* This procedure can only differentiate nodes with different names.
182 * Architectures where different physical nodes have identical names,
183 * such as IBM Blue Gene, should use an architecture specific solution.
185 MPI_Get_processor_name(mpi_hostname, &resultlen);
187 /* The string hash function returns an unsigned int. We cast to an int.
188 * Negative numbers are converted to positive by setting the sign bit to 0.
189 * This makes the hash one bit smaller.
190 * A 63-bit hash (with 64-bit int) should be enough for unique node hashes,
191 * even on a million node machine. 31 bits might not be enough though!
194 (int)gmx_string_fullhash_func(mpi_hostname, gmx_string_hash_init);
204 /* TODO: this function should be fully replaced by gmx_physicalnode_id_hash */
205 int gmx_hostname_num()
210 #ifdef GMX_THREAD_MPI
211 /* thread-MPI currently puts the thread number in the process name,
212 * we might want to change this, as this is inconsistent with what
213 * most MPI implementations would do when running on a single node.
217 int resultlen, hostnum, i, j;
218 char mpi_hostname[MPI_MAX_PROCESSOR_NAME], hostnum_str[MPI_MAX_PROCESSOR_NAME];
220 MPI_Get_processor_name(mpi_hostname, &resultlen);
221 #ifdef GMX_TARGET_BGQ
222 Personality_t personality;
223 Kernel_GetPersonality(&personality, sizeof(personality));
224 /* Each MPI rank has a unique coordinate in a 6-dimensional space
225 (A,B,C,D,E,T), with dimensions A-E corresponding to different
226 physical nodes, and T within each node. Each node has sixteen
227 physical cores, each of which can have up to four hardware
228 threads, so 0 <= T <= 63 (but the maximum value of T depends on
229 the confituration of ranks and OpenMP threads per
230 node). However, T is irrelevant for computing a suitable return
231 value for gmx_hostname_num().
233 hostnum = personality.Network_Config.Acoord;
234 hostnum *= personality.Network_Config.Bnodes;
235 hostnum += personality.Network_Config.Bcoord;
236 hostnum *= personality.Network_Config.Cnodes;
237 hostnum += personality.Network_Config.Ccoord;
238 hostnum *= personality.Network_Config.Dnodes;
239 hostnum += personality.Network_Config.Dcoord;
240 hostnum *= personality.Network_Config.Enodes;
241 hostnum += personality.Network_Config.Ecoord;
243 /* This procedure can only differentiate nodes with host names
244 * that end on unique numbers.
248 /* Only parse the host name up to the first dot */
249 while (i < resultlen && mpi_hostname[i] != '.')
251 if (isdigit(mpi_hostname[i]))
253 hostnum_str[j++] = mpi_hostname[i];
257 hostnum_str[j] = '\0';
264 /* Use only the last 9 decimals, so we don't overflow an int */
265 hostnum = strtol(hostnum_str + max(0, j-9), NULL, 10);
271 fprintf(debug, "In gmx_hostname_num: hostname '%s', hostnum %d\n",
272 mpi_hostname, hostnum);
273 #ifdef GMX_TARGET_BGQ
275 "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",
276 personality.Network_Config.Acoord,
277 personality.Network_Config.Anodes,
278 personality.Network_Config.Bcoord,
279 personality.Network_Config.Bnodes,
280 personality.Network_Config.Ccoord,
281 personality.Network_Config.Cnodes,
282 personality.Network_Config.Dcoord,
283 personality.Network_Config.Dnodes,
284 personality.Network_Config.Ecoord,
285 personality.Network_Config.Enodes,
286 Kernel_ProcessorCoreID(),
288 Kernel_ProcessorID(),
290 Kernel_ProcessorThreadID(),
299 void gmx_setup_nodecomm(FILE gmx_unused *fplog, t_commrec *cr)
302 int n, rank, hostnum, ng, ni;
304 /* Many MPI implementations do not optimize MPI_Allreduce
305 * (and probably also other global communication calls)
306 * for multi-core nodes connected by a network.
307 * We can optimize such communication by using one MPI call
308 * within each node and one between the nodes.
309 * For MVAPICH2 and Intel MPI this reduces the time for
310 * the global_stat communication by 25%
311 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
312 * B. Hess, November 2007
318 #ifndef GMX_THREAD_MPI
320 MPI_Comm_size(cr->mpi_comm_mygroup, &n);
321 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
323 hostnum = gmx_hostname_num();
327 fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
331 /* The intra-node communicator, split on node number */
332 MPI_Comm_split(cr->mpi_comm_mygroup, hostnum, rank, &nc->comm_intra);
333 MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
336 fprintf(debug, "In gmx_setup_nodecomm: node rank %d rank_intra %d\n",
337 rank, nc->rank_intra);
339 /* The inter-node communicator, split on rank_intra.
340 * We actually only need the one for rank=0,
341 * but it is easier to create them all.
343 MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
344 /* Check if this really created two step communication */
345 MPI_Comm_size(nc->comm_inter, &ng);
346 MPI_Comm_size(nc->comm_intra, &ni);
349 fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
353 if (getenv("GMX_NO_NODECOMM") == NULL &&
354 ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
359 fprintf(fplog, "Using two step summing over %d groups of on average %.1f processes\n\n",
360 ng, (real)n/(real)ng);
362 if (nc->rank_intra > 0)
364 MPI_Comm_free(&nc->comm_inter);
369 /* One group or all processes in a separate group, use normal summing */
370 MPI_Comm_free(&nc->comm_inter);
371 MPI_Comm_free(&nc->comm_intra);
374 fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
379 /* tMPI runs only on a single node so just use the nodeid */
380 nc->rank_intra = cr->nodeid;
384 void gmx_init_intranode_counters(t_commrec *cr)
386 /* counters for PP+PME and PP-only processes on my physical node */
387 int nrank_intranode, rank_intranode;
388 int nrank_pp_intranode, rank_pp_intranode;
389 /* thread-MPI is not initialized when not running in parallel */
390 #if defined GMX_MPI && !defined GMX_THREAD_MPI
391 int nrank_world, rank_world;
392 int i, mynum, *num, *num_s, *num_pp, *num_pp_s;
394 MPI_Comm_size(MPI_COMM_WORLD, &nrank_world);
395 MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
397 /* Get the node number from the hostname to identify the nodes */
398 mynum = gmx_hostname_num();
400 /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
401 snew(num, nrank_world);
402 snew(num_s, nrank_world);
403 snew(num_pp, nrank_world);
404 snew(num_pp_s, nrank_world);
406 num_s[rank_world] = mynum;
407 num_pp_s[rank_world] = (cr->duty & DUTY_PP) ? mynum : -1;
409 MPI_Allreduce(num_s, num, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
410 MPI_Allreduce(num_pp_s, num_pp, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
414 nrank_pp_intranode = 0;
415 rank_pp_intranode = 0;
416 for (i = 0; i < nrank_world; i++)
426 if ((cr->duty & DUTY_PP) && num_pp[i] == mynum)
428 nrank_pp_intranode++;
440 /* Serial or thread-MPI code: we run within a single physical node */
441 nrank_intranode = cr->nnodes;
442 rank_intranode = cr->sim_nodeid;
443 nrank_pp_intranode = cr->nnodes - cr->npmenodes;
444 rank_pp_intranode = cr->nodeid;
450 if (cr->duty & DUTY_PP && cr->duty & DUTY_PME)
452 sprintf(sbuf, "PP+PME");
456 sprintf(sbuf, "%s", cr->duty & DUTY_PP ? "PP" : "PME");
458 fprintf(debug, "On %3s node %d: nrank_intranode=%d, rank_intranode=%d, "
459 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
460 sbuf, cr->sim_nodeid,
461 nrank_intranode, rank_intranode,
462 nrank_pp_intranode, rank_pp_intranode);
465 cr->nrank_intranode = nrank_intranode;
466 cr->rank_intranode = rank_intranode;
467 cr->nrank_pp_intranode = nrank_pp_intranode;
468 cr->rank_pp_intranode = rank_pp_intranode;
472 void gmx_barrier(const t_commrec gmx_unused *cr)
475 gmx_call("gmx_barrier");
477 MPI_Barrier(cr->mpi_comm_mygroup);
481 void gmx_abort(int gmx_unused noderank, int gmx_unused nnodes, int gmx_unused errorno)
484 gmx_call("gmx_abort");
486 #ifdef GMX_THREAD_MPI
487 fprintf(stderr, "Halting program %s\n", ShortProgram());
493 fprintf(stderr, "Halting parallel program %s on CPU %d out of %d\n",
494 ShortProgram(), noderank, nnodes);
498 fprintf(stderr, "Halting program %s\n", ShortProgram());
502 MPI_Abort(MPI_COMM_WORLD, errorno);
508 void gmx_bcast(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
511 gmx_call("gmx_bast");
513 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
517 void gmx_bcast_sim(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
520 gmx_call("gmx_bast");
522 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
526 void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unused *cr)
529 gmx_call("gmx_sumd");
531 #if defined(MPI_IN_PLACE_EXISTS)
534 if (cr->nc.rank_intra == 0)
536 /* Use two step summing. */
537 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
539 /* Sum the roots of the internal (intra) buffers. */
540 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
545 /* This is here because of the silly MPI specification
546 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
547 MPI_Reduce(r, NULL, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
549 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
553 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
554 cr->mpi_comm_mygroup);
559 if (nr > cr->mpb->dbuf_alloc)
561 cr->mpb->dbuf_alloc = nr;
562 srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
566 /* Use two step summing */
567 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
568 if (cr->nc.rank_intra == 0)
570 /* Sum with the buffers reversed */
571 MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
574 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
578 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
579 cr->mpi_comm_mygroup);
580 for (i = 0; i < nr; i++)
582 r[i] = cr->mpb->dbuf[i];
589 void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unused *cr)
592 gmx_call("gmx_sumf");
594 #if defined(MPI_IN_PLACE_EXISTS)
597 /* Use two step summing. */
598 if (cr->nc.rank_intra == 0)
600 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
602 /* Sum the roots of the internal (intra) buffers */
603 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
608 /* This is here because of the silly MPI specification
609 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
610 MPI_Reduce(r, NULL, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
612 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
616 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
621 if (nr > cr->mpb->fbuf_alloc)
623 cr->mpb->fbuf_alloc = nr;
624 srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
628 /* Use two step summing */
629 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
630 if (cr->nc.rank_intra == 0)
632 /* Sum with the buffers reversed */
633 MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
636 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
640 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
641 cr->mpi_comm_mygroup);
642 for (i = 0; i < nr; i++)
644 r[i] = cr->mpb->fbuf[i];
651 void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused *cr)
654 gmx_call("gmx_sumi");
656 #if defined(MPI_IN_PLACE_EXISTS)
659 /* Use two step summing */
660 if (cr->nc.rank_intra == 0)
662 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
663 /* Sum with the buffers reversed */
664 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
668 /* This is here because of the silly MPI specification
669 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
670 MPI_Reduce(r, NULL, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
672 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
676 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
681 if (nr > cr->mpb->ibuf_alloc)
683 cr->mpb->ibuf_alloc = nr;
684 srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
688 /* Use two step summing */
689 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
690 if (cr->nc.rank_intra == 0)
692 /* Sum with the buffers reversed */
693 MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
695 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
699 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
700 for (i = 0; i < nr; i++)
702 r[i] = cr->mpb->ibuf[i];
709 void gmx_sumli(int gmx_unused nr, gmx_int64_t gmx_unused r[], const t_commrec gmx_unused *cr)
712 gmx_call("gmx_sumli");
714 #if defined(MPI_IN_PLACE_EXISTS)
717 /* Use two step summing */
718 if (cr->nc.rank_intra == 0)
720 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, 0,
722 /* Sum with the buffers reversed */
723 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
728 /* This is here because of the silly MPI specification
729 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
730 MPI_Reduce(r, NULL, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
732 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
736 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
741 if (nr > cr->mpb->libuf_alloc)
743 cr->mpb->libuf_alloc = nr;
744 srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
748 /* Use two step summing */
749 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
751 if (cr->nc.rank_intra == 0)
753 /* Sum with the buffers reversed */
754 MPI_Allreduce(cr->mpb->libuf, r, nr, MPI_INT64_T, MPI_SUM,
757 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
761 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
762 cr->mpi_comm_mygroup);
763 for (i = 0; i < nr; i++)
765 r[i] = cr->mpb->libuf[i];
775 void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
777 #if defined(MPI_IN_PLACE_EXISTS)
778 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
780 /* this function is only used in code that is not performance critical,
781 (during setup, when comm_rec is not the appropriate communication
782 structure), so this isn't as bad as it looks. */
787 MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
788 for (i = 0; i < nr; i++)
798 void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
800 #if defined(MPI_IN_PLACE_EXISTS)
801 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
803 /* this function is only used in code that is not performance critical,
804 (during setup, when comm_rec is not the appropriate communication
805 structure), so this isn't as bad as it looks. */
810 MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
811 for (i = 0; i < nr; i++)
820 void gmx_sumd_sim(int gmx_unused nr, double gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
823 gmx_call("gmx_sumd_sim");
825 gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
829 void gmx_sumf_sim(int gmx_unused nr, float gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
832 gmx_call("gmx_sumf_sim");
834 gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
838 void gmx_sumi_sim(int gmx_unused nr, int gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
841 gmx_call("gmx_sumi_sim");
843 #if defined(MPI_IN_PLACE_EXISTS)
844 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
846 /* this is thread-unsafe, but it will do for now: */
849 if (nr > ms->mpb->ibuf_alloc)
851 ms->mpb->ibuf_alloc = nr;
852 srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
854 MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
855 for (i = 0; i < nr; i++)
857 r[i] = ms->mpb->ibuf[i];
863 void gmx_sumli_sim(int gmx_unused nr, gmx_int64_t gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
866 gmx_call("gmx_sumli_sim");
868 #if defined(MPI_IN_PLACE_EXISTS)
869 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
870 ms->mpi_comm_masters);
872 /* this is thread-unsafe, but it will do for now: */
875 if (nr > ms->mpb->libuf_alloc)
877 ms->mpb->libuf_alloc = nr;
878 srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
880 MPI_Allreduce(r, ms->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
881 ms->mpi_comm_masters);
882 for (i = 0; i < nr; i++)
884 r[i] = ms->mpb->libuf[i];