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,2015,2016,2017, 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.
48 #include "gromacs/commandline/filenm.h"
49 #include "gromacs/mdtypes/commrec.h"
50 #include "gromacs/utility/basenetwork.h"
51 #include "gromacs/utility/cstringutil.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/futil.h"
54 #include "gromacs/utility/gmxmpi.h"
55 #include "gromacs/utility/smalloc.h"
57 /* The source code in this file should be thread-safe.
58 Please keep it that way. */
60 void gmx_fill_commrec_from_mpi(t_commrec gmx_unused
*cr
)
63 gmx_call("gmx_fill_commrec_from_mpi");
65 if (!gmx_mpi_initialized())
67 gmx_comm("MPI has not been initialized properly");
70 cr
->nnodes
= gmx_node_num();
71 cr
->nodeid
= gmx_node_rank();
72 // TODO This communicator should be always available. Currently we
73 // make it multiple times, and keep it only when relevant. But the
74 // cost of an extra communicator is negligible in single-node
75 // cases (both thread-MPI and real MPI) case, and we need it in
76 // all multi-node MPI cases with more than one PP rank per node,
77 // with and without GPUs. By always having it available, we also
78 // don't need to protect calls to mpi_comm_physicalnode, etc.
79 if (PAR(cr
) || MULTISIM(cr
))
81 MPI_Comm_split(MPI_COMM_WORLD
, gmx_physicalnode_id_hash(), cr
->nodeid
, &cr
->mpi_comm_physicalnode
);
83 cr
->sim_nodeid
= cr
->nodeid
;
84 cr
->mpi_comm_mysim
= MPI_COMM_WORLD
;
85 cr
->mpi_comm_mygroup
= MPI_COMM_WORLD
;
90 t_commrec
*init_commrec()
97 gmx_fill_commrec_from_mpi(cr
);
99 cr
->mpi_comm_mysim
= nullptr;
100 cr
->mpi_comm_mygroup
= nullptr;
103 cr
->nodeid
= cr
->sim_nodeid
;
106 // TODO cr->duty should not be initialized here
107 cr
->duty
= (DUTY_PP
| DUTY_PME
);
109 #if GMX_MPI && !MPI_IN_PLACE_EXISTS
110 /* initialize the MPI_IN_PLACE replacement buffers */
112 cr
->mpb
->ibuf
= NULL
;
113 cr
->mpb
->libuf
= NULL
;
114 cr
->mpb
->fbuf
= NULL
;
115 cr
->mpb
->dbuf
= NULL
;
116 cr
->mpb
->ibuf_alloc
= 0;
117 cr
->mpb
->libuf_alloc
= 0;
118 cr
->mpb
->fbuf_alloc
= 0;
119 cr
->mpb
->dbuf_alloc
= 0;
125 static void done_mpi_in_place_buf(mpi_in_place_buf_t
*buf
)
137 void done_commrec(t_commrec
*cr
)
140 if (PAR(cr
) || MULTISIM(cr
))
142 MPI_Comm_free(&cr
->mpi_comm_physicalnode
);
145 if (nullptr != cr
->dd
)
148 // done_domdec(cr->dd);
150 if (nullptr != cr
->ms
)
152 done_mpi_in_place_buf(cr
->ms
->mpb
);
155 done_mpi_in_place_buf(cr
->mpb
);
159 t_commrec
*reinitialize_commrec_for_this_thread(const t_commrec gmx_unused
*cro
)
164 /* make a thread-specific commrec */
166 /* now copy the whole thing, so settings like the number of PME nodes
170 /* and we start setting our own thread-specific values for things */
171 gmx_fill_commrec_from_mpi(cr
);
173 // TODO cr->duty should not be initialized here
174 cr
->duty
= (DUTY_PP
| DUTY_PME
);
182 void gmx_setup_nodecomm(FILE gmx_unused
*fplog
, t_commrec
*cr
)
186 /* Many MPI implementations do not optimize MPI_Allreduce
187 * (and probably also other global communication calls)
188 * for multi-core nodes connected by a network.
189 * We can optimize such communication by using one MPI call
190 * within each node and one between the nodes.
191 * For MVAPICH2 and Intel MPI this reduces the time for
192 * the global_stat communication by 25%
193 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
194 * B. Hess, November 2007
204 MPI_Comm_size(cr
->mpi_comm_mygroup
, &n
);
205 MPI_Comm_rank(cr
->mpi_comm_mygroup
, &rank
);
207 int nodehash
= gmx_physicalnode_id_hash();
211 fprintf(debug
, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n
);
215 /* The intra-node communicator, split on node number */
216 MPI_Comm_split(cr
->mpi_comm_mygroup
, nodehash
, rank
, &nc
->comm_intra
);
217 MPI_Comm_rank(nc
->comm_intra
, &nc
->rank_intra
);
220 fprintf(debug
, "In gmx_setup_nodecomm: node ID %d rank within node %d\n",
221 rank
, nc
->rank_intra
);
223 /* The inter-node communicator, split on rank_intra.
224 * We actually only need the one for rank=0,
225 * but it is easier to create them all.
227 MPI_Comm_split(cr
->mpi_comm_mygroup
, nc
->rank_intra
, rank
, &nc
->comm_inter
);
228 /* Check if this really created two step communication */
231 MPI_Comm_size(nc
->comm_inter
, &ng
);
232 MPI_Comm_size(nc
->comm_intra
, &ni
);
235 fprintf(debug
, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
239 if (getenv("GMX_NO_NODECOMM") == nullptr &&
240 ((ng
> 1 && ng
< n
) || (ni
> 1 && ni
< n
)))
245 fprintf(fplog
, "Using two step summing over %d groups of on average %.1f ranks\n\n",
246 ng
, (real
)n
/(real
)ng
);
248 if (nc
->rank_intra
> 0)
250 MPI_Comm_free(&nc
->comm_inter
);
255 /* One group or all processes in a separate group, use normal summing */
256 MPI_Comm_free(&nc
->comm_inter
);
257 MPI_Comm_free(&nc
->comm_intra
);
260 fprintf(debug
, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
265 /* tMPI runs only on a single node so just use the nodeid */
266 nc
->rank_intra
= cr
->nodeid
;
270 void gmx_init_intranode_counters(t_commrec
*cr
)
272 /* counters for PP+PME and PP-only processes on my physical node */
273 int nrank_intranode
, rank_intranode
;
274 int nrank_pp_intranode
, rank_pp_intranode
;
275 /* thread-MPI is not initialized when not running in parallel */
276 #if GMX_MPI && !GMX_THREAD_MPI
277 int nrank_world
, rank_world
;
278 int i
, myhash
, *hash
, *hash_s
, *hash_pp
, *hash_pp_s
;
280 MPI_Comm_size(MPI_COMM_WORLD
, &nrank_world
);
281 MPI_Comm_rank(MPI_COMM_WORLD
, &rank_world
);
283 /* Get a (hopefully unique) hash that identifies our physical node */
284 myhash
= gmx_physicalnode_id_hash();
286 /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
287 snew(hash
, nrank_world
);
288 snew(hash_s
, nrank_world
);
289 snew(hash_pp
, nrank_world
);
290 snew(hash_pp_s
, nrank_world
);
292 hash_s
[rank_world
] = myhash
;
293 hash_pp_s
[rank_world
] = thisRankHasDuty(cr
, DUTY_PP
) ? myhash
: -1;
295 MPI_Allreduce(hash_s
, hash
, nrank_world
, MPI_INT
, MPI_SUM
, MPI_COMM_WORLD
);
296 MPI_Allreduce(hash_pp_s
, hash_pp
, nrank_world
, MPI_INT
, MPI_SUM
, MPI_COMM_WORLD
);
300 nrank_pp_intranode
= 0;
301 rank_pp_intranode
= 0;
302 for (i
= 0; i
< nrank_world
; i
++)
304 if (hash
[i
] == myhash
)
312 if (hash_pp
[i
] == myhash
)
314 nrank_pp_intranode
++;
315 if (thisRankHasDuty(cr
, DUTY_PP
) && i
< rank_world
)
326 /* Serial or thread-MPI code: we run within a single physical node */
327 nrank_intranode
= cr
->nnodes
;
328 rank_intranode
= cr
->sim_nodeid
;
329 nrank_pp_intranode
= cr
->nnodes
- cr
->npmenodes
;
330 rank_pp_intranode
= cr
->nodeid
;
336 if (thisRankHasDuty(cr
, DUTY_PP
) && thisRankHasDuty(cr
, DUTY_PME
))
338 sprintf(sbuf
, "PP+PME");
342 sprintf(sbuf
, "%s", thisRankHasDuty(cr
, DUTY_PP
) ? "PP" : "PME");
344 fprintf(debug
, "On %3s rank %d: nrank_intranode=%d, rank_intranode=%d, "
345 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
346 sbuf
, cr
->sim_nodeid
,
347 nrank_intranode
, rank_intranode
,
348 nrank_pp_intranode
, rank_pp_intranode
);
351 cr
->nrank_intranode
= nrank_intranode
;
352 cr
->rank_intranode
= rank_intranode
;
353 cr
->nrank_pp_intranode
= nrank_pp_intranode
;
354 cr
->rank_pp_intranode
= rank_pp_intranode
;
358 void gmx_barrier(const t_commrec gmx_unused
*cr
)
361 gmx_call("gmx_barrier");
363 MPI_Barrier(cr
->mpi_comm_mygroup
);
367 void gmx_barrier_physical_node(const t_commrec gmx_unused
*cr
)
370 gmx_call("gmx_barrier_physical_node");
372 MPI_Barrier(cr
->mpi_comm_physicalnode
);
376 void gmx_bcast(int gmx_unused nbytes
, void gmx_unused
*b
, const t_commrec gmx_unused
*cr
)
379 gmx_call("gmx_bast");
381 MPI_Bcast(b
, nbytes
, MPI_BYTE
, MASTERRANK(cr
), cr
->mpi_comm_mygroup
);
385 void gmx_bcast_sim(int gmx_unused nbytes
, void gmx_unused
*b
, const t_commrec gmx_unused
*cr
)
388 gmx_call("gmx_bast");
390 MPI_Bcast(b
, nbytes
, MPI_BYTE
, MASTERRANK(cr
), cr
->mpi_comm_mysim
);
394 void gmx_sumd(int gmx_unused nr
, double gmx_unused r
[], const t_commrec gmx_unused
*cr
)
397 gmx_call("gmx_sumd");
399 #if MPI_IN_PLACE_EXISTS
402 if (cr
->nc
.rank_intra
== 0)
404 /* Use two step summing. */
405 MPI_Reduce(MPI_IN_PLACE
, r
, nr
, MPI_DOUBLE
, MPI_SUM
, 0,
407 /* Sum the roots of the internal (intra) buffers. */
408 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_DOUBLE
, MPI_SUM
,
413 /* This is here because of the silly MPI specification
414 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
415 MPI_Reduce(r
, nullptr, nr
, MPI_DOUBLE
, MPI_SUM
, 0, cr
->nc
.comm_intra
);
417 MPI_Bcast(r
, nr
, MPI_DOUBLE
, 0, cr
->nc
.comm_intra
);
421 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_DOUBLE
, MPI_SUM
,
422 cr
->mpi_comm_mygroup
);
427 if (nr
> cr
->mpb
->dbuf_alloc
)
429 cr
->mpb
->dbuf_alloc
= nr
;
430 srenew(cr
->mpb
->dbuf
, cr
->mpb
->dbuf_alloc
);
434 /* Use two step summing */
435 MPI_Allreduce(r
, cr
->mpb
->dbuf
, nr
, MPI_DOUBLE
, MPI_SUM
, cr
->nc
.comm_intra
);
436 if (cr
->nc
.rank_intra
== 0)
438 /* Sum with the buffers reversed */
439 MPI_Allreduce(cr
->mpb
->dbuf
, r
, nr
, MPI_DOUBLE
, MPI_SUM
,
442 MPI_Bcast(r
, nr
, MPI_DOUBLE
, 0, cr
->nc
.comm_intra
);
446 MPI_Allreduce(r
, cr
->mpb
->dbuf
, nr
, MPI_DOUBLE
, MPI_SUM
,
447 cr
->mpi_comm_mygroup
);
448 for (i
= 0; i
< nr
; i
++)
450 r
[i
] = cr
->mpb
->dbuf
[i
];
457 void gmx_sumf(int gmx_unused nr
, float gmx_unused r
[], const t_commrec gmx_unused
*cr
)
460 gmx_call("gmx_sumf");
462 #if MPI_IN_PLACE_EXISTS
465 /* Use two step summing. */
466 if (cr
->nc
.rank_intra
== 0)
468 MPI_Reduce(MPI_IN_PLACE
, r
, nr
, MPI_FLOAT
, MPI_SUM
, 0,
470 /* Sum the roots of the internal (intra) buffers */
471 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_FLOAT
, MPI_SUM
,
476 /* This is here because of the silly MPI specification
477 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
478 MPI_Reduce(r
, nullptr, nr
, MPI_FLOAT
, MPI_SUM
, 0, cr
->nc
.comm_intra
);
480 MPI_Bcast(r
, nr
, MPI_FLOAT
, 0, cr
->nc
.comm_intra
);
484 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_FLOAT
, MPI_SUM
, cr
->mpi_comm_mygroup
);
489 if (nr
> cr
->mpb
->fbuf_alloc
)
491 cr
->mpb
->fbuf_alloc
= nr
;
492 srenew(cr
->mpb
->fbuf
, cr
->mpb
->fbuf_alloc
);
496 /* Use two step summing */
497 MPI_Allreduce(r
, cr
->mpb
->fbuf
, nr
, MPI_FLOAT
, MPI_SUM
, cr
->nc
.comm_intra
);
498 if (cr
->nc
.rank_intra
== 0)
500 /* Sum with the buffers reversed */
501 MPI_Allreduce(cr
->mpb
->fbuf
, r
, nr
, MPI_FLOAT
, MPI_SUM
,
504 MPI_Bcast(r
, nr
, MPI_FLOAT
, 0, cr
->nc
.comm_intra
);
508 MPI_Allreduce(r
, cr
->mpb
->fbuf
, nr
, MPI_FLOAT
, MPI_SUM
,
509 cr
->mpi_comm_mygroup
);
510 for (i
= 0; i
< nr
; i
++)
512 r
[i
] = cr
->mpb
->fbuf
[i
];
519 void gmx_sumi(int gmx_unused nr
, int gmx_unused r
[], const t_commrec gmx_unused
*cr
)
522 gmx_call("gmx_sumi");
524 #if MPI_IN_PLACE_EXISTS
527 /* Use two step summing */
528 if (cr
->nc
.rank_intra
== 0)
530 MPI_Reduce(MPI_IN_PLACE
, r
, nr
, MPI_INT
, MPI_SUM
, 0, cr
->nc
.comm_intra
);
531 /* Sum with the buffers reversed */
532 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_INT
, MPI_SUM
, cr
->nc
.comm_inter
);
536 /* This is here because of the silly MPI specification
537 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
538 MPI_Reduce(r
, nullptr, nr
, MPI_INT
, MPI_SUM
, 0, cr
->nc
.comm_intra
);
540 MPI_Bcast(r
, nr
, MPI_INT
, 0, cr
->nc
.comm_intra
);
544 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_INT
, MPI_SUM
, cr
->mpi_comm_mygroup
);
549 if (nr
> cr
->mpb
->ibuf_alloc
)
551 cr
->mpb
->ibuf_alloc
= nr
;
552 srenew(cr
->mpb
->ibuf
, cr
->mpb
->ibuf_alloc
);
556 /* Use two step summing */
557 MPI_Allreduce(r
, cr
->mpb
->ibuf
, nr
, MPI_INT
, MPI_SUM
, cr
->nc
.comm_intra
);
558 if (cr
->nc
.rank_intra
== 0)
560 /* Sum with the buffers reversed */
561 MPI_Allreduce(cr
->mpb
->ibuf
, r
, nr
, MPI_INT
, MPI_SUM
, cr
->nc
.comm_inter
);
563 MPI_Bcast(r
, nr
, MPI_INT
, 0, cr
->nc
.comm_intra
);
567 MPI_Allreduce(r
, cr
->mpb
->ibuf
, nr
, MPI_INT
, MPI_SUM
, cr
->mpi_comm_mygroup
);
568 for (i
= 0; i
< nr
; i
++)
570 r
[i
] = cr
->mpb
->ibuf
[i
];
577 void gmx_sumli(int gmx_unused nr
, gmx_int64_t gmx_unused r
[], const t_commrec gmx_unused
*cr
)
580 gmx_call("gmx_sumli");
582 #if MPI_IN_PLACE_EXISTS
585 /* Use two step summing */
586 if (cr
->nc
.rank_intra
== 0)
588 MPI_Reduce(MPI_IN_PLACE
, r
, nr
, MPI_INT64_T
, MPI_SUM
, 0,
590 /* Sum with the buffers reversed */
591 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_INT64_T
, MPI_SUM
,
596 /* This is here because of the silly MPI specification
597 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
598 MPI_Reduce(r
, nullptr, nr
, MPI_INT64_T
, MPI_SUM
, 0, cr
->nc
.comm_intra
);
600 MPI_Bcast(r
, nr
, MPI_INT64_T
, 0, cr
->nc
.comm_intra
);
604 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_INT64_T
, MPI_SUM
, cr
->mpi_comm_mygroup
);
609 if (nr
> cr
->mpb
->libuf_alloc
)
611 cr
->mpb
->libuf_alloc
= nr
;
612 srenew(cr
->mpb
->libuf
, cr
->mpb
->libuf_alloc
);
616 /* Use two step summing */
617 MPI_Allreduce(r
, cr
->mpb
->libuf
, nr
, MPI_INT64_T
, MPI_SUM
,
619 if (cr
->nc
.rank_intra
== 0)
621 /* Sum with the buffers reversed */
622 MPI_Allreduce(cr
->mpb
->libuf
, r
, nr
, MPI_INT64_T
, MPI_SUM
,
625 MPI_Bcast(r
, nr
, MPI_INT64_T
, 0, cr
->nc
.comm_intra
);
629 MPI_Allreduce(r
, cr
->mpb
->libuf
, nr
, MPI_INT64_T
, MPI_SUM
,
630 cr
->mpi_comm_mygroup
);
631 for (i
= 0; i
< nr
; i
++)
633 r
[i
] = cr
->mpb
->libuf
[i
];
643 static void gmx_sumd_comm(int nr
, double r
[], MPI_Comm mpi_comm
)
645 #if MPI_IN_PLACE_EXISTS
646 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_DOUBLE
, MPI_SUM
, mpi_comm
);
648 /* this function is only used in code that is not performance critical,
649 (during setup, when comm_rec is not the appropriate communication
650 structure), so this isn't as bad as it looks. */
655 MPI_Allreduce(r
, buf
, nr
, MPI_DOUBLE
, MPI_SUM
, mpi_comm
);
656 for (i
= 0; i
< nr
; i
++)
666 static void gmx_sumf_comm(int nr
, float r
[], MPI_Comm mpi_comm
)
668 #if MPI_IN_PLACE_EXISTS
669 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_FLOAT
, MPI_SUM
, mpi_comm
);
671 /* this function is only used in code that is not performance critical,
672 (during setup, when comm_rec is not the appropriate communication
673 structure), so this isn't as bad as it looks. */
678 MPI_Allreduce(r
, buf
, nr
, MPI_FLOAT
, MPI_SUM
, mpi_comm
);
679 for (i
= 0; i
< nr
; i
++)
688 void gmx_sumd_sim(int gmx_unused nr
, double gmx_unused r
[], const gmx_multisim_t gmx_unused
*ms
)
691 gmx_call("gmx_sumd_sim");
693 gmx_sumd_comm(nr
, r
, ms
->mpi_comm_masters
);
697 void gmx_sumf_sim(int gmx_unused nr
, float gmx_unused r
[], const gmx_multisim_t gmx_unused
*ms
)
700 gmx_call("gmx_sumf_sim");
702 gmx_sumf_comm(nr
, r
, ms
->mpi_comm_masters
);
706 void gmx_sumi_sim(int gmx_unused nr
, int gmx_unused r
[], const gmx_multisim_t gmx_unused
*ms
)
709 gmx_call("gmx_sumi_sim");
711 #if MPI_IN_PLACE_EXISTS
712 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_INT
, MPI_SUM
, ms
->mpi_comm_masters
);
714 /* this is thread-unsafe, but it will do for now: */
717 if (nr
> ms
->mpb
->ibuf_alloc
)
719 ms
->mpb
->ibuf_alloc
= nr
;
720 srenew(ms
->mpb
->ibuf
, ms
->mpb
->ibuf_alloc
);
722 MPI_Allreduce(r
, ms
->mpb
->ibuf
, nr
, MPI_INT
, MPI_SUM
, ms
->mpi_comm_masters
);
723 for (i
= 0; i
< nr
; i
++)
725 r
[i
] = ms
->mpb
->ibuf
[i
];
731 void gmx_sumli_sim(int gmx_unused nr
, gmx_int64_t gmx_unused r
[], const gmx_multisim_t gmx_unused
*ms
)
734 gmx_call("gmx_sumli_sim");
736 #if MPI_IN_PLACE_EXISTS
737 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_INT64_T
, MPI_SUM
,
738 ms
->mpi_comm_masters
);
740 /* this is thread-unsafe, but it will do for now: */
743 if (nr
> ms
->mpb
->libuf_alloc
)
745 ms
->mpb
->libuf_alloc
= nr
;
746 srenew(ms
->mpb
->libuf
, ms
->mpb
->libuf_alloc
);
748 MPI_Allreduce(r
, ms
->mpb
->libuf
, nr
, MPI_INT64_T
, MPI_SUM
,
749 ms
->mpi_comm_masters
);
750 for (i
= 0; i
< nr
; i
++)
752 r
[i
] = ms
->mpb
->libuf
[i
];
758 const char *opt2fn_master(const char *opt
, int nfile
, const t_filenm fnm
[],
761 return SIMMASTER(cr
) ? opt2fn(opt
, nfile
, fnm
) : nullptr;
764 void gmx_fatal_collective(int f_errno
, const char *file
, int line
,
765 MPI_Comm comm
, gmx_bool bMaster
,
766 const char *fmt
, ...)
772 /* Check if we are calling on all processes in MPI_COMM_WORLD */
773 MPI_Comm_compare(comm
, MPI_COMM_WORLD
, &result
);
774 /* Any result except MPI_UNEQUAL allows us to call MPI_Finalize */
775 bFinalize
= (result
!= MPI_UNEQUAL
);
777 GMX_UNUSED_VALUE(comm
);
782 gmx_fatal_mpi_va(f_errno
, file
, line
, bMaster
, bFinalize
, fmt
, ap
);