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"
50 #include "gromacs/utility/gmxmpi.h"
53 /* The source code in this file should be thread-safe.
54 Please keep it that way. */
56 gmx_bool gmx_mpi_initialized(void)
68 void gmx_fill_commrec_from_mpi(t_commrec *cr)
71 gmx_call("gmx_fill_commrec_from_mpi");
73 if (!gmx_mpi_initialized())
75 gmx_comm("MPI has not been initialized properly");
78 cr->nnodes = gmx_node_num();
79 cr->nodeid = gmx_node_rank();
80 cr->sim_nodeid = cr->nodeid;
81 cr->mpi_comm_mysim = MPI_COMM_WORLD;
82 cr->mpi_comm_mygroup = MPI_COMM_WORLD;
87 t_commrec *init_commrec()
94 gmx_fill_commrec_from_mpi(cr);
96 cr->mpi_comm_mysim = NULL;
97 cr->mpi_comm_mygroup = NULL;
100 cr->nodeid = cr->sim_nodeid;
103 // TODO cr->duty should not be initialized here
104 cr->duty = (DUTY_PP | DUTY_PME);
106 #if defined GMX_LIB_MPI && !defined MPI_IN_PLACE_EXISTS
107 /* initialize the MPI_IN_PLACE replacement buffers */
109 cr->mpb->ibuf = NULL;
110 cr->mpb->libuf = NULL;
111 cr->mpb->fbuf = NULL;
112 cr->mpb->dbuf = NULL;
113 cr->mpb->ibuf_alloc = 0;
114 cr->mpb->libuf_alloc = 0;
115 cr->mpb->fbuf_alloc = 0;
116 cr->mpb->dbuf_alloc = 0;
122 t_commrec *reinitialize_commrec_for_this_thread(const t_commrec *cro)
124 #ifdef GMX_THREAD_MPI
127 /* make a thread-specific commrec */
129 /* now copy the whole thing, so settings like the number of PME nodes
133 /* and we start setting our own thread-specific values for things */
134 gmx_fill_commrec_from_mpi(cr);
136 // TODO cr->duty should not be initialized here
137 cr->duty = (DUTY_PP | DUTY_PME);
145 int gmx_node_num(void)
151 (void) MPI_Comm_size(MPI_COMM_WORLD, &i);
156 int gmx_node_rank(void)
162 (void) MPI_Comm_rank(MPI_COMM_WORLD, &i);
167 #if defined GMX_LIB_MPI && defined GMX_TARGET_BGQ
168 #include <spi/include/kernel/location.h>
171 int gmx_physicalnode_id_hash(void)
176 /* We have a single physical node */
180 char mpi_hostname[MPI_MAX_PROCESSOR_NAME];
182 /* This procedure can only differentiate nodes with different names.
183 * Architectures where different physical nodes have identical names,
184 * such as IBM Blue Gene, should use an architecture specific solution.
186 MPI_Get_processor_name(mpi_hostname, &resultlen);
188 /* The string hash function returns an unsigned int. We cast to an int.
189 * Negative numbers are converted to positive by setting the sign bit to 0.
190 * This makes the hash one bit smaller.
191 * A 63-bit hash (with 64-bit int) should be enough for unique node hashes,
192 * even on a million node machine. 31 bits might not be enough though!
195 (int)gmx_string_fullhash_func(mpi_hostname, gmx_string_hash_init);
205 /* TODO: this function should be fully replaced by gmx_physicalnode_id_hash */
206 int gmx_hostname_num()
211 #ifdef GMX_THREAD_MPI
212 /* thread-MPI currently puts the thread number in the process name,
213 * we might want to change this, as this is inconsistent with what
214 * most MPI implementations would do when running on a single node.
218 int resultlen, hostnum, i, j;
219 char mpi_hostname[MPI_MAX_PROCESSOR_NAME], hostnum_str[MPI_MAX_PROCESSOR_NAME];
221 MPI_Get_processor_name(mpi_hostname, &resultlen);
222 #ifdef GMX_TARGET_BGQ
223 Personality_t personality;
224 Kernel_GetPersonality(&personality, sizeof(personality));
225 /* Each MPI rank has a unique coordinate in a 6-dimensional space
226 (A,B,C,D,E,T), with dimensions A-E corresponding to different
227 physical nodes, and T within each node. Each node has sixteen
228 physical cores, each of which can have up to four hardware
229 threads, so 0 <= T <= 63 (but the maximum value of T depends on
230 the confituration of ranks and OpenMP threads per
231 node). However, T is irrelevant for computing a suitable return
232 value for gmx_hostname_num().
234 hostnum = personality.Network_Config.Acoord;
235 hostnum *= personality.Network_Config.Bnodes;
236 hostnum += personality.Network_Config.Bcoord;
237 hostnum *= personality.Network_Config.Cnodes;
238 hostnum += personality.Network_Config.Ccoord;
239 hostnum *= personality.Network_Config.Dnodes;
240 hostnum += personality.Network_Config.Dcoord;
241 hostnum *= personality.Network_Config.Enodes;
242 hostnum += personality.Network_Config.Ecoord;
244 /* This procedure can only differentiate nodes with host names
245 * that end on unique numbers.
249 /* Only parse the host name up to the first dot */
250 while (i < resultlen && mpi_hostname[i] != '.')
252 if (isdigit(mpi_hostname[i]))
254 hostnum_str[j++] = mpi_hostname[i];
258 hostnum_str[j] = '\0';
265 /* Use only the last 9 decimals, so we don't overflow an int */
266 hostnum = strtol(hostnum_str + max(0, j-9), NULL, 10);
272 fprintf(debug, "In gmx_hostname_num: hostname '%s', hostnum %d\n",
273 mpi_hostname, hostnum);
274 #ifdef GMX_TARGET_BGQ
276 "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",
277 personality.Network_Config.Acoord,
278 personality.Network_Config.Anodes,
279 personality.Network_Config.Bcoord,
280 personality.Network_Config.Bnodes,
281 personality.Network_Config.Ccoord,
282 personality.Network_Config.Cnodes,
283 personality.Network_Config.Dcoord,
284 personality.Network_Config.Dnodes,
285 personality.Network_Config.Ecoord,
286 personality.Network_Config.Enodes,
287 Kernel_ProcessorCoreID(),
289 Kernel_ProcessorID(),
291 Kernel_ProcessorThreadID(),
300 void gmx_setup_nodecomm(FILE gmx_unused *fplog, t_commrec *cr)
303 int n, rank, hostnum, ng, ni;
305 /* Many MPI implementations do not optimize MPI_Allreduce
306 * (and probably also other global communication calls)
307 * for multi-core nodes connected by a network.
308 * We can optimize such communication by using one MPI call
309 * within each node and one between the nodes.
310 * For MVAPICH2 and Intel MPI this reduces the time for
311 * the global_stat communication by 25%
312 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
313 * B. Hess, November 2007
319 #ifndef GMX_THREAD_MPI
321 MPI_Comm_size(cr->mpi_comm_mygroup, &n);
322 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
324 hostnum = gmx_hostname_num();
328 fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
332 /* The intra-node communicator, split on node number */
333 MPI_Comm_split(cr->mpi_comm_mygroup, hostnum, rank, &nc->comm_intra);
334 MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
337 fprintf(debug, "In gmx_setup_nodecomm: node rank %d rank_intra %d\n",
338 rank, nc->rank_intra);
340 /* The inter-node communicator, split on rank_intra.
341 * We actually only need the one for rank=0,
342 * but it is easier to create them all.
344 MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
345 /* Check if this really created two step communication */
346 MPI_Comm_size(nc->comm_inter, &ng);
347 MPI_Comm_size(nc->comm_intra, &ni);
350 fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
354 if (getenv("GMX_NO_NODECOMM") == NULL &&
355 ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
360 fprintf(fplog, "Using two step summing over %d groups of on average %.1f processes\n\n",
361 ng, (real)n/(real)ng);
363 if (nc->rank_intra > 0)
365 MPI_Comm_free(&nc->comm_inter);
370 /* One group or all processes in a separate group, use normal summing */
371 MPI_Comm_free(&nc->comm_inter);
372 MPI_Comm_free(&nc->comm_intra);
375 fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
380 /* tMPI runs only on a single node so just use the nodeid */
381 nc->rank_intra = cr->nodeid;
385 void gmx_init_intranode_counters(t_commrec *cr)
387 /* counters for PP+PME and PP-only processes on my physical node */
388 int nrank_intranode, rank_intranode;
389 int nrank_pp_intranode, rank_pp_intranode;
390 /* thread-MPI is not initialized when not running in parallel */
391 #if defined GMX_MPI && !defined GMX_THREAD_MPI
392 int nrank_world, rank_world;
393 int i, mynum, *num, *num_s, *num_pp, *num_pp_s;
395 MPI_Comm_size(MPI_COMM_WORLD, &nrank_world);
396 MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
398 /* Get the node number from the hostname to identify the nodes */
399 mynum = gmx_hostname_num();
401 /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
402 snew(num, nrank_world);
403 snew(num_s, nrank_world);
404 snew(num_pp, nrank_world);
405 snew(num_pp_s, nrank_world);
407 num_s[rank_world] = mynum;
408 num_pp_s[rank_world] = (cr->duty & DUTY_PP) ? mynum : -1;
410 MPI_Allreduce(num_s, num, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
411 MPI_Allreduce(num_pp_s, num_pp, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
415 nrank_pp_intranode = 0;
416 rank_pp_intranode = 0;
417 for (i = 0; i < nrank_world; i++)
427 if ((cr->duty & DUTY_PP) && num_pp[i] == mynum)
429 nrank_pp_intranode++;
441 /* Serial or thread-MPI code: we run within a single physical node */
442 nrank_intranode = cr->nnodes;
443 rank_intranode = cr->sim_nodeid;
444 nrank_pp_intranode = cr->nnodes - cr->npmenodes;
445 rank_pp_intranode = cr->nodeid;
451 if (cr->duty & DUTY_PP && cr->duty & DUTY_PME)
453 sprintf(sbuf, "PP+PME");
457 sprintf(sbuf, "%s", cr->duty & DUTY_PP ? "PP" : "PME");
459 fprintf(debug, "On %3s node %d: nrank_intranode=%d, rank_intranode=%d, "
460 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
461 sbuf, cr->sim_nodeid,
462 nrank_intranode, rank_intranode,
463 nrank_pp_intranode, rank_pp_intranode);
466 cr->nrank_intranode = nrank_intranode;
467 cr->rank_intranode = rank_intranode;
468 cr->nrank_pp_intranode = nrank_pp_intranode;
469 cr->rank_pp_intranode = rank_pp_intranode;
473 void gmx_barrier(const t_commrec *cr)
476 gmx_call("gmx_barrier");
478 MPI_Barrier(cr->mpi_comm_mygroup);
482 void gmx_abort(int gmx_unused noderank, int gmx_unused nnodes, int gmx_unused errorno)
485 gmx_call("gmx_abort");
487 #ifdef GMX_THREAD_MPI
488 fprintf(stderr, "Halting program %s\n", ShortProgram());
494 fprintf(stderr, "Halting parallel program %s on CPU %d out of %d\n",
495 ShortProgram(), noderank, nnodes);
499 fprintf(stderr, "Halting program %s\n", ShortProgram());
503 MPI_Abort(MPI_COMM_WORLD, errorno);
509 void gmx_bcast(int nbytes, void *b, const t_commrec *cr)
512 gmx_call("gmx_bast");
514 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
518 void gmx_bcast_sim(int nbytes, void *b, const t_commrec *cr)
521 gmx_call("gmx_bast");
523 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
527 void gmx_sumd(int nr, double r[], const t_commrec *cr)
530 gmx_call("gmx_sumd");
532 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
535 if (cr->nc.rank_intra == 0)
537 /* Use two step summing. */
538 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
540 /* Sum the roots of the internal (intra) buffers. */
541 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
546 /* This is here because of the silly MPI specification
547 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
548 MPI_Reduce(r, NULL, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
550 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
554 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
555 cr->mpi_comm_mygroup);
560 if (nr > cr->mpb->dbuf_alloc)
562 cr->mpb->dbuf_alloc = nr;
563 srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
567 /* Use two step summing */
568 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
569 if (cr->nc.rank_intra == 0)
571 /* Sum with the buffers reversed */
572 MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
575 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
579 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
580 cr->mpi_comm_mygroup);
581 for (i = 0; i < nr; i++)
583 r[i] = cr->mpb->dbuf[i];
590 void gmx_sumf(int nr, float r[], const t_commrec *cr)
593 gmx_call("gmx_sumf");
595 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
598 /* Use two step summing. */
599 if (cr->nc.rank_intra == 0)
601 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
603 /* Sum the roots of the internal (intra) buffers */
604 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
609 /* This is here because of the silly MPI specification
610 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
611 MPI_Reduce(r, NULL, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
613 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
617 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
622 if (nr > cr->mpb->fbuf_alloc)
624 cr->mpb->fbuf_alloc = nr;
625 srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
629 /* Use two step summing */
630 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
631 if (cr->nc.rank_intra == 0)
633 /* Sum with the buffers reversed */
634 MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
637 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
641 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
642 cr->mpi_comm_mygroup);
643 for (i = 0; i < nr; i++)
645 r[i] = cr->mpb->fbuf[i];
652 void gmx_sumi(int nr, int r[], const t_commrec *cr)
655 gmx_call("gmx_sumi");
657 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
660 /* Use two step summing */
661 if (cr->nc.rank_intra == 0)
663 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
664 /* Sum with the buffers reversed */
665 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
669 /* This is here because of the silly MPI specification
670 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
671 MPI_Reduce(r, NULL, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
673 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
677 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
682 if (nr > cr->mpb->ibuf_alloc)
684 cr->mpb->ibuf_alloc = nr;
685 srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
689 /* Use two step summing */
690 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
691 if (cr->nc.rank_intra == 0)
693 /* Sum with the buffers reversed */
694 MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
696 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
700 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
701 for (i = 0; i < nr; i++)
703 r[i] = cr->mpb->ibuf[i];
710 void gmx_sumli(int nr, gmx_large_int_t r[], const t_commrec *cr)
713 gmx_call("gmx_sumli");
715 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
718 /* Use two step summing */
719 if (cr->nc.rank_intra == 0)
721 MPI_Reduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM, 0,
723 /* Sum with the buffers reversed */
724 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
729 /* This is here because of the silly MPI specification
730 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
731 MPI_Reduce(r, NULL, nr, GMX_MPI_LARGE_INT, MPI_SUM, 0, cr->nc.comm_intra);
733 MPI_Bcast(r, nr, GMX_MPI_LARGE_INT, 0, cr->nc.comm_intra);
737 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM, cr->mpi_comm_mygroup);
742 if (nr > cr->mpb->libuf_alloc)
744 cr->mpb->libuf_alloc = nr;
745 srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
749 /* Use two step summing */
750 MPI_Allreduce(r, cr->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
752 if (cr->nc.rank_intra == 0)
754 /* Sum with the buffers reversed */
755 MPI_Allreduce(cr->mpb->libuf, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
758 MPI_Bcast(r, nr, GMX_MPI_LARGE_INT, 0, cr->nc.comm_intra);
762 MPI_Allreduce(r, cr->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
763 cr->mpi_comm_mygroup);
764 for (i = 0; i < nr; i++)
766 r[i] = cr->mpb->libuf[i];
776 void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
778 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
779 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
781 /* this function is only used in code that is not performance critical,
782 (during setup, when comm_rec is not the appropriate communication
783 structure), so this isn't as bad as it looks. */
788 MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
789 for (i = 0; i < nr; i++)
799 void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
801 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
802 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
804 /* this function is only used in code that is not performance critical,
805 (during setup, when comm_rec is not the appropriate communication
806 structure), so this isn't as bad as it looks. */
811 MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
812 for (i = 0; i < nr; i++)
821 void gmx_sumd_sim(int nr, double r[], const gmx_multisim_t *ms)
824 gmx_call("gmx_sumd_sim");
826 gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
830 void gmx_sumf_sim(int nr, float r[], const gmx_multisim_t *ms)
833 gmx_call("gmx_sumf_sim");
835 gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
839 void gmx_sumi_sim(int nr, int r[], const gmx_multisim_t *ms)
842 gmx_call("gmx_sumi_sim");
844 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
845 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
847 /* this is thread-unsafe, but it will do for now: */
850 if (nr > ms->mpb->ibuf_alloc)
852 ms->mpb->ibuf_alloc = nr;
853 srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
855 MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
856 for (i = 0; i < nr; i++)
858 r[i] = ms->mpb->ibuf[i];
864 void gmx_sumli_sim(int nr, gmx_large_int_t r[], const gmx_multisim_t *ms)
867 gmx_call("gmx_sumli_sim");
869 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
870 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
871 ms->mpi_comm_masters);
873 /* this is thread-unsafe, but it will do for now: */
876 if (nr > ms->mpb->libuf_alloc)
878 ms->mpb->libuf_alloc = nr;
879 srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
881 MPI_Allreduce(r, ms->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
882 ms->mpi_comm_masters);
883 for (i = 0; i < nr; i++)
885 r[i] = ms->mpb->libuf[i];