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 if (PAR(cr
) || MULTISIM(cr
))
74 MPI_Comm_split(MPI_COMM_WORLD
, gmx_physicalnode_id_hash(), cr
->nodeid
, &cr
->mpi_comm_physicalnode
);
76 cr
->sim_nodeid
= cr
->nodeid
;
77 cr
->mpi_comm_mysim
= MPI_COMM_WORLD
;
78 cr
->mpi_comm_mygroup
= MPI_COMM_WORLD
;
83 t_commrec
*init_commrec()
90 gmx_fill_commrec_from_mpi(cr
);
92 cr
->mpi_comm_mysim
= nullptr;
93 cr
->mpi_comm_mygroup
= nullptr;
96 cr
->nodeid
= cr
->sim_nodeid
;
99 // TODO cr->duty should not be initialized here
100 cr
->duty
= (DUTY_PP
| DUTY_PME
);
102 #if GMX_MPI && !MPI_IN_PLACE_EXISTS
103 /* initialize the MPI_IN_PLACE replacement buffers */
105 cr
->mpb
->ibuf
= NULL
;
106 cr
->mpb
->libuf
= NULL
;
107 cr
->mpb
->fbuf
= NULL
;
108 cr
->mpb
->dbuf
= NULL
;
109 cr
->mpb
->ibuf_alloc
= 0;
110 cr
->mpb
->libuf_alloc
= 0;
111 cr
->mpb
->fbuf_alloc
= 0;
112 cr
->mpb
->dbuf_alloc
= 0;
118 static void done_mpi_in_place_buf(mpi_in_place_buf_t
*buf
)
130 void done_commrec(t_commrec
*cr
)
133 if (PAR(cr
) || MULTISIM(cr
))
135 MPI_Comm_free(&cr
->mpi_comm_physicalnode
);
138 if (nullptr != cr
->dd
)
141 // done_domdec(cr->dd);
143 if (nullptr != cr
->ms
)
145 done_mpi_in_place_buf(cr
->ms
->mpb
);
148 done_mpi_in_place_buf(cr
->mpb
);
152 t_commrec
*reinitialize_commrec_for_this_thread(const t_commrec gmx_unused
*cro
)
157 /* make a thread-specific commrec */
159 /* now copy the whole thing, so settings like the number of PME nodes
163 /* and we start setting our own thread-specific values for things */
164 gmx_fill_commrec_from_mpi(cr
);
166 // TODO cr->duty should not be initialized here
167 cr
->duty
= (DUTY_PP
| DUTY_PME
);
175 void gmx_setup_nodecomm(FILE gmx_unused
*fplog
, t_commrec
*cr
)
179 /* Many MPI implementations do not optimize MPI_Allreduce
180 * (and probably also other global communication calls)
181 * for multi-core nodes connected by a network.
182 * We can optimize such communication by using one MPI call
183 * within each node and one between the nodes.
184 * For MVAPICH2 and Intel MPI this reduces the time for
185 * the global_stat communication by 25%
186 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
187 * B. Hess, November 2007
197 MPI_Comm_size(cr
->mpi_comm_mygroup
, &n
);
198 MPI_Comm_rank(cr
->mpi_comm_mygroup
, &rank
);
200 int nodehash
= gmx_physicalnode_id_hash();
204 fprintf(debug
, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n
);
208 /* The intra-node communicator, split on node number */
209 MPI_Comm_split(cr
->mpi_comm_mygroup
, nodehash
, rank
, &nc
->comm_intra
);
210 MPI_Comm_rank(nc
->comm_intra
, &nc
->rank_intra
);
213 fprintf(debug
, "In gmx_setup_nodecomm: node ID %d rank within node %d\n",
214 rank
, nc
->rank_intra
);
216 /* The inter-node communicator, split on rank_intra.
217 * We actually only need the one for rank=0,
218 * but it is easier to create them all.
220 MPI_Comm_split(cr
->mpi_comm_mygroup
, nc
->rank_intra
, rank
, &nc
->comm_inter
);
221 /* Check if this really created two step communication */
224 MPI_Comm_size(nc
->comm_inter
, &ng
);
225 MPI_Comm_size(nc
->comm_intra
, &ni
);
228 fprintf(debug
, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
232 if (getenv("GMX_NO_NODECOMM") == nullptr &&
233 ((ng
> 1 && ng
< n
) || (ni
> 1 && ni
< n
)))
238 fprintf(fplog
, "Using two step summing over %d groups of on average %.1f ranks\n\n",
239 ng
, (real
)n
/(real
)ng
);
241 if (nc
->rank_intra
> 0)
243 MPI_Comm_free(&nc
->comm_inter
);
248 /* One group or all processes in a separate group, use normal summing */
249 MPI_Comm_free(&nc
->comm_inter
);
250 MPI_Comm_free(&nc
->comm_intra
);
253 fprintf(debug
, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
258 /* tMPI runs only on a single node so just use the nodeid */
259 nc
->rank_intra
= cr
->nodeid
;
263 void gmx_init_intranode_counters(t_commrec
*cr
)
265 /* counters for PP+PME and PP-only processes on my physical node */
266 int nrank_intranode
, rank_intranode
;
267 int nrank_pp_intranode
, rank_pp_intranode
;
268 /* thread-MPI is not initialized when not running in parallel */
269 #if GMX_MPI && !GMX_THREAD_MPI
270 int nrank_world
, rank_world
;
271 int i
, myhash
, *hash
, *hash_s
, *hash_pp
, *hash_pp_s
;
273 MPI_Comm_size(MPI_COMM_WORLD
, &nrank_world
);
274 MPI_Comm_rank(MPI_COMM_WORLD
, &rank_world
);
276 /* Get a (hopefully unique) hash that identifies our physical node */
277 myhash
= gmx_physicalnode_id_hash();
279 /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
280 snew(hash
, nrank_world
);
281 snew(hash_s
, nrank_world
);
282 snew(hash_pp
, nrank_world
);
283 snew(hash_pp_s
, nrank_world
);
285 hash_s
[rank_world
] = myhash
;
286 hash_pp_s
[rank_world
] = (cr
->duty
& DUTY_PP
) ? myhash
: -1;
288 MPI_Allreduce(hash_s
, hash
, nrank_world
, MPI_INT
, MPI_SUM
, MPI_COMM_WORLD
);
289 MPI_Allreduce(hash_pp_s
, hash_pp
, nrank_world
, MPI_INT
, MPI_SUM
, MPI_COMM_WORLD
);
293 nrank_pp_intranode
= 0;
294 rank_pp_intranode
= 0;
295 for (i
= 0; i
< nrank_world
; i
++)
297 if (hash
[i
] == myhash
)
305 if (hash_pp
[i
] == myhash
)
307 nrank_pp_intranode
++;
308 if ((cr
->duty
& DUTY_PP
) && i
< rank_world
)
319 /* Serial or thread-MPI code: we run within a single physical node */
320 nrank_intranode
= cr
->nnodes
;
321 rank_intranode
= cr
->sim_nodeid
;
322 nrank_pp_intranode
= cr
->nnodes
- cr
->npmenodes
;
323 rank_pp_intranode
= cr
->nodeid
;
329 if ((cr
->duty
& DUTY_PP
) && (cr
->duty
& DUTY_PME
))
331 sprintf(sbuf
, "PP+PME");
335 sprintf(sbuf
, "%s", (cr
->duty
& DUTY_PP
) ? "PP" : "PME");
337 fprintf(debug
, "On %3s rank %d: nrank_intranode=%d, rank_intranode=%d, "
338 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
339 sbuf
, cr
->sim_nodeid
,
340 nrank_intranode
, rank_intranode
,
341 nrank_pp_intranode
, rank_pp_intranode
);
344 cr
->nrank_intranode
= nrank_intranode
;
345 cr
->rank_intranode
= rank_intranode
;
346 cr
->nrank_pp_intranode
= nrank_pp_intranode
;
347 cr
->rank_pp_intranode
= rank_pp_intranode
;
351 void gmx_barrier(const t_commrec gmx_unused
*cr
)
354 gmx_call("gmx_barrier");
356 MPI_Barrier(cr
->mpi_comm_mygroup
);
360 void gmx_barrier_physical_node(const t_commrec gmx_unused
*cr
)
363 gmx_call("gmx_barrier_physical_node");
365 MPI_Barrier(cr
->mpi_comm_physicalnode
);
369 void gmx_bcast(int gmx_unused nbytes
, void gmx_unused
*b
, const t_commrec gmx_unused
*cr
)
372 gmx_call("gmx_bast");
374 MPI_Bcast(b
, nbytes
, MPI_BYTE
, MASTERRANK(cr
), cr
->mpi_comm_mygroup
);
378 void gmx_bcast_sim(int gmx_unused nbytes
, void gmx_unused
*b
, const t_commrec gmx_unused
*cr
)
381 gmx_call("gmx_bast");
383 MPI_Bcast(b
, nbytes
, MPI_BYTE
, MASTERRANK(cr
), cr
->mpi_comm_mysim
);
387 void gmx_sumd(int gmx_unused nr
, double gmx_unused r
[], const t_commrec gmx_unused
*cr
)
390 gmx_call("gmx_sumd");
392 #if MPI_IN_PLACE_EXISTS
395 if (cr
->nc
.rank_intra
== 0)
397 /* Use two step summing. */
398 MPI_Reduce(MPI_IN_PLACE
, r
, nr
, MPI_DOUBLE
, MPI_SUM
, 0,
400 /* Sum the roots of the internal (intra) buffers. */
401 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_DOUBLE
, MPI_SUM
,
406 /* This is here because of the silly MPI specification
407 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
408 MPI_Reduce(r
, nullptr, nr
, MPI_DOUBLE
, MPI_SUM
, 0, cr
->nc
.comm_intra
);
410 MPI_Bcast(r
, nr
, MPI_DOUBLE
, 0, cr
->nc
.comm_intra
);
414 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_DOUBLE
, MPI_SUM
,
415 cr
->mpi_comm_mygroup
);
420 if (nr
> cr
->mpb
->dbuf_alloc
)
422 cr
->mpb
->dbuf_alloc
= nr
;
423 srenew(cr
->mpb
->dbuf
, cr
->mpb
->dbuf_alloc
);
427 /* Use two step summing */
428 MPI_Allreduce(r
, cr
->mpb
->dbuf
, nr
, MPI_DOUBLE
, MPI_SUM
, cr
->nc
.comm_intra
);
429 if (cr
->nc
.rank_intra
== 0)
431 /* Sum with the buffers reversed */
432 MPI_Allreduce(cr
->mpb
->dbuf
, r
, nr
, MPI_DOUBLE
, MPI_SUM
,
435 MPI_Bcast(r
, nr
, MPI_DOUBLE
, 0, cr
->nc
.comm_intra
);
439 MPI_Allreduce(r
, cr
->mpb
->dbuf
, nr
, MPI_DOUBLE
, MPI_SUM
,
440 cr
->mpi_comm_mygroup
);
441 for (i
= 0; i
< nr
; i
++)
443 r
[i
] = cr
->mpb
->dbuf
[i
];
450 void gmx_sumf(int gmx_unused nr
, float gmx_unused r
[], const t_commrec gmx_unused
*cr
)
453 gmx_call("gmx_sumf");
455 #if MPI_IN_PLACE_EXISTS
458 /* Use two step summing. */
459 if (cr
->nc
.rank_intra
== 0)
461 MPI_Reduce(MPI_IN_PLACE
, r
, nr
, MPI_FLOAT
, MPI_SUM
, 0,
463 /* Sum the roots of the internal (intra) buffers */
464 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_FLOAT
, MPI_SUM
,
469 /* This is here because of the silly MPI specification
470 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
471 MPI_Reduce(r
, nullptr, nr
, MPI_FLOAT
, MPI_SUM
, 0, cr
->nc
.comm_intra
);
473 MPI_Bcast(r
, nr
, MPI_FLOAT
, 0, cr
->nc
.comm_intra
);
477 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_FLOAT
, MPI_SUM
, cr
->mpi_comm_mygroup
);
482 if (nr
> cr
->mpb
->fbuf_alloc
)
484 cr
->mpb
->fbuf_alloc
= nr
;
485 srenew(cr
->mpb
->fbuf
, cr
->mpb
->fbuf_alloc
);
489 /* Use two step summing */
490 MPI_Allreduce(r
, cr
->mpb
->fbuf
, nr
, MPI_FLOAT
, MPI_SUM
, cr
->nc
.comm_intra
);
491 if (cr
->nc
.rank_intra
== 0)
493 /* Sum with the buffers reversed */
494 MPI_Allreduce(cr
->mpb
->fbuf
, r
, nr
, MPI_FLOAT
, MPI_SUM
,
497 MPI_Bcast(r
, nr
, MPI_FLOAT
, 0, cr
->nc
.comm_intra
);
501 MPI_Allreduce(r
, cr
->mpb
->fbuf
, nr
, MPI_FLOAT
, MPI_SUM
,
502 cr
->mpi_comm_mygroup
);
503 for (i
= 0; i
< nr
; i
++)
505 r
[i
] = cr
->mpb
->fbuf
[i
];
512 void gmx_sumi(int gmx_unused nr
, int gmx_unused r
[], const t_commrec gmx_unused
*cr
)
515 gmx_call("gmx_sumi");
517 #if MPI_IN_PLACE_EXISTS
520 /* Use two step summing */
521 if (cr
->nc
.rank_intra
== 0)
523 MPI_Reduce(MPI_IN_PLACE
, r
, nr
, MPI_INT
, MPI_SUM
, 0, cr
->nc
.comm_intra
);
524 /* Sum with the buffers reversed */
525 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_INT
, MPI_SUM
, cr
->nc
.comm_inter
);
529 /* This is here because of the silly MPI specification
530 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
531 MPI_Reduce(r
, nullptr, nr
, MPI_INT
, MPI_SUM
, 0, cr
->nc
.comm_intra
);
533 MPI_Bcast(r
, nr
, MPI_INT
, 0, cr
->nc
.comm_intra
);
537 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_INT
, MPI_SUM
, cr
->mpi_comm_mygroup
);
542 if (nr
> cr
->mpb
->ibuf_alloc
)
544 cr
->mpb
->ibuf_alloc
= nr
;
545 srenew(cr
->mpb
->ibuf
, cr
->mpb
->ibuf_alloc
);
549 /* Use two step summing */
550 MPI_Allreduce(r
, cr
->mpb
->ibuf
, nr
, MPI_INT
, MPI_SUM
, cr
->nc
.comm_intra
);
551 if (cr
->nc
.rank_intra
== 0)
553 /* Sum with the buffers reversed */
554 MPI_Allreduce(cr
->mpb
->ibuf
, r
, nr
, MPI_INT
, MPI_SUM
, cr
->nc
.comm_inter
);
556 MPI_Bcast(r
, nr
, MPI_INT
, 0, cr
->nc
.comm_intra
);
560 MPI_Allreduce(r
, cr
->mpb
->ibuf
, nr
, MPI_INT
, MPI_SUM
, cr
->mpi_comm_mygroup
);
561 for (i
= 0; i
< nr
; i
++)
563 r
[i
] = cr
->mpb
->ibuf
[i
];
570 void gmx_sumli(int gmx_unused nr
, gmx_int64_t gmx_unused r
[], const t_commrec gmx_unused
*cr
)
573 gmx_call("gmx_sumli");
575 #if MPI_IN_PLACE_EXISTS
578 /* Use two step summing */
579 if (cr
->nc
.rank_intra
== 0)
581 MPI_Reduce(MPI_IN_PLACE
, r
, nr
, MPI_INT64_T
, MPI_SUM
, 0,
583 /* Sum with the buffers reversed */
584 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_INT64_T
, MPI_SUM
,
589 /* This is here because of the silly MPI specification
590 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
591 MPI_Reduce(r
, nullptr, nr
, MPI_INT64_T
, MPI_SUM
, 0, cr
->nc
.comm_intra
);
593 MPI_Bcast(r
, nr
, MPI_INT64_T
, 0, cr
->nc
.comm_intra
);
597 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_INT64_T
, MPI_SUM
, cr
->mpi_comm_mygroup
);
602 if (nr
> cr
->mpb
->libuf_alloc
)
604 cr
->mpb
->libuf_alloc
= nr
;
605 srenew(cr
->mpb
->libuf
, cr
->mpb
->libuf_alloc
);
609 /* Use two step summing */
610 MPI_Allreduce(r
, cr
->mpb
->libuf
, nr
, MPI_INT64_T
, MPI_SUM
,
612 if (cr
->nc
.rank_intra
== 0)
614 /* Sum with the buffers reversed */
615 MPI_Allreduce(cr
->mpb
->libuf
, r
, nr
, MPI_INT64_T
, MPI_SUM
,
618 MPI_Bcast(r
, nr
, MPI_INT64_T
, 0, cr
->nc
.comm_intra
);
622 MPI_Allreduce(r
, cr
->mpb
->libuf
, nr
, MPI_INT64_T
, MPI_SUM
,
623 cr
->mpi_comm_mygroup
);
624 for (i
= 0; i
< nr
; i
++)
626 r
[i
] = cr
->mpb
->libuf
[i
];
636 static void gmx_sumd_comm(int nr
, double r
[], MPI_Comm mpi_comm
)
638 #if MPI_IN_PLACE_EXISTS
639 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_DOUBLE
, MPI_SUM
, mpi_comm
);
641 /* this function is only used in code that is not performance critical,
642 (during setup, when comm_rec is not the appropriate communication
643 structure), so this isn't as bad as it looks. */
648 MPI_Allreduce(r
, buf
, nr
, MPI_DOUBLE
, MPI_SUM
, mpi_comm
);
649 for (i
= 0; i
< nr
; i
++)
659 static void gmx_sumf_comm(int nr
, float r
[], MPI_Comm mpi_comm
)
661 #if MPI_IN_PLACE_EXISTS
662 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_FLOAT
, MPI_SUM
, mpi_comm
);
664 /* this function is only used in code that is not performance critical,
665 (during setup, when comm_rec is not the appropriate communication
666 structure), so this isn't as bad as it looks. */
671 MPI_Allreduce(r
, buf
, nr
, MPI_FLOAT
, MPI_SUM
, mpi_comm
);
672 for (i
= 0; i
< nr
; i
++)
681 void gmx_sumd_sim(int gmx_unused nr
, double gmx_unused r
[], const gmx_multisim_t gmx_unused
*ms
)
684 gmx_call("gmx_sumd_sim");
686 gmx_sumd_comm(nr
, r
, ms
->mpi_comm_masters
);
690 void gmx_sumf_sim(int gmx_unused nr
, float gmx_unused r
[], const gmx_multisim_t gmx_unused
*ms
)
693 gmx_call("gmx_sumf_sim");
695 gmx_sumf_comm(nr
, r
, ms
->mpi_comm_masters
);
699 void gmx_sumi_sim(int gmx_unused nr
, int gmx_unused r
[], const gmx_multisim_t gmx_unused
*ms
)
702 gmx_call("gmx_sumi_sim");
704 #if MPI_IN_PLACE_EXISTS
705 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_INT
, MPI_SUM
, ms
->mpi_comm_masters
);
707 /* this is thread-unsafe, but it will do for now: */
710 if (nr
> ms
->mpb
->ibuf_alloc
)
712 ms
->mpb
->ibuf_alloc
= nr
;
713 srenew(ms
->mpb
->ibuf
, ms
->mpb
->ibuf_alloc
);
715 MPI_Allreduce(r
, ms
->mpb
->ibuf
, nr
, MPI_INT
, MPI_SUM
, ms
->mpi_comm_masters
);
716 for (i
= 0; i
< nr
; i
++)
718 r
[i
] = ms
->mpb
->ibuf
[i
];
724 void gmx_sumli_sim(int gmx_unused nr
, gmx_int64_t gmx_unused r
[], const gmx_multisim_t gmx_unused
*ms
)
727 gmx_call("gmx_sumli_sim");
729 #if MPI_IN_PLACE_EXISTS
730 MPI_Allreduce(MPI_IN_PLACE
, r
, nr
, MPI_INT64_T
, MPI_SUM
,
731 ms
->mpi_comm_masters
);
733 /* this is thread-unsafe, but it will do for now: */
736 if (nr
> ms
->mpb
->libuf_alloc
)
738 ms
->mpb
->libuf_alloc
= nr
;
739 srenew(ms
->mpb
->libuf
, ms
->mpb
->libuf_alloc
);
741 MPI_Allreduce(r
, ms
->mpb
->libuf
, nr
, MPI_INT64_T
, MPI_SUM
,
742 ms
->mpi_comm_masters
);
743 for (i
= 0; i
< nr
; i
++)
745 r
[i
] = ms
->mpb
->libuf
[i
];
751 const char *opt2fn_master(const char *opt
, int nfile
, const t_filenm fnm
[],
754 return SIMMASTER(cr
) ? opt2fn(opt
, nfile
, fnm
) : nullptr;
757 void gmx_fatal_collective(int f_errno
, const char *file
, int line
,
758 MPI_Comm comm
, gmx_bool bMaster
,
759 const char *fmt
, ...)
765 /* Check if we are calling on all processes in MPI_COMM_WORLD */
766 MPI_Comm_compare(comm
, MPI_COMM_WORLD
, &result
);
767 /* Any result except MPI_UNEQUAL allows us to call MPI_Finalize */
768 bFinalize
= (result
!= MPI_UNEQUAL
);
770 GMX_UNUSED_VALUE(comm
);
775 gmx_fatal_mpi_va(f_errno
, file
, line
, bMaster
, bFinalize
, fmt
, ap
);