Remove gmx custom fixed int (e.g. gmx_int64_t) types
[alexxy/gromacs.git] / src / gromacs / gmxlib / network.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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,2018, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "network.h"
40
41 #include "config.h"
42
43 #include <cctype>
44 #include <cstdarg>
45 #include <cstdlib>
46 #include <cstring>
47
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"
56
57 /* The source code in this file should be thread-safe.
58       Please keep it that way. */
59
60 void gmx_fill_commrec_from_mpi(t_commrec *cr)
61 {
62 #if !GMX_MPI
63     gmx_call("gmx_fill_commrec_from_mpi");
64     GMX_UNUSED_VALUE(cr);
65 #else
66     if (!gmx_mpi_initialized())
67     {
68         gmx_comm("MPI has not been initialized properly");
69     }
70
71     cr->nnodes           = gmx_node_num();
72     cr->nodeid           = gmx_node_rank();
73     cr->sim_nodeid       = cr->nodeid;
74     cr->mpi_comm_mysim   = MPI_COMM_WORLD;
75     cr->mpi_comm_mygroup = MPI_COMM_WORLD;
76 #endif
77 }
78
79 t_commrec *init_commrec()
80 {
81     t_commrec    *cr;
82
83     snew(cr, 1);
84
85 #if GMX_LIB_MPI
86     gmx_fill_commrec_from_mpi(cr);
87 #else
88     cr->mpi_comm_mysim   = MPI_COMM_NULL;
89     cr->mpi_comm_mygroup = MPI_COMM_NULL;
90     cr->nnodes           = 1;
91     cr->sim_nodeid       = 0;
92     cr->nodeid           = cr->sim_nodeid;
93 #endif
94
95     // TODO cr->duty should not be initialized here
96     cr->duty = (DUTY_PP | DUTY_PME);
97
98 #if GMX_MPI && !MPI_IN_PLACE_EXISTS
99     /* initialize the MPI_IN_PLACE replacement buffers */
100     snew(cr->mpb, 1);
101     cr->mpb->ibuf        = NULL;
102     cr->mpb->libuf       = NULL;
103     cr->mpb->fbuf        = NULL;
104     cr->mpb->dbuf        = NULL;
105     cr->mpb->ibuf_alloc  = 0;
106     cr->mpb->libuf_alloc = 0;
107     cr->mpb->fbuf_alloc  = 0;
108     cr->mpb->dbuf_alloc  = 0;
109 #endif
110
111     return cr;
112 }
113
114 void done_mpi_in_place_buf(mpi_in_place_buf_t *buf)
115 {
116     if (nullptr != buf)
117     {
118         sfree(buf->ibuf);
119         sfree(buf->libuf);
120         sfree(buf->fbuf);
121         sfree(buf->dbuf);
122         sfree(buf);
123     }
124 }
125
126 void done_commrec(t_commrec *cr)
127 {
128     if (nullptr != cr->dd)
129     {
130         // TODO: implement
131         // done_domdec(cr->dd);
132     }
133     done_mpi_in_place_buf(cr->mpb);
134     sfree(cr);
135 }
136
137 t_commrec *reinitialize_commrec_for_this_thread(const t_commrec *cro)
138 {
139 #if GMX_THREAD_MPI
140     t_commrec *cr;
141
142     /* make a thread-specific commrec */
143     snew(cr, 1);
144     /* now copy the whole thing, so settings like the number of PME nodes
145        get propagated. */
146     *cr = *cro;
147
148     /* and we start setting our own thread-specific values for things */
149     gmx_fill_commrec_from_mpi(cr);
150
151     // TODO cr->duty should not be initialized here
152     cr->duty             = (DUTY_PP | DUTY_PME);
153
154     return cr;
155 #else
156     GMX_UNUSED_VALUE(cro);
157     return nullptr;
158 #endif
159 }
160
161 void gmx_setup_nodecomm(FILE gmx_unused *fplog, t_commrec *cr)
162 {
163     gmx_nodecomm_t *nc;
164
165     /* Many MPI implementations do not optimize MPI_Allreduce
166      * (and probably also other global communication calls)
167      * for multi-core nodes connected by a network.
168      * We can optimize such communication by using one MPI call
169      * within each node and one between the nodes.
170      * For MVAPICH2 and Intel MPI this reduces the time for
171      * the global_stat communication by 25%
172      * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
173      * B. Hess, November 2007
174      */
175
176     nc = &cr->nc;
177
178     nc->bUse = FALSE;
179 #if !GMX_THREAD_MPI
180 #if GMX_MPI
181     int n, rank;
182
183     // TODO PhysicalNodeCommunicator could be extended/used to handle
184     // the need for per-node per-group communicators.
185     MPI_Comm_size(cr->mpi_comm_mygroup, &n);
186     MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
187
188     int nodehash = gmx_physicalnode_id_hash();
189
190     if (debug)
191     {
192         fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
193     }
194
195
196     /* The intra-node communicator, split on node number */
197     MPI_Comm_split(cr->mpi_comm_mygroup, nodehash, rank, &nc->comm_intra);
198     MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
199     if (debug)
200     {
201         fprintf(debug, "In gmx_setup_nodecomm: node ID %d rank within node %d\n",
202                 rank, nc->rank_intra);
203     }
204     /* The inter-node communicator, split on rank_intra.
205      * We actually only need the one for rank=0,
206      * but it is easier to create them all.
207      */
208     MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
209     /* Check if this really created two step communication */
210     int ng, ni;
211
212     MPI_Comm_size(nc->comm_inter, &ng);
213     MPI_Comm_size(nc->comm_intra, &ni);
214     if (debug)
215     {
216         fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
217                 ng, ni);
218     }
219
220     if (getenv("GMX_NO_NODECOMM") == nullptr &&
221         ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
222     {
223         nc->bUse = TRUE;
224         if (fplog)
225         {
226             fprintf(fplog, "Using two step summing over %d groups of on average %.1f ranks\n\n",
227                     ng, (real)n/(real)ng);
228         }
229         if (nc->rank_intra > 0)
230         {
231             MPI_Comm_free(&nc->comm_inter);
232         }
233     }
234     else
235     {
236         /* One group or all processes in a separate group, use normal summing */
237         MPI_Comm_free(&nc->comm_inter);
238         MPI_Comm_free(&nc->comm_intra);
239         if (debug)
240         {
241             fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
242         }
243     }
244 #endif
245 #else
246     /* tMPI runs only on a single node so just use the nodeid */
247     nc->rank_intra = cr->nodeid;
248 #endif
249 }
250
251 void gmx_barrier(const t_commrec gmx_unused *cr)
252 {
253 #if !GMX_MPI
254     gmx_call("gmx_barrier");
255 #else
256     MPI_Barrier(cr->mpi_comm_mygroup);
257 #endif
258 }
259
260 void gmx_bcast(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
261 {
262 #if !GMX_MPI
263     gmx_call("gmx_bast");
264 #else
265     MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
266 #endif
267 }
268
269 void gmx_bcast_sim(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
270 {
271 #if !GMX_MPI
272     gmx_call("gmx_bast");
273 #else
274     MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
275 #endif
276 }
277
278 void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unused *cr)
279 {
280 #if !GMX_MPI
281     gmx_call("gmx_sumd");
282 #else
283 #if MPI_IN_PLACE_EXISTS
284     if (cr->nc.bUse)
285     {
286         if (cr->nc.rank_intra == 0)
287         {
288             /* Use two step summing. */
289             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
290                        cr->nc.comm_intra);
291             /* Sum the roots of the internal (intra) buffers. */
292             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
293                           cr->nc.comm_inter);
294         }
295         else
296         {
297             /* This is here because of the silly MPI specification
298                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
299             MPI_Reduce(r, nullptr, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
300         }
301         MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
302     }
303     else
304     {
305         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
306                       cr->mpi_comm_mygroup);
307     }
308 #else
309     int i;
310
311     if (nr > cr->mpb->dbuf_alloc)
312     {
313         cr->mpb->dbuf_alloc = nr;
314         srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
315     }
316     if (cr->nc.bUse)
317     {
318         /* Use two step summing */
319         MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
320         if (cr->nc.rank_intra == 0)
321         {
322             /* Sum with the buffers reversed */
323             MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
324                           cr->nc.comm_inter);
325         }
326         MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
327     }
328     else
329     {
330         MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
331                       cr->mpi_comm_mygroup);
332         for (i = 0; i < nr; i++)
333         {
334             r[i] = cr->mpb->dbuf[i];
335         }
336     }
337 #endif
338 #endif
339 }
340
341 void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unused *cr)
342 {
343 #if !GMX_MPI
344     gmx_call("gmx_sumf");
345 #else
346 #if MPI_IN_PLACE_EXISTS
347     if (cr->nc.bUse)
348     {
349         /* Use two step summing.  */
350         if (cr->nc.rank_intra == 0)
351         {
352             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
353                        cr->nc.comm_intra);
354             /* Sum the roots of the internal (intra) buffers */
355             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
356                           cr->nc.comm_inter);
357         }
358         else
359         {
360             /* This is here because of the silly MPI specification
361                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
362             MPI_Reduce(r, nullptr, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
363         }
364         MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
365     }
366     else
367     {
368         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
369     }
370 #else
371     int i;
372
373     if (nr > cr->mpb->fbuf_alloc)
374     {
375         cr->mpb->fbuf_alloc = nr;
376         srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
377     }
378     if (cr->nc.bUse)
379     {
380         /* Use two step summing */
381         MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
382         if (cr->nc.rank_intra == 0)
383         {
384             /* Sum with the buffers reversed */
385             MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
386                           cr->nc.comm_inter);
387         }
388         MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
389     }
390     else
391     {
392         MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
393                       cr->mpi_comm_mygroup);
394         for (i = 0; i < nr; i++)
395         {
396             r[i] = cr->mpb->fbuf[i];
397         }
398     }
399 #endif
400 #endif
401 }
402
403 void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused *cr)
404 {
405 #if !GMX_MPI
406     gmx_call("gmx_sumi");
407 #else
408 #if MPI_IN_PLACE_EXISTS
409     if (cr->nc.bUse)
410     {
411         /* Use two step summing */
412         if (cr->nc.rank_intra == 0)
413         {
414             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
415             /* Sum with the buffers reversed */
416             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
417         }
418         else
419         {
420             /* This is here because of the silly MPI specification
421                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
422             MPI_Reduce(r, nullptr, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
423         }
424         MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
425     }
426     else
427     {
428         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
429     }
430 #else
431     int i;
432
433     if (nr > cr->mpb->ibuf_alloc)
434     {
435         cr->mpb->ibuf_alloc = nr;
436         srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
437     }
438     if (cr->nc.bUse)
439     {
440         /* Use two step summing */
441         MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
442         if (cr->nc.rank_intra == 0)
443         {
444             /* Sum with the buffers reversed */
445             MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
446         }
447         MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
448     }
449     else
450     {
451         MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
452         for (i = 0; i < nr; i++)
453         {
454             r[i] = cr->mpb->ibuf[i];
455         }
456     }
457 #endif
458 #endif
459 }
460
461 void gmx_sumli(int gmx_unused nr, int64_t gmx_unused r[], const t_commrec gmx_unused *cr)
462 {
463 #if !GMX_MPI
464     gmx_call("gmx_sumli");
465 #else
466 #if MPI_IN_PLACE_EXISTS
467     if (cr->nc.bUse)
468     {
469         /* Use two step summing */
470         if (cr->nc.rank_intra == 0)
471         {
472             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, 0,
473                        cr->nc.comm_intra);
474             /* Sum with the buffers reversed */
475             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
476                           cr->nc.comm_inter);
477         }
478         else
479         {
480             /* This is here because of the silly MPI specification
481                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
482             MPI_Reduce(r, nullptr, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
483         }
484         MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
485     }
486     else
487     {
488         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
489     }
490 #else
491     int i;
492
493     if (nr > cr->mpb->libuf_alloc)
494     {
495         cr->mpb->libuf_alloc = nr;
496         srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
497     }
498     if (cr->nc.bUse)
499     {
500         /* Use two step summing */
501         MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
502                       cr->nc.comm_intra);
503         if (cr->nc.rank_intra == 0)
504         {
505             /* Sum with the buffers reversed */
506             MPI_Allreduce(cr->mpb->libuf, r, nr, MPI_INT64_T, MPI_SUM,
507                           cr->nc.comm_inter);
508         }
509         MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
510     }
511     else
512     {
513         MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
514                       cr->mpi_comm_mygroup);
515         for (i = 0; i < nr; i++)
516         {
517             r[i] = cr->mpb->libuf[i];
518         }
519     }
520 #endif
521 #endif
522 }
523
524
525
526 #if GMX_MPI
527 static void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
528 {
529 #if MPI_IN_PLACE_EXISTS
530     MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
531 #else
532     /* this function is only used in code that is not performance critical,
533        (during setup, when comm_rec is not the appropriate communication
534        structure), so this isn't as bad as it looks. */
535     double *buf;
536     int     i;
537
538     snew(buf, nr);
539     MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
540     for (i = 0; i < nr; i++)
541     {
542         r[i] = buf[i];
543     }
544     sfree(buf);
545 #endif
546 }
547 #endif
548
549 #if GMX_MPI
550 static void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
551 {
552 #if MPI_IN_PLACE_EXISTS
553     MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
554 #else
555     /* this function is only used in code that is not performance critical,
556        (during setup, when comm_rec is not the appropriate communication
557        structure), so this isn't as bad as it looks. */
558     float *buf;
559     int    i;
560
561     snew(buf, nr);
562     MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
563     for (i = 0; i < nr; i++)
564     {
565         r[i] = buf[i];
566     }
567     sfree(buf);
568 #endif
569 }
570 #endif
571
572 void gmx_sumd_sim(int gmx_unused nr, double gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
573 {
574 #if !GMX_MPI
575     gmx_call("gmx_sumd_sim");
576 #else
577     gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
578 #endif
579 }
580
581 void gmx_sumf_sim(int gmx_unused nr, float gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
582 {
583 #if !GMX_MPI
584     gmx_call("gmx_sumf_sim");
585 #else
586     gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
587 #endif
588 }
589
590 void gmx_sumi_sim(int gmx_unused nr, int gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
591 {
592 #if !GMX_MPI
593     gmx_call("gmx_sumi_sim");
594 #else
595 #if MPI_IN_PLACE_EXISTS
596     MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
597 #else
598     /* this is thread-unsafe, but it will do for now: */
599     int i;
600
601     if (nr > ms->mpb->ibuf_alloc)
602     {
603         ms->mpb->ibuf_alloc = nr;
604         srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
605     }
606     MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
607     for (i = 0; i < nr; i++)
608     {
609         r[i] = ms->mpb->ibuf[i];
610     }
611 #endif
612 #endif
613 }
614
615 void gmx_sumli_sim(int gmx_unused nr, int64_t gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
616 {
617 #if !GMX_MPI
618     gmx_call("gmx_sumli_sim");
619 #else
620 #if MPI_IN_PLACE_EXISTS
621     MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
622                   ms->mpi_comm_masters);
623 #else
624     /* this is thread-unsafe, but it will do for now: */
625     int i;
626
627     if (nr > ms->mpb->libuf_alloc)
628     {
629         ms->mpb->libuf_alloc = nr;
630         srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
631     }
632     MPI_Allreduce(r, ms->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
633                   ms->mpi_comm_masters);
634     for (i = 0; i < nr; i++)
635     {
636         r[i] = ms->mpb->libuf[i];
637     }
638 #endif
639 #endif
640 }
641
642 const char *opt2fn_master(const char *opt, int nfile, const t_filenm fnm[],
643                           t_commrec *cr)
644 {
645     return SIMMASTER(cr) ? opt2fn(opt, nfile, fnm) : nullptr;
646 }
647
648 void gmx_fatal_collective(int f_errno, const char *file, int line,
649                           MPI_Comm comm, gmx_bool bMaster,
650                           const char *fmt, ...)
651 {
652     va_list  ap;
653     gmx_bool bFinalize;
654 #if GMX_MPI
655     int      result;
656     /* Check if we are calling on all processes in MPI_COMM_WORLD */
657     MPI_Comm_compare(comm, MPI_COMM_WORLD, &result);
658     /* Any result except MPI_UNEQUAL allows us to call MPI_Finalize */
659     bFinalize = (result != MPI_UNEQUAL);
660 #else
661     GMX_UNUSED_VALUE(comm);
662     bFinalize = TRUE;
663 #endif
664
665     va_start(ap, fmt);
666     gmx_fatal_mpi_va(f_errno, file, line, bMaster, bFinalize, fmt, ap);
667     va_end(ap);
668 }