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 *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_LIB_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 *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_IS_BGQ
167 #include <spi/include/kernel/location.h>
170 int gmx_hostname_num()
175 #ifdef GMX_THREAD_MPI
176 /* thread-MPI currently puts the thread number in the process name,
177 * we might want to change this, as this is inconsistent with what
178 * most MPI implementations would do when running on a single node.
182 int resultlen, hostnum, i, j;
183 char mpi_hostname[MPI_MAX_PROCESSOR_NAME], hostnum_str[MPI_MAX_PROCESSOR_NAME];
185 MPI_Get_processor_name(mpi_hostname, &resultlen);
187 Personality_t personality;
188 Kernel_GetPersonality(&personality, sizeof(personality));
189 /* Each MPI rank has a unique coordinate in a 6-dimensional space
190 (A,B,C,D,E,T), with dimensions A-E corresponding to different
191 physical nodes, and T within each node. Each node has sixteen
192 physical cores, each of which can have up to four hardware
193 threads, so 0 <= T <= 63 (but the maximum value of T depends on
194 the confituration of ranks and OpenMP threads per
195 node). However, T is irrelevant for computing a suitable return
196 value for gmx_hostname_num().
198 hostnum = personality.Network_Config.Acoord;
199 hostnum *= personality.Network_Config.Bnodes;
200 hostnum += personality.Network_Config.Bcoord;
201 hostnum *= personality.Network_Config.Cnodes;
202 hostnum += personality.Network_Config.Ccoord;
203 hostnum *= personality.Network_Config.Dnodes;
204 hostnum += personality.Network_Config.Dcoord;
205 hostnum *= personality.Network_Config.Enodes;
206 hostnum += personality.Network_Config.Ecoord;
208 /* This procedure can only differentiate nodes with host names
209 * that end on unique numbers.
213 /* Only parse the host name up to the first dot */
214 while (i < resultlen && mpi_hostname[i] != '.')
216 if (isdigit(mpi_hostname[i]))
218 hostnum_str[j++] = mpi_hostname[i];
222 hostnum_str[j] = '\0';
229 /* Use only the last 9 decimals, so we don't overflow an int */
230 hostnum = strtol(hostnum_str + max(0, j-9), NULL, 10);
236 fprintf(debug, "In gmx_hostname_num: hostname '%s', hostnum %d\n",
237 mpi_hostname, hostnum);
240 "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",
241 personality.Network_Config.Acoord,
242 personality.Network_Config.Anodes,
243 personality.Network_Config.Bcoord,
244 personality.Network_Config.Bnodes,
245 personality.Network_Config.Ccoord,
246 personality.Network_Config.Cnodes,
247 personality.Network_Config.Dcoord,
248 personality.Network_Config.Dnodes,
249 personality.Network_Config.Ecoord,
250 personality.Network_Config.Enodes,
251 Kernel_ProcessorCoreID(),
253 Kernel_ProcessorID(),
255 Kernel_ProcessorThreadID(),
264 void gmx_setup_nodecomm(FILE gmx_unused *fplog, t_commrec *cr)
267 int n, rank, hostnum, ng, ni;
269 /* Many MPI implementations do not optimize MPI_Allreduce
270 * (and probably also other global communication calls)
271 * for multi-core nodes connected by a network.
272 * We can optimize such communication by using one MPI call
273 * within each node and one between the nodes.
274 * For MVAPICH2 and Intel MPI this reduces the time for
275 * the global_stat communication by 25%
276 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
277 * B. Hess, November 2007
283 #ifndef GMX_THREAD_MPI
285 MPI_Comm_size(cr->mpi_comm_mygroup, &n);
286 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
288 hostnum = gmx_hostname_num();
292 fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
296 /* The intra-node communicator, split on node number */
297 MPI_Comm_split(cr->mpi_comm_mygroup, hostnum, rank, &nc->comm_intra);
298 MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
301 fprintf(debug, "In gmx_setup_nodecomm: node rank %d rank_intra %d\n",
302 rank, nc->rank_intra);
304 /* The inter-node communicator, split on rank_intra.
305 * We actually only need the one for rank=0,
306 * but it is easier to create them all.
308 MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
309 /* Check if this really created two step communication */
310 MPI_Comm_size(nc->comm_inter, &ng);
311 MPI_Comm_size(nc->comm_intra, &ni);
314 fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
318 if (getenv("GMX_NO_NODECOMM") == NULL &&
319 ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
324 fprintf(fplog, "Using two step summing over %d groups of on average %.1f processes\n\n",
325 ng, (real)n/(real)ng);
327 if (nc->rank_intra > 0)
329 MPI_Comm_free(&nc->comm_inter);
334 /* One group or all processes in a separate group, use normal summing */
335 MPI_Comm_free(&nc->comm_inter);
336 MPI_Comm_free(&nc->comm_intra);
339 fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
344 /* tMPI runs only on a single node so just use the nodeid */
345 nc->rank_intra = cr->nodeid;
349 void gmx_init_intranode_counters(t_commrec *cr)
351 /* counters for PP+PME and PP-only processes on my physical node */
352 int nrank_intranode, rank_intranode;
353 int nrank_pp_intranode, rank_pp_intranode;
354 /* thread-MPI is not initialized when not running in parallel */
355 #if defined GMX_MPI && !defined GMX_THREAD_MPI
356 int nrank_world, rank_world;
357 int i, mynum, *num, *num_s, *num_pp, *num_pp_s;
359 MPI_Comm_size(MPI_COMM_WORLD, &nrank_world);
360 MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
362 /* Get the node number from the hostname to identify the nodes */
363 mynum = gmx_hostname_num();
365 /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
366 snew(num, nrank_world);
367 snew(num_s, nrank_world);
368 snew(num_pp, nrank_world);
369 snew(num_pp_s, nrank_world);
371 num_s[rank_world] = mynum;
372 num_pp_s[rank_world] = (cr->duty & DUTY_PP) ? mynum : -1;
374 MPI_Allreduce(num_s, num, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
375 MPI_Allreduce(num_pp_s, num_pp, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
379 nrank_pp_intranode = 0;
380 rank_pp_intranode = 0;
381 for (i = 0; i < nrank_world; i++)
391 if ((cr->duty & DUTY_PP) && num_pp[i] == mynum)
393 nrank_pp_intranode++;
405 /* Serial or thread-MPI code: we run within a single physical node */
406 nrank_intranode = cr->nnodes;
407 rank_intranode = cr->sim_nodeid;
408 nrank_pp_intranode = cr->nnodes - cr->npmenodes;
409 rank_pp_intranode = cr->nodeid;
415 if (cr->duty & DUTY_PP && cr->duty & DUTY_PME)
417 sprintf(sbuf, "PP+PME");
421 sprintf(sbuf, "%s", cr->duty & DUTY_PP ? "PP" : "PME");
423 fprintf(debug, "On %3s node %d: nrank_intranode=%d, rank_intranode=%d, "
424 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
425 sbuf, cr->sim_nodeid,
426 nrank_intranode, rank_intranode,
427 nrank_pp_intranode, rank_pp_intranode);
430 cr->nrank_intranode = nrank_intranode;
431 cr->rank_intranode = rank_intranode;
432 cr->nrank_pp_intranode = nrank_pp_intranode;
433 cr->rank_pp_intranode = rank_pp_intranode;
437 void gmx_barrier(const t_commrec *cr)
440 gmx_call("gmx_barrier");
442 MPI_Barrier(cr->mpi_comm_mygroup);
446 void gmx_abort(int gmx_unused noderank, int gmx_unused nnodes, int gmx_unused errorno)
449 gmx_call("gmx_abort");
451 #ifdef GMX_THREAD_MPI
452 fprintf(stderr, "Halting program %s\n", ShortProgram());
458 fprintf(stderr, "Halting parallel program %s on CPU %d out of %d\n",
459 ShortProgram(), noderank, nnodes);
463 fprintf(stderr, "Halting program %s\n", ShortProgram());
467 MPI_Abort(MPI_COMM_WORLD, errorno);
473 void gmx_bcast(int nbytes, void *b, const t_commrec *cr)
476 gmx_call("gmx_bast");
478 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
482 void gmx_bcast_sim(int nbytes, void *b, const t_commrec *cr)
485 gmx_call("gmx_bast");
487 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
491 void gmx_sumd(int nr, double r[], const t_commrec *cr)
494 gmx_call("gmx_sumd");
496 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
499 if (cr->nc.rank_intra == 0)
501 /* Use two step summing. */
502 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
504 /* Sum the roots of the internal (intra) buffers. */
505 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
510 /* This is here because of the silly MPI specification
511 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
512 MPI_Reduce(r, NULL, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
514 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
518 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
519 cr->mpi_comm_mygroup);
524 if (nr > cr->mpb->dbuf_alloc)
526 cr->mpb->dbuf_alloc = nr;
527 srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
531 /* Use two step summing */
532 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
533 if (cr->nc.rank_intra == 0)
535 /* Sum with the buffers reversed */
536 MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
539 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
543 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
544 cr->mpi_comm_mygroup);
545 for (i = 0; i < nr; i++)
547 r[i] = cr->mpb->dbuf[i];
554 void gmx_sumf(int nr, float r[], const t_commrec *cr)
557 gmx_call("gmx_sumf");
559 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
562 /* Use two step summing. */
563 if (cr->nc.rank_intra == 0)
565 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
567 /* Sum the roots of the internal (intra) buffers */
568 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
573 /* This is here because of the silly MPI specification
574 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
575 MPI_Reduce(r, NULL, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
577 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
581 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
586 if (nr > cr->mpb->fbuf_alloc)
588 cr->mpb->fbuf_alloc = nr;
589 srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
593 /* Use two step summing */
594 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
595 if (cr->nc.rank_intra == 0)
597 /* Sum with the buffers reversed */
598 MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
601 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
605 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
606 cr->mpi_comm_mygroup);
607 for (i = 0; i < nr; i++)
609 r[i] = cr->mpb->fbuf[i];
616 void gmx_sumi(int nr, int r[], const t_commrec *cr)
619 gmx_call("gmx_sumi");
621 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
624 /* Use two step summing */
625 if (cr->nc.rank_intra == 0)
627 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
628 /* Sum with the buffers reversed */
629 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
633 /* This is here because of the silly MPI specification
634 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
635 MPI_Reduce(r, NULL, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
637 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
641 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
646 if (nr > cr->mpb->ibuf_alloc)
648 cr->mpb->ibuf_alloc = nr;
649 srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
653 /* Use two step summing */
654 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
655 if (cr->nc.rank_intra == 0)
657 /* Sum with the buffers reversed */
658 MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
660 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
664 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
665 for (i = 0; i < nr; i++)
667 r[i] = cr->mpb->ibuf[i];
674 void gmx_sumli(int nr, gmx_large_int_t r[], const t_commrec *cr)
677 gmx_call("gmx_sumli");
679 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
682 /* Use two step summing */
683 if (cr->nc.rank_intra == 0)
685 MPI_Reduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM, 0,
687 /* Sum with the buffers reversed */
688 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
693 /* This is here because of the silly MPI specification
694 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
695 MPI_Reduce(r, NULL, nr, GMX_MPI_LARGE_INT, MPI_SUM, 0, cr->nc.comm_intra);
697 MPI_Bcast(r, nr, GMX_MPI_LARGE_INT, 0, cr->nc.comm_intra);
701 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM, cr->mpi_comm_mygroup);
706 if (nr > cr->mpb->libuf_alloc)
708 cr->mpb->libuf_alloc = nr;
709 srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
713 /* Use two step summing */
714 MPI_Allreduce(r, cr->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
716 if (cr->nc.rank_intra == 0)
718 /* Sum with the buffers reversed */
719 MPI_Allreduce(cr->mpb->libuf, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
722 MPI_Bcast(r, nr, GMX_MPI_LARGE_INT, 0, cr->nc.comm_intra);
726 MPI_Allreduce(r, cr->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
727 cr->mpi_comm_mygroup);
728 for (i = 0; i < nr; i++)
730 r[i] = cr->mpb->libuf[i];
740 void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
742 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
743 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
745 /* this function is only used in code that is not performance critical,
746 (during setup, when comm_rec is not the appropriate communication
747 structure), so this isn't as bad as it looks. */
752 MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
753 for (i = 0; i < nr; i++)
763 void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
765 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
766 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
768 /* this function is only used in code that is not performance critical,
769 (during setup, when comm_rec is not the appropriate communication
770 structure), so this isn't as bad as it looks. */
775 MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
776 for (i = 0; i < nr; i++)
785 void gmx_sumd_sim(int nr, double r[], const gmx_multisim_t *ms)
788 gmx_call("gmx_sumd_sim");
790 gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
794 void gmx_sumf_sim(int nr, float r[], const gmx_multisim_t *ms)
797 gmx_call("gmx_sumf_sim");
799 gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
803 void gmx_sumi_sim(int nr, int r[], const gmx_multisim_t *ms)
806 gmx_call("gmx_sumi_sim");
808 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
809 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
811 /* this is thread-unsafe, but it will do for now: */
814 if (nr > ms->mpb->ibuf_alloc)
816 ms->mpb->ibuf_alloc = nr;
817 srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
819 MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
820 for (i = 0; i < nr; i++)
822 r[i] = ms->mpb->ibuf[i];
828 void gmx_sumli_sim(int nr, gmx_large_int_t r[], const gmx_multisim_t *ms)
831 gmx_call("gmx_sumli_sim");
833 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
834 MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
835 ms->mpi_comm_masters);
837 /* this is thread-unsafe, but it will do for now: */
840 if (nr > ms->mpb->libuf_alloc)
842 ms->mpb->libuf_alloc = nr;
843 srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
845 MPI_Allreduce(r, ms->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
846 ms->mpi_comm_masters);
847 for (i = 0; i < nr; i++)
849 r[i] = ms->mpb->libuf[i];