More config.h macros to 0/1
[gromacs/AngularHB.git] / src / gromacs / gmxlib / network.cpp
blob485bfe8154dad330868cdda3db2b517423b3e50f
1 /*
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, 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.
37 #include "gmxpre.h"
39 #include "network.h"
41 #include "config.h"
43 #include <cctype>
44 #include <cstdarg>
45 #include <cstdlib>
46 #include <cstring>
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)
62 #ifndef GMX_MPI
63 gmx_call("gmx_fill_commrec_from_mpi");
64 #else
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 cr->sim_nodeid = cr->nodeid;
73 cr->mpi_comm_mysim = MPI_COMM_WORLD;
74 cr->mpi_comm_mygroup = MPI_COMM_WORLD;
76 #endif
79 t_commrec *init_commrec()
81 t_commrec *cr;
83 snew(cr, 1);
85 #ifdef GMX_LIB_MPI
86 gmx_fill_commrec_from_mpi(cr);
87 #else
88 cr->mpi_comm_mysim = NULL;
89 cr->mpi_comm_mygroup = NULL;
90 cr->nnodes = 1;
91 cr->sim_nodeid = 0;
92 cr->nodeid = cr->sim_nodeid;
93 #endif
95 // TODO cr->duty should not be initialized here
96 cr->duty = (DUTY_PP | DUTY_PME);
98 #if defined 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
111 return cr;
114 t_commrec *reinitialize_commrec_for_this_thread(const t_commrec gmx_unused *cro)
116 #ifdef GMX_THREAD_MPI
117 t_commrec *cr;
119 /* make a thread-specific commrec */
120 snew(cr, 1);
121 /* now copy the whole thing, so settings like the number of PME nodes
122 get propagated. */
123 *cr = *cro;
125 /* and we start setting our own thread-specific values for things */
126 gmx_fill_commrec_from_mpi(cr);
128 // TODO cr->duty should not be initialized here
129 cr->duty = (DUTY_PP | DUTY_PME);
131 return cr;
132 #else
133 return NULL;
134 #endif
137 void gmx_setup_nodecomm(FILE gmx_unused *fplog, t_commrec *cr)
139 gmx_nodecomm_t *nc;
141 /* Many MPI implementations do not optimize MPI_Allreduce
142 * (and probably also other global communication calls)
143 * for multi-core nodes connected by a network.
144 * We can optimize such communication by using one MPI call
145 * within each node and one between the nodes.
146 * For MVAPICH2 and Intel MPI this reduces the time for
147 * the global_stat communication by 25%
148 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
149 * B. Hess, November 2007
152 nc = &cr->nc;
154 nc->bUse = FALSE;
155 #ifndef GMX_THREAD_MPI
156 #ifdef GMX_MPI
157 int n, rank;
159 MPI_Comm_size(cr->mpi_comm_mygroup, &n);
160 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
162 int nodehash = gmx_physicalnode_id_hash();
164 if (debug)
166 fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
170 /* The intra-node communicator, split on node number */
171 MPI_Comm_split(cr->mpi_comm_mygroup, nodehash, rank, &nc->comm_intra);
172 MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
173 if (debug)
175 fprintf(debug, "In gmx_setup_nodecomm: node ID %d rank within node %d\n",
176 rank, nc->rank_intra);
178 /* The inter-node communicator, split on rank_intra.
179 * We actually only need the one for rank=0,
180 * but it is easier to create them all.
182 MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
183 /* Check if this really created two step communication */
184 int ng, ni;
186 MPI_Comm_size(nc->comm_inter, &ng);
187 MPI_Comm_size(nc->comm_intra, &ni);
188 if (debug)
190 fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
191 ng, ni);
194 if (getenv("GMX_NO_NODECOMM") == NULL &&
195 ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
197 nc->bUse = TRUE;
198 if (fplog)
200 fprintf(fplog, "Using two step summing over %d groups of on average %.1f ranks\n\n",
201 ng, (real)n/(real)ng);
203 if (nc->rank_intra > 0)
205 MPI_Comm_free(&nc->comm_inter);
208 else
210 /* One group or all processes in a separate group, use normal summing */
211 MPI_Comm_free(&nc->comm_inter);
212 MPI_Comm_free(&nc->comm_intra);
213 if (debug)
215 fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
218 #endif
219 #else
220 /* tMPI runs only on a single node so just use the nodeid */
221 nc->rank_intra = cr->nodeid;
222 #endif
225 void gmx_init_intranode_counters(t_commrec *cr)
227 /* counters for PP+PME and PP-only processes on my physical node */
228 int nrank_intranode, rank_intranode;
229 int nrank_pp_intranode, rank_pp_intranode;
230 /* thread-MPI is not initialized when not running in parallel */
231 #if defined GMX_MPI && !defined GMX_THREAD_MPI
232 int nrank_world, rank_world;
233 int i, myhash, *hash, *hash_s, *hash_pp, *hash_pp_s;
235 MPI_Comm_size(MPI_COMM_WORLD, &nrank_world);
236 MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
238 /* Get a (hopefully unique) hash that identifies our physical node */
239 myhash = gmx_physicalnode_id_hash();
241 /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
242 snew(hash, nrank_world);
243 snew(hash_s, nrank_world);
244 snew(hash_pp, nrank_world);
245 snew(hash_pp_s, nrank_world);
247 hash_s[rank_world] = myhash;
248 hash_pp_s[rank_world] = (cr->duty & DUTY_PP) ? myhash : -1;
250 MPI_Allreduce(hash_s, hash, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
251 MPI_Allreduce(hash_pp_s, hash_pp, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
253 nrank_intranode = 0;
254 rank_intranode = 0;
255 nrank_pp_intranode = 0;
256 rank_pp_intranode = 0;
257 for (i = 0; i < nrank_world; i++)
259 if (hash[i] == myhash)
261 nrank_intranode++;
262 if (i < rank_world)
264 rank_intranode++;
267 if (hash_pp[i] == myhash)
269 nrank_pp_intranode++;
270 if ((cr->duty & DUTY_PP) && i < rank_world)
272 rank_pp_intranode++;
276 sfree(hash);
277 sfree(hash_s);
278 sfree(hash_pp);
279 sfree(hash_pp_s);
280 #else
281 /* Serial or thread-MPI code: we run within a single physical node */
282 nrank_intranode = cr->nnodes;
283 rank_intranode = cr->sim_nodeid;
284 nrank_pp_intranode = cr->nnodes - cr->npmenodes;
285 rank_pp_intranode = cr->nodeid;
286 #endif
288 if (debug)
290 char sbuf[STRLEN];
291 if (cr->duty & DUTY_PP && cr->duty & DUTY_PME)
293 sprintf(sbuf, "PP+PME");
295 else
297 sprintf(sbuf, "%s", cr->duty & DUTY_PP ? "PP" : "PME");
299 fprintf(debug, "On %3s rank %d: nrank_intranode=%d, rank_intranode=%d, "
300 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
301 sbuf, cr->sim_nodeid,
302 nrank_intranode, rank_intranode,
303 nrank_pp_intranode, rank_pp_intranode);
306 cr->nrank_intranode = nrank_intranode;
307 cr->rank_intranode = rank_intranode;
308 cr->nrank_pp_intranode = nrank_pp_intranode;
309 cr->rank_pp_intranode = rank_pp_intranode;
313 void gmx_barrier(const t_commrec gmx_unused *cr)
315 #ifndef GMX_MPI
316 gmx_call("gmx_barrier");
317 #else
318 MPI_Barrier(cr->mpi_comm_mygroup);
319 #endif
322 void gmx_bcast(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
324 #ifndef GMX_MPI
325 gmx_call("gmx_bast");
326 #else
327 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
328 #endif
331 void gmx_bcast_sim(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
333 #ifndef GMX_MPI
334 gmx_call("gmx_bast");
335 #else
336 MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
337 #endif
340 void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unused *cr)
342 #ifndef GMX_MPI
343 gmx_call("gmx_sumd");
344 #else
345 #if MPI_IN_PLACE_EXISTS
346 if (cr->nc.bUse)
348 if (cr->nc.rank_intra == 0)
350 /* Use two step summing. */
351 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
352 cr->nc.comm_intra);
353 /* Sum the roots of the internal (intra) buffers. */
354 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
355 cr->nc.comm_inter);
357 else
359 /* This is here because of the silly MPI specification
360 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
361 MPI_Reduce(r, NULL, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
363 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
365 else
367 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
368 cr->mpi_comm_mygroup);
370 #else
371 int i;
373 if (nr > cr->mpb->dbuf_alloc)
375 cr->mpb->dbuf_alloc = nr;
376 srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
378 if (cr->nc.bUse)
380 /* Use two step summing */
381 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
382 if (cr->nc.rank_intra == 0)
384 /* Sum with the buffers reversed */
385 MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
386 cr->nc.comm_inter);
388 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
390 else
392 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
393 cr->mpi_comm_mygroup);
394 for (i = 0; i < nr; i++)
396 r[i] = cr->mpb->dbuf[i];
399 #endif
400 #endif
403 void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unused *cr)
405 #ifndef GMX_MPI
406 gmx_call("gmx_sumf");
407 #else
408 #if MPI_IN_PLACE_EXISTS
409 if (cr->nc.bUse)
411 /* Use two step summing. */
412 if (cr->nc.rank_intra == 0)
414 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
415 cr->nc.comm_intra);
416 /* Sum the roots of the internal (intra) buffers */
417 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
418 cr->nc.comm_inter);
420 else
422 /* This is here because of the silly MPI specification
423 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
424 MPI_Reduce(r, NULL, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
426 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
428 else
430 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
432 #else
433 int i;
435 if (nr > cr->mpb->fbuf_alloc)
437 cr->mpb->fbuf_alloc = nr;
438 srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
440 if (cr->nc.bUse)
442 /* Use two step summing */
443 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
444 if (cr->nc.rank_intra == 0)
446 /* Sum with the buffers reversed */
447 MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
448 cr->nc.comm_inter);
450 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
452 else
454 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
455 cr->mpi_comm_mygroup);
456 for (i = 0; i < nr; i++)
458 r[i] = cr->mpb->fbuf[i];
461 #endif
462 #endif
465 void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused *cr)
467 #ifndef GMX_MPI
468 gmx_call("gmx_sumi");
469 #else
470 #if MPI_IN_PLACE_EXISTS
471 if (cr->nc.bUse)
473 /* Use two step summing */
474 if (cr->nc.rank_intra == 0)
476 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
477 /* Sum with the buffers reversed */
478 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
480 else
482 /* This is here because of the silly MPI specification
483 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
484 MPI_Reduce(r, NULL, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
486 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
488 else
490 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
492 #else
493 int i;
495 if (nr > cr->mpb->ibuf_alloc)
497 cr->mpb->ibuf_alloc = nr;
498 srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
500 if (cr->nc.bUse)
502 /* Use two step summing */
503 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
504 if (cr->nc.rank_intra == 0)
506 /* Sum with the buffers reversed */
507 MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
509 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
511 else
513 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
514 for (i = 0; i < nr; i++)
516 r[i] = cr->mpb->ibuf[i];
519 #endif
520 #endif
523 void gmx_sumli(int gmx_unused nr, gmx_int64_t gmx_unused r[], const t_commrec gmx_unused *cr)
525 #ifndef GMX_MPI
526 gmx_call("gmx_sumli");
527 #else
528 #if MPI_IN_PLACE_EXISTS
529 if (cr->nc.bUse)
531 /* Use two step summing */
532 if (cr->nc.rank_intra == 0)
534 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, 0,
535 cr->nc.comm_intra);
536 /* Sum with the buffers reversed */
537 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
538 cr->nc.comm_inter);
540 else
542 /* This is here because of the silly MPI specification
543 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
544 MPI_Reduce(r, NULL, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
546 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
548 else
550 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
552 #else
553 int i;
555 if (nr > cr->mpb->libuf_alloc)
557 cr->mpb->libuf_alloc = nr;
558 srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
560 if (cr->nc.bUse)
562 /* Use two step summing */
563 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
564 cr->nc.comm_intra);
565 if (cr->nc.rank_intra == 0)
567 /* Sum with the buffers reversed */
568 MPI_Allreduce(cr->mpb->libuf, r, nr, MPI_INT64_T, MPI_SUM,
569 cr->nc.comm_inter);
571 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
573 else
575 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
576 cr->mpi_comm_mygroup);
577 for (i = 0; i < nr; i++)
579 r[i] = cr->mpb->libuf[i];
582 #endif
583 #endif
588 #ifdef GMX_MPI
589 static void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
591 #if MPI_IN_PLACE_EXISTS
592 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
593 #else
594 /* this function is only used in code that is not performance critical,
595 (during setup, when comm_rec is not the appropriate communication
596 structure), so this isn't as bad as it looks. */
597 double *buf;
598 int i;
600 snew(buf, nr);
601 MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
602 for (i = 0; i < nr; i++)
604 r[i] = buf[i];
606 sfree(buf);
607 #endif
609 #endif
611 #ifdef GMX_MPI
612 static void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
614 #if MPI_IN_PLACE_EXISTS
615 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
616 #else
617 /* this function is only used in code that is not performance critical,
618 (during setup, when comm_rec is not the appropriate communication
619 structure), so this isn't as bad as it looks. */
620 float *buf;
621 int i;
623 snew(buf, nr);
624 MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
625 for (i = 0; i < nr; i++)
627 r[i] = buf[i];
629 sfree(buf);
630 #endif
632 #endif
634 void gmx_sumd_sim(int gmx_unused nr, double gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
636 #ifndef GMX_MPI
637 gmx_call("gmx_sumd_sim");
638 #else
639 gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
640 #endif
643 void gmx_sumf_sim(int gmx_unused nr, float gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
645 #ifndef GMX_MPI
646 gmx_call("gmx_sumf_sim");
647 #else
648 gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
649 #endif
652 void gmx_sumi_sim(int gmx_unused nr, int gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
654 #ifndef GMX_MPI
655 gmx_call("gmx_sumi_sim");
656 #else
657 #if MPI_IN_PLACE_EXISTS
658 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
659 #else
660 /* this is thread-unsafe, but it will do for now: */
661 int i;
663 if (nr > ms->mpb->ibuf_alloc)
665 ms->mpb->ibuf_alloc = nr;
666 srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
668 MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
669 for (i = 0; i < nr; i++)
671 r[i] = ms->mpb->ibuf[i];
673 #endif
674 #endif
677 void gmx_sumli_sim(int gmx_unused nr, gmx_int64_t gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
679 #ifndef GMX_MPI
680 gmx_call("gmx_sumli_sim");
681 #else
682 #if MPI_IN_PLACE_EXISTS
683 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
684 ms->mpi_comm_masters);
685 #else
686 /* this is thread-unsafe, but it will do for now: */
687 int i;
689 if (nr > ms->mpb->libuf_alloc)
691 ms->mpb->libuf_alloc = nr;
692 srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
694 MPI_Allreduce(r, ms->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
695 ms->mpi_comm_masters);
696 for (i = 0; i < nr; i++)
698 r[i] = ms->mpb->libuf[i];
700 #endif
701 #endif
704 const char *opt2fn_master(const char *opt, int nfile, const t_filenm fnm[],
705 t_commrec *cr)
707 return SIMMASTER(cr) ? opt2fn(opt, nfile, fnm) : NULL;
710 void gmx_fatal_collective(int f_errno, const char *file, int line,
711 MPI_Comm comm, gmx_bool bMaster,
712 const char *fmt, ...)
714 va_list ap;
715 gmx_bool bFinalize;
716 #ifdef GMX_MPI
717 int result;
718 /* Check if we are calling on all processes in MPI_COMM_WORLD */
719 MPI_Comm_compare(comm, MPI_COMM_WORLD, &result);
720 /* Any result except MPI_UNEQUAL allows us to call MPI_Finalize */
721 bFinalize = (result != MPI_UNEQUAL);
722 #else
723 bFinalize = TRUE;
724 #endif
726 va_start(ap, fmt);
727 gmx_fatal_mpi_va(f_errno, file, line, bMaster, bFinalize, fmt, ap);
728 va_end(ap);