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