f6ddd9a1f0de240b7fb01db408f6145e40219a79
[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, 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 gmx_unused *cr)
61 {
62 #if !GMX_MPI
63     gmx_call("gmx_fill_commrec_from_mpi");
64 #else
65     if (!gmx_mpi_initialized())
66     {
67         gmx_comm("MPI has not been initialized properly");
68     }
69
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;
75
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   = NULL;
89     cr->mpi_comm_mygroup = 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 static void done_mpi_in_place_buf(mpi_in_place_buf_t *buf)
115 {
116     if (NULL != 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 (NULL != cr->dd)
129     {
130         // TODO: implement
131         // done_domdec(cr->dd);
132     }
133     if (NULL != cr->ms)
134     {
135         done_mpi_in_place_buf(cr->ms->mpb);
136         sfree(cr->ms);
137     }
138     done_mpi_in_place_buf(cr->mpb);
139     sfree(cr);
140 }
141
142 t_commrec *reinitialize_commrec_for_this_thread(const t_commrec gmx_unused *cro)
143 {
144 #if GMX_THREAD_MPI
145     t_commrec *cr;
146
147     /* make a thread-specific commrec */
148     snew(cr, 1);
149     /* now copy the whole thing, so settings like the number of PME nodes
150        get propagated. */
151     *cr = *cro;
152
153     /* and we start setting our own thread-specific values for things */
154     gmx_fill_commrec_from_mpi(cr);
155
156     // TODO cr->duty should not be initialized here
157     cr->duty             = (DUTY_PP | DUTY_PME);
158
159     return cr;
160 #else
161     return NULL;
162 #endif
163 }
164
165 void gmx_setup_nodecomm(FILE gmx_unused *fplog, t_commrec *cr)
166 {
167     gmx_nodecomm_t *nc;
168
169     /* Many MPI implementations do not optimize MPI_Allreduce
170      * (and probably also other global communication calls)
171      * for multi-core nodes connected by a network.
172      * We can optimize such communication by using one MPI call
173      * within each node and one between the nodes.
174      * For MVAPICH2 and Intel MPI this reduces the time for
175      * the global_stat communication by 25%
176      * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
177      * B. Hess, November 2007
178      */
179
180     nc = &cr->nc;
181
182     nc->bUse = FALSE;
183 #if !GMX_THREAD_MPI
184 #if GMX_MPI
185     int n, rank;
186
187     MPI_Comm_size(cr->mpi_comm_mygroup, &n);
188     MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
189
190     int nodehash = gmx_physicalnode_id_hash();
191
192     if (debug)
193     {
194         fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
195     }
196
197
198     /* The intra-node communicator, split on node number */
199     MPI_Comm_split(cr->mpi_comm_mygroup, nodehash, rank, &nc->comm_intra);
200     MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
201     if (debug)
202     {
203         fprintf(debug, "In gmx_setup_nodecomm: node ID %d rank within node %d\n",
204                 rank, nc->rank_intra);
205     }
206     /* The inter-node communicator, split on rank_intra.
207      * We actually only need the one for rank=0,
208      * but it is easier to create them all.
209      */
210     MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
211     /* Check if this really created two step communication */
212     int ng, ni;
213
214     MPI_Comm_size(nc->comm_inter, &ng);
215     MPI_Comm_size(nc->comm_intra, &ni);
216     if (debug)
217     {
218         fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
219                 ng, ni);
220     }
221
222     if (getenv("GMX_NO_NODECOMM") == NULL &&
223         ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
224     {
225         nc->bUse = TRUE;
226         if (fplog)
227         {
228             fprintf(fplog, "Using two step summing over %d groups of on average %.1f ranks\n\n",
229                     ng, (real)n/(real)ng);
230         }
231         if (nc->rank_intra > 0)
232         {
233             MPI_Comm_free(&nc->comm_inter);
234         }
235     }
236     else
237     {
238         /* One group or all processes in a separate group, use normal summing */
239         MPI_Comm_free(&nc->comm_inter);
240         MPI_Comm_free(&nc->comm_intra);
241         if (debug)
242         {
243             fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
244         }
245     }
246 #endif
247 #else
248     /* tMPI runs only on a single node so just use the nodeid */
249     nc->rank_intra = cr->nodeid;
250 #endif
251 }
252
253 void gmx_init_intranode_counters(t_commrec *cr)
254 {
255     /* counters for PP+PME and PP-only processes on my physical node */
256     int nrank_intranode, rank_intranode;
257     int nrank_pp_intranode, rank_pp_intranode;
258     /* thread-MPI is not initialized when not running in parallel */
259 #if GMX_MPI && !GMX_THREAD_MPI
260     int nrank_world, rank_world;
261     int i, myhash, *hash, *hash_s, *hash_pp, *hash_pp_s;
262
263     MPI_Comm_size(MPI_COMM_WORLD, &nrank_world);
264     MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
265
266     /* Get a (hopefully unique) hash that identifies our physical node */
267     myhash = gmx_physicalnode_id_hash();
268
269     /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
270     snew(hash,   nrank_world);
271     snew(hash_s, nrank_world);
272     snew(hash_pp,   nrank_world);
273     snew(hash_pp_s, nrank_world);
274
275     hash_s[rank_world]    = myhash;
276     hash_pp_s[rank_world] = (cr->duty & DUTY_PP) ? myhash : -1;
277
278     MPI_Allreduce(hash_s,    hash,    nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
279     MPI_Allreduce(hash_pp_s, hash_pp, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
280
281     nrank_intranode    = 0;
282     rank_intranode     = 0;
283     nrank_pp_intranode = 0;
284     rank_pp_intranode  = 0;
285     for (i = 0; i < nrank_world; i++)
286     {
287         if (hash[i] == myhash)
288         {
289             nrank_intranode++;
290             if (i < rank_world)
291             {
292                 rank_intranode++;
293             }
294         }
295         if (hash_pp[i] == myhash)
296         {
297             nrank_pp_intranode++;
298             if ((cr->duty & DUTY_PP) && i < rank_world)
299             {
300                 rank_pp_intranode++;
301             }
302         }
303     }
304     sfree(hash);
305     sfree(hash_s);
306     sfree(hash_pp);
307     sfree(hash_pp_s);
308 #else
309     /* Serial or thread-MPI code: we run within a single physical node */
310     nrank_intranode    = cr->nnodes;
311     rank_intranode     = cr->sim_nodeid;
312     nrank_pp_intranode = cr->nnodes - cr->npmenodes;
313     rank_pp_intranode  = cr->nodeid;
314 #endif
315
316     if (debug)
317     {
318         char sbuf[STRLEN];
319         if ((cr->duty & DUTY_PP) && (cr->duty & DUTY_PME))
320         {
321             sprintf(sbuf, "PP+PME");
322         }
323         else
324         {
325             sprintf(sbuf, "%s", (cr->duty & DUTY_PP) ? "PP" : "PME");
326         }
327         fprintf(debug, "On %3s rank %d: nrank_intranode=%d, rank_intranode=%d, "
328                 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
329                 sbuf, cr->sim_nodeid,
330                 nrank_intranode, rank_intranode,
331                 nrank_pp_intranode, rank_pp_intranode);
332     }
333
334     cr->nrank_intranode    = nrank_intranode;
335     cr->rank_intranode     = rank_intranode;
336     cr->nrank_pp_intranode = nrank_pp_intranode;
337     cr->rank_pp_intranode  = rank_pp_intranode;
338 }
339
340
341 void gmx_barrier(const t_commrec gmx_unused *cr)
342 {
343 #if !GMX_MPI
344     gmx_call("gmx_barrier");
345 #else
346     MPI_Barrier(cr->mpi_comm_mygroup);
347 #endif
348 }
349
350 void gmx_bcast(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
351 {
352 #if !GMX_MPI
353     gmx_call("gmx_bast");
354 #else
355     MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
356 #endif
357 }
358
359 void gmx_bcast_sim(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
360 {
361 #if !GMX_MPI
362     gmx_call("gmx_bast");
363 #else
364     MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
365 #endif
366 }
367
368 void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unused *cr)
369 {
370 #if !GMX_MPI
371     gmx_call("gmx_sumd");
372 #else
373 #if MPI_IN_PLACE_EXISTS
374     if (cr->nc.bUse)
375     {
376         if (cr->nc.rank_intra == 0)
377         {
378             /* Use two step summing. */
379             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
380                        cr->nc.comm_intra);
381             /* Sum the roots of the internal (intra) buffers. */
382             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
383                           cr->nc.comm_inter);
384         }
385         else
386         {
387             /* This is here because of the silly MPI specification
388                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
389             MPI_Reduce(r, NULL, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
390         }
391         MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
392     }
393     else
394     {
395         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
396                       cr->mpi_comm_mygroup);
397     }
398 #else
399     int i;
400
401     if (nr > cr->mpb->dbuf_alloc)
402     {
403         cr->mpb->dbuf_alloc = nr;
404         srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
405     }
406     if (cr->nc.bUse)
407     {
408         /* Use two step summing */
409         MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
410         if (cr->nc.rank_intra == 0)
411         {
412             /* Sum with the buffers reversed */
413             MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
414                           cr->nc.comm_inter);
415         }
416         MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
417     }
418     else
419     {
420         MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
421                       cr->mpi_comm_mygroup);
422         for (i = 0; i < nr; i++)
423         {
424             r[i] = cr->mpb->dbuf[i];
425         }
426     }
427 #endif
428 #endif
429 }
430
431 void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unused *cr)
432 {
433 #if !GMX_MPI
434     gmx_call("gmx_sumf");
435 #else
436 #if MPI_IN_PLACE_EXISTS
437     if (cr->nc.bUse)
438     {
439         /* Use two step summing.  */
440         if (cr->nc.rank_intra == 0)
441         {
442             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
443                        cr->nc.comm_intra);
444             /* Sum the roots of the internal (intra) buffers */
445             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
446                           cr->nc.comm_inter);
447         }
448         else
449         {
450             /* This is here because of the silly MPI specification
451                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
452             MPI_Reduce(r, NULL, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
453         }
454         MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
455     }
456     else
457     {
458         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
459     }
460 #else
461     int i;
462
463     if (nr > cr->mpb->fbuf_alloc)
464     {
465         cr->mpb->fbuf_alloc = nr;
466         srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
467     }
468     if (cr->nc.bUse)
469     {
470         /* Use two step summing */
471         MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
472         if (cr->nc.rank_intra == 0)
473         {
474             /* Sum with the buffers reversed */
475             MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
476                           cr->nc.comm_inter);
477         }
478         MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
479     }
480     else
481     {
482         MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
483                       cr->mpi_comm_mygroup);
484         for (i = 0; i < nr; i++)
485         {
486             r[i] = cr->mpb->fbuf[i];
487         }
488     }
489 #endif
490 #endif
491 }
492
493 void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused *cr)
494 {
495 #if !GMX_MPI
496     gmx_call("gmx_sumi");
497 #else
498 #if MPI_IN_PLACE_EXISTS
499     if (cr->nc.bUse)
500     {
501         /* Use two step summing */
502         if (cr->nc.rank_intra == 0)
503         {
504             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
505             /* Sum with the buffers reversed */
506             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
507         }
508         else
509         {
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_INT, MPI_SUM, 0, cr->nc.comm_intra);
513         }
514         MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
515     }
516     else
517     {
518         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
519     }
520 #else
521     int i;
522
523     if (nr > cr->mpb->ibuf_alloc)
524     {
525         cr->mpb->ibuf_alloc = nr;
526         srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
527     }
528     if (cr->nc.bUse)
529     {
530         /* Use two step summing */
531         MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
532         if (cr->nc.rank_intra == 0)
533         {
534             /* Sum with the buffers reversed */
535             MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
536         }
537         MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
538     }
539     else
540     {
541         MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
542         for (i = 0; i < nr; i++)
543         {
544             r[i] = cr->mpb->ibuf[i];
545         }
546     }
547 #endif
548 #endif
549 }
550
551 void gmx_sumli(int gmx_unused nr, gmx_int64_t gmx_unused r[], const t_commrec gmx_unused *cr)
552 {
553 #if !GMX_MPI
554     gmx_call("gmx_sumli");
555 #else
556 #if MPI_IN_PLACE_EXISTS
557     if (cr->nc.bUse)
558     {
559         /* Use two step summing */
560         if (cr->nc.rank_intra == 0)
561         {
562             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, 0,
563                        cr->nc.comm_intra);
564             /* Sum with the buffers reversed */
565             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
566                           cr->nc.comm_inter);
567         }
568         else
569         {
570             /* This is here because of the silly MPI specification
571                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
572             MPI_Reduce(r, NULL, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
573         }
574         MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
575     }
576     else
577     {
578         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
579     }
580 #else
581     int i;
582
583     if (nr > cr->mpb->libuf_alloc)
584     {
585         cr->mpb->libuf_alloc = nr;
586         srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
587     }
588     if (cr->nc.bUse)
589     {
590         /* Use two step summing */
591         MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
592                       cr->nc.comm_intra);
593         if (cr->nc.rank_intra == 0)
594         {
595             /* Sum with the buffers reversed */
596             MPI_Allreduce(cr->mpb->libuf, r, nr, MPI_INT64_T, MPI_SUM,
597                           cr->nc.comm_inter);
598         }
599         MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
600     }
601     else
602     {
603         MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
604                       cr->mpi_comm_mygroup);
605         for (i = 0; i < nr; i++)
606         {
607             r[i] = cr->mpb->libuf[i];
608         }
609     }
610 #endif
611 #endif
612 }
613
614
615
616 #if GMX_MPI
617 static void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
618 {
619 #if MPI_IN_PLACE_EXISTS
620     MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
621 #else
622     /* this function is only used in code that is not performance critical,
623        (during setup, when comm_rec is not the appropriate communication
624        structure), so this isn't as bad as it looks. */
625     double *buf;
626     int     i;
627
628     snew(buf, nr);
629     MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
630     for (i = 0; i < nr; i++)
631     {
632         r[i] = buf[i];
633     }
634     sfree(buf);
635 #endif
636 }
637 #endif
638
639 #if GMX_MPI
640 static void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
641 {
642 #if MPI_IN_PLACE_EXISTS
643     MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
644 #else
645     /* this function is only used in code that is not performance critical,
646        (during setup, when comm_rec is not the appropriate communication
647        structure), so this isn't as bad as it looks. */
648     float *buf;
649     int    i;
650
651     snew(buf, nr);
652     MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
653     for (i = 0; i < nr; i++)
654     {
655         r[i] = buf[i];
656     }
657     sfree(buf);
658 #endif
659 }
660 #endif
661
662 void gmx_sumd_sim(int gmx_unused nr, double gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
663 {
664 #if !GMX_MPI
665     gmx_call("gmx_sumd_sim");
666 #else
667     gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
668 #endif
669 }
670
671 void gmx_sumf_sim(int gmx_unused nr, float gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
672 {
673 #if !GMX_MPI
674     gmx_call("gmx_sumf_sim");
675 #else
676     gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
677 #endif
678 }
679
680 void gmx_sumi_sim(int gmx_unused nr, int gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
681 {
682 #if !GMX_MPI
683     gmx_call("gmx_sumi_sim");
684 #else
685 #if MPI_IN_PLACE_EXISTS
686     MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
687 #else
688     /* this is thread-unsafe, but it will do for now: */
689     int i;
690
691     if (nr > ms->mpb->ibuf_alloc)
692     {
693         ms->mpb->ibuf_alloc = nr;
694         srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
695     }
696     MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
697     for (i = 0; i < nr; i++)
698     {
699         r[i] = ms->mpb->ibuf[i];
700     }
701 #endif
702 #endif
703 }
704
705 void gmx_sumli_sim(int gmx_unused nr, gmx_int64_t gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
706 {
707 #if !GMX_MPI
708     gmx_call("gmx_sumli_sim");
709 #else
710 #if MPI_IN_PLACE_EXISTS
711     MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
712                   ms->mpi_comm_masters);
713 #else
714     /* this is thread-unsafe, but it will do for now: */
715     int i;
716
717     if (nr > ms->mpb->libuf_alloc)
718     {
719         ms->mpb->libuf_alloc = nr;
720         srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
721     }
722     MPI_Allreduce(r, ms->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
723                   ms->mpi_comm_masters);
724     for (i = 0; i < nr; i++)
725     {
726         r[i] = ms->mpb->libuf[i];
727     }
728 #endif
729 #endif
730 }
731
732 const char *opt2fn_master(const char *opt, int nfile, const t_filenm fnm[],
733                           t_commrec *cr)
734 {
735     return SIMMASTER(cr) ? opt2fn(opt, nfile, fnm) : NULL;
736 }
737
738 void gmx_fatal_collective(int f_errno, const char *file, int line,
739                           MPI_Comm comm, gmx_bool bMaster,
740                           const char *fmt, ...)
741 {
742     va_list  ap;
743     gmx_bool bFinalize;
744 #if GMX_MPI
745     int      result;
746     /* Check if we are calling on all processes in MPI_COMM_WORLD */
747     MPI_Comm_compare(comm, MPI_COMM_WORLD, &result);
748     /* Any result except MPI_UNEQUAL allows us to call MPI_Finalize */
749     bFinalize = (result != MPI_UNEQUAL);
750 #else
751     GMX_UNUSED_VALUE(comm);
752     bFinalize = TRUE;
753 #endif
754
755     va_start(ap, fmt);
756     gmx_fatal_mpi_va(f_errno, file, line, bMaster, bFinalize, fmt, ap);
757     va_end(ap);
758 }