More consolidation of MPI initialization
[alexxy/gromacs.git] / src / gromacs / gmxlib / network.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * GROningen Mixture of Alchemy and Childrens' Stories
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <string.h>
40 #include "gmx_fatal.h"
41 #include "main.h"
42 #include "smalloc.h"
43 #include "network.h"
44 #include "copyrite.h"
45 #include "statutil.h"
46 #include <ctype.h>
47 #include "macros.h"
48
49 #include "gromacs/utility/gmxmpi.h"
50
51
52 /* The source code in this file should be thread-safe.
53       Please keep it that way. */
54
55 gmx_bool gmx_mpi_initialized(void)
56 {
57     int n;
58 #ifndef GMX_MPI
59     return 0;
60 #else
61     MPI_Initialized(&n);
62
63     return n;
64 #endif
65 }
66
67 void gmx_fill_commrec_from_mpi(t_commrec *cr)
68 {
69 #ifndef GMX_MPI
70     gmx_call("gmx_fill_commrec_from_mpi");
71 #else
72     cr->nnodes           = gmx_node_num();
73     cr->nodeid           = gmx_node_rank();
74     cr->sim_nodeid       = cr->nodeid;
75     cr->mpi_comm_mysim   = MPI_COMM_WORLD;
76     cr->mpi_comm_mygroup = MPI_COMM_WORLD;
77
78 #endif
79 }
80
81 t_commrec *init_commrec()
82 {
83     t_commrec    *cr;
84
85     snew(cr, 1);
86
87 #ifdef GMX_LIB_MPI
88     if (!gmx_mpi_initialized())
89     {
90         gmx_comm("MPI has not been initialized properly");
91     }
92
93     gmx_fill_commrec_from_mpi(cr);
94 #else
95     /* These should never be accessed */
96     cr->mpi_comm_mysim   = NULL;
97     cr->mpi_comm_mygroup = NULL;
98     cr->nnodes           = 1;
99     cr->sim_nodeid       = 0;
100     cr->nodeid           = cr->sim_nodeid;
101 #endif
102
103     // TODO this should be initialized elsewhere
104     cr->duty = (DUTY_PP | DUTY_PME);
105
106 #if defined GMX_LIB_MPI && !defined MPI_IN_PLACE_EXISTS
107     /* initialize the MPI_IN_PLACE replacement buffers */
108     snew(cr->mpb, 1);
109     cr->mpb->ibuf        = NULL;
110     cr->mpb->libuf       = NULL;
111     cr->mpb->fbuf        = NULL;
112     cr->mpb->dbuf        = NULL;
113     cr->mpb->ibuf_alloc  = 0;
114     cr->mpb->libuf_alloc = 0;
115     cr->mpb->fbuf_alloc  = 0;
116     cr->mpb->dbuf_alloc  = 0;
117 #endif
118
119     return cr;
120 }
121
122 t_commrec *init_par_threads(const t_commrec *cro)
123 {
124 #ifdef GMX_THREAD_MPI
125     int        initialized;
126     t_commrec *cr;
127
128     /* make a thread-specific commrec */
129     snew(cr, 1);
130     /* now copy the whole thing, so settings like the number of PME nodes
131        get propagated. */
132     *cr = *cro;
133
134     /* and we start setting our own thread-specific values for things */
135     MPI_Initialized(&initialized);
136     if (!initialized)
137     {
138         gmx_comm("Initializing threads without comm");
139     }
140
141     /* No need to do MPI_Init() for thread-MPI. */
142     gmx_fill_commrec_from_mpi(cr);
143
144     // TODO cr->duty should not be initialized here
145     cr->duty             = (DUTY_PP | DUTY_PME);
146
147     return cr;
148 #else
149     return NULL;
150 #endif
151 }
152
153 int  gmx_node_num(void)
154 {
155 #ifndef GMX_MPI
156     return 1;
157 #else
158     int i;
159     (void) MPI_Comm_size(MPI_COMM_WORLD, &i);
160     return i;
161 #endif
162 }
163
164 int gmx_node_rank(void)
165 {
166 #ifndef GMX_MPI
167     return 0;
168 #else
169     int i;
170     (void) MPI_Comm_rank(MPI_COMM_WORLD, &i);
171     return i;
172 #endif
173 }
174
175
176 int gmx_hostname_num()
177 {
178 #ifndef GMX_MPI
179     return 0;
180 #else
181 #ifdef GMX_THREAD_MPI
182     /* thread-MPI currently puts the thread number in the process name,
183      * we might want to change this, as this is inconsistent with what
184      * most MPI implementations would do when running on a single node.
185      */
186     return 0;
187 #else
188     int  resultlen, hostnum, i, j;
189     char mpi_hostname[MPI_MAX_PROCESSOR_NAME], hostnum_str[MPI_MAX_PROCESSOR_NAME];
190
191     MPI_Get_processor_name(mpi_hostname, &resultlen);
192     /* This procedure can only differentiate nodes with host names
193      * that end on unique numbers.
194      */
195     i = 0;
196     j = 0;
197     /* Only parse the host name up to the first dot */
198     while (i < resultlen && mpi_hostname[i] != '.')
199     {
200         if (isdigit(mpi_hostname[i]))
201         {
202             hostnum_str[j++] = mpi_hostname[i];
203         }
204         i++;
205     }
206     hostnum_str[j] = '\0';
207     if (j == 0)
208     {
209         hostnum = 0;
210     }
211     else
212     {
213         /* Use only the last 9 decimals, so we don't overflow an int */
214         hostnum = strtol(hostnum_str + max(0, j-9), NULL, 10);
215     }
216
217     if (debug)
218     {
219         fprintf(debug, "In gmx_setup_nodecomm: hostname '%s', hostnum %d\n",
220                 mpi_hostname, hostnum);
221     }
222     return hostnum;
223 #endif
224 #endif
225 }
226
227 void gmx_setup_nodecomm(FILE gmx_unused *fplog, t_commrec *cr)
228 {
229     gmx_nodecomm_t *nc;
230     int             n, rank, hostnum, ng, ni;
231
232     /* Many MPI implementations do not optimize MPI_Allreduce
233      * (and probably also other global communication calls)
234      * for multi-core nodes connected by a network.
235      * We can optimize such communication by using one MPI call
236      * within each node and one between the nodes.
237      * For MVAPICH2 and Intel MPI this reduces the time for
238      * the global_stat communication by 25%
239      * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
240      * B. Hess, November 2007
241      */
242
243     nc = &cr->nc;
244
245     nc->bUse = FALSE;
246 #ifndef GMX_THREAD_MPI
247 #ifdef GMX_MPI
248     MPI_Comm_size(cr->mpi_comm_mygroup, &n);
249     MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
250
251     hostnum = gmx_hostname_num();
252
253     if (debug)
254     {
255         fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
256     }
257
258
259     /* The intra-node communicator, split on node number */
260     MPI_Comm_split(cr->mpi_comm_mygroup, hostnum, rank, &nc->comm_intra);
261     MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
262     if (debug)
263     {
264         fprintf(debug, "In gmx_setup_nodecomm: node rank %d rank_intra %d\n",
265                 rank, nc->rank_intra);
266     }
267     /* The inter-node communicator, split on rank_intra.
268      * We actually only need the one for rank=0,
269      * but it is easier to create them all.
270      */
271     MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
272     /* Check if this really created two step communication */
273     MPI_Comm_size(nc->comm_inter, &ng);
274     MPI_Comm_size(nc->comm_intra, &ni);
275     if (debug)
276     {
277         fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
278                 ng, ni);
279     }
280
281     if (getenv("GMX_NO_NODECOMM") == NULL &&
282         ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
283     {
284         nc->bUse = TRUE;
285         if (fplog)
286         {
287             fprintf(fplog, "Using two step summing over %d groups of on average %.1f processes\n\n",
288                     ng, (real)n/(real)ng);
289         }
290         if (nc->rank_intra > 0)
291         {
292             MPI_Comm_free(&nc->comm_inter);
293         }
294     }
295     else
296     {
297         /* One group or all processes in a separate group, use normal summing */
298         MPI_Comm_free(&nc->comm_inter);
299         MPI_Comm_free(&nc->comm_intra);
300         if (debug)
301         {
302             fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
303         }
304     }
305 #endif
306 #else
307     /* tMPI runs only on a single node so just use the nodeid */
308     nc->rank_intra = cr->nodeid;
309 #endif
310 }
311
312 void gmx_init_intranode_counters(t_commrec *cr)
313 {
314     /* counters for PP+PME and PP-only processes on my physical node */
315     int nrank_intranode, rank_intranode;
316     int nrank_pp_intranode, rank_pp_intranode;
317     /* thread-MPI is not initialized when not running in parallel */
318 #if defined GMX_MPI && !defined GMX_THREAD_MPI
319     int nrank_world, rank_world;
320     int i, mynum, *num, *num_s, *num_pp, *num_pp_s;
321
322     MPI_Comm_size(MPI_COMM_WORLD, &nrank_world);
323     MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
324
325     /* Get the node number from the hostname to identify the nodes */
326     mynum = gmx_hostname_num();
327
328     /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
329     snew(num,   nrank_world);
330     snew(num_s, nrank_world);
331     snew(num_pp,   nrank_world);
332     snew(num_pp_s, nrank_world);
333
334     num_s[rank_world]    = mynum;
335     num_pp_s[rank_world] = (cr->duty & DUTY_PP) ? mynum : -1;
336
337     MPI_Allreduce(num_s,    num,    nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
338     MPI_Allreduce(num_pp_s, num_pp, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
339
340     nrank_intranode    = 0;
341     rank_intranode     = 0;
342     nrank_pp_intranode = 0;
343     rank_pp_intranode  = 0;
344     for (i = 0; i < nrank_world; i++)
345     {
346         if (num[i] == mynum)
347         {
348             nrank_intranode++;
349             if (i < rank_world)
350             {
351                 rank_intranode++;
352             }
353         }
354         if ((cr->duty & DUTY_PP) && num_pp[i] == mynum)
355         {
356             nrank_pp_intranode++;
357             if (i < rank_world)
358             {
359                 rank_pp_intranode++;
360             }
361         }
362     }
363     sfree(num);
364     sfree(num_s);
365     sfree(num_pp);
366     sfree(num_pp_s);
367 #else
368     /* Serial or thread-MPI code: we run within a single physical node */
369     nrank_intranode    = cr->nnodes;
370     rank_intranode     = cr->sim_nodeid;
371     nrank_pp_intranode = cr->nnodes - cr->npmenodes;
372     rank_pp_intranode  = cr->nodeid;
373 #endif
374
375     if (debug)
376     {
377         char sbuf[STRLEN];
378         if (cr->duty & DUTY_PP && cr->duty & DUTY_PME)
379         {
380             sprintf(sbuf, "PP+PME");
381         }
382         else
383         {
384             sprintf(sbuf, "%s", cr->duty & DUTY_PP ? "PP" : "PME");
385         }
386         fprintf(debug, "On %3s node %d: nrank_intranode=%d, rank_intranode=%d, "
387                 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
388                 sbuf, cr->sim_nodeid,
389                 nrank_intranode, rank_intranode,
390                 nrank_pp_intranode, rank_pp_intranode);
391     }
392
393     cr->nrank_intranode    = nrank_intranode;
394     cr->rank_intranode     = rank_intranode;
395     cr->nrank_pp_intranode = nrank_pp_intranode;
396     cr->rank_pp_intranode  = rank_pp_intranode;
397 }
398
399
400 void gmx_barrier(const t_commrec *cr)
401 {
402 #ifndef GMX_MPI
403     gmx_call("gmx_barrier");
404 #else
405     MPI_Barrier(cr->mpi_comm_mygroup);
406 #endif
407 }
408
409 void gmx_abort(int gmx_unused noderank, int gmx_unused nnodes, int gmx_unused errorno)
410 {
411 #ifndef GMX_MPI
412     gmx_call("gmx_abort");
413 #else
414 #ifdef GMX_THREAD_MPI
415     fprintf(stderr, "Halting program %s\n", ShortProgram());
416     gmx_thanx(stderr);
417     exit(1);
418 #else
419     if (nnodes > 1)
420     {
421         fprintf(stderr, "Halting parallel program %s on CPU %d out of %d\n",
422                 ShortProgram(), noderank, nnodes);
423     }
424     else
425     {
426         fprintf(stderr, "Halting program %s\n", ShortProgram());
427     }
428
429     gmx_thanx(stderr);
430     MPI_Abort(MPI_COMM_WORLD, errorno);
431     exit(1);
432 #endif
433 #endif
434 }
435
436 void gmx_bcast(int nbytes, void *b, const t_commrec *cr)
437 {
438 #ifndef GMX_MPI
439     gmx_call("gmx_bast");
440 #else
441     MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
442 #endif
443 }
444
445 void gmx_bcast_sim(int nbytes, void *b, const t_commrec *cr)
446 {
447 #ifndef GMX_MPI
448     gmx_call("gmx_bast");
449 #else
450     MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
451 #endif
452 }
453
454 void gmx_sumd(int nr, double r[], const t_commrec *cr)
455 {
456 #ifndef GMX_MPI
457     gmx_call("gmx_sumd");
458 #else
459 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
460     if (cr->nc.bUse)
461     {
462         if (cr->nc.rank_intra == 0)
463         {
464             /* Use two step summing. */
465             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
466                        cr->nc.comm_intra);
467             /* Sum the roots of the internal (intra) buffers. */
468             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
469                           cr->nc.comm_inter);
470         }
471         else
472         {
473             /* This is here because of the silly MPI specification
474                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
475             MPI_Reduce(r, NULL, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
476         }
477         MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
478     }
479     else
480     {
481         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
482                       cr->mpi_comm_mygroup);
483     }
484 #else
485     int i;
486
487     if (nr > cr->mpb->dbuf_alloc)
488     {
489         cr->mpb->dbuf_alloc = nr;
490         srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
491     }
492     if (cr->nc.bUse)
493     {
494         /* Use two step summing */
495         MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
496         if (cr->nc.rank_intra == 0)
497         {
498             /* Sum with the buffers reversed */
499             MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
500                           cr->nc.comm_inter);
501         }
502         MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
503     }
504     else
505     {
506         MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
507                       cr->mpi_comm_mygroup);
508         for (i = 0; i < nr; i++)
509         {
510             r[i] = cr->mpb->dbuf[i];
511         }
512     }
513 #endif
514 #endif
515 }
516
517 void gmx_sumf(int nr, float r[], const t_commrec *cr)
518 {
519 #ifndef GMX_MPI
520     gmx_call("gmx_sumf");
521 #else
522 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
523     if (cr->nc.bUse)
524     {
525         /* Use two step summing.  */
526         if (cr->nc.rank_intra == 0)
527         {
528             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
529                        cr->nc.comm_intra);
530             /* Sum the roots of the internal (intra) buffers */
531             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
532                           cr->nc.comm_inter);
533         }
534         else
535         {
536             /* This is here because of the silly MPI specification
537                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
538             MPI_Reduce(r, NULL, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
539         }
540         MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
541     }
542     else
543     {
544         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
545     }
546 #else
547     int i;
548
549     if (nr > cr->mpb->fbuf_alloc)
550     {
551         cr->mpb->fbuf_alloc = nr;
552         srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
553     }
554     if (cr->nc.bUse)
555     {
556         /* Use two step summing */
557         MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
558         if (cr->nc.rank_intra == 0)
559         {
560             /* Sum with the buffers reversed */
561             MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
562                           cr->nc.comm_inter);
563         }
564         MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
565     }
566     else
567     {
568         MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
569                       cr->mpi_comm_mygroup);
570         for (i = 0; i < nr; i++)
571         {
572             r[i] = cr->mpb->fbuf[i];
573         }
574     }
575 #endif
576 #endif
577 }
578
579 void gmx_sumi(int nr, int r[], const t_commrec *cr)
580 {
581 #ifndef GMX_MPI
582     gmx_call("gmx_sumi");
583 #else
584 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
585     if (cr->nc.bUse)
586     {
587         /* Use two step summing */
588         if (cr->nc.rank_intra == 0)
589         {
590             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
591             /* Sum with the buffers reversed */
592             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
593         }
594         else
595         {
596             /* This is here because of the silly MPI specification
597                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
598             MPI_Reduce(r, NULL, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
599         }
600         MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
601     }
602     else
603     {
604         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
605     }
606 #else
607     int i;
608
609     if (nr > cr->mpb->ibuf_alloc)
610     {
611         cr->mpb->ibuf_alloc = nr;
612         srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
613     }
614     if (cr->nc.bUse)
615     {
616         /* Use two step summing */
617         MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
618         if (cr->nc.rank_intra == 0)
619         {
620             /* Sum with the buffers reversed */
621             MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
622         }
623         MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
624     }
625     else
626     {
627         MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
628         for (i = 0; i < nr; i++)
629         {
630             r[i] = cr->mpb->ibuf[i];
631         }
632     }
633 #endif
634 #endif
635 }
636
637 void gmx_sumli(int nr, gmx_large_int_t r[], const t_commrec *cr)
638 {
639 #ifndef GMX_MPI
640     gmx_call("gmx_sumli");
641 #else
642 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
643     if (cr->nc.bUse)
644     {
645         /* Use two step summing */
646         if (cr->nc.rank_intra == 0)
647         {
648             MPI_Reduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM, 0,
649                        cr->nc.comm_intra);
650             /* Sum with the buffers reversed */
651             MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
652                           cr->nc.comm_inter);
653         }
654         else
655         {
656             /* This is here because of the silly MPI specification
657                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
658             MPI_Reduce(r, NULL, nr, GMX_MPI_LARGE_INT, MPI_SUM, 0, cr->nc.comm_intra);
659         }
660         MPI_Bcast(r, nr, GMX_MPI_LARGE_INT, 0, cr->nc.comm_intra);
661     }
662     else
663     {
664         MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM, cr->mpi_comm_mygroup);
665     }
666 #else
667     int i;
668
669     if (nr > cr->mpb->libuf_alloc)
670     {
671         cr->mpb->libuf_alloc = nr;
672         srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
673     }
674     if (cr->nc.bUse)
675     {
676         /* Use two step summing */
677         MPI_Allreduce(r, cr->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
678                       cr->nc.comm_intra);
679         if (cr->nc.rank_intra == 0)
680         {
681             /* Sum with the buffers reversed */
682             MPI_Allreduce(cr->mpb->libuf, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
683                           cr->nc.comm_inter);
684         }
685         MPI_Bcast(r, nr, GMX_MPI_LARGE_INT, 0, cr->nc.comm_intra);
686     }
687     else
688     {
689         MPI_Allreduce(r, cr->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
690                       cr->mpi_comm_mygroup);
691         for (i = 0; i < nr; i++)
692         {
693             r[i] = cr->mpb->libuf[i];
694         }
695     }
696 #endif
697 #endif
698 }
699
700
701
702 #ifdef GMX_MPI
703 void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
704 {
705 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
706     MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
707 #else
708     /* this function is only used in code that is not performance critical,
709        (during setup, when comm_rec is not the appropriate communication
710        structure), so this isn't as bad as it looks. */
711     double *buf;
712     int     i;
713
714     snew(buf, nr);
715     MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
716     for (i = 0; i < nr; i++)
717     {
718         r[i] = buf[i];
719     }
720     sfree(buf);
721 #endif
722 }
723 #endif
724
725 #ifdef GMX_MPI
726 void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
727 {
728 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
729     MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
730 #else
731     /* this function is only used in code that is not performance critical,
732        (during setup, when comm_rec is not the appropriate communication
733        structure), so this isn't as bad as it looks. */
734     float *buf;
735     int    i;
736
737     snew(buf, nr);
738     MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
739     for (i = 0; i < nr; i++)
740     {
741         r[i] = buf[i];
742     }
743     sfree(buf);
744 #endif
745 }
746 #endif
747
748 void gmx_sumd_sim(int nr, double r[], const gmx_multisim_t *ms)
749 {
750 #ifndef GMX_MPI
751     gmx_call("gmx_sumd_sim");
752 #else
753     gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
754 #endif
755 }
756
757 void gmx_sumf_sim(int nr, float r[], const gmx_multisim_t *ms)
758 {
759 #ifndef GMX_MPI
760     gmx_call("gmx_sumf_sim");
761 #else
762     gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
763 #endif
764 }
765
766 void gmx_sumi_sim(int nr, int r[], const gmx_multisim_t *ms)
767 {
768 #ifndef GMX_MPI
769     gmx_call("gmx_sumi_sim");
770 #else
771 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
772     MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
773 #else
774     /* this is thread-unsafe, but it will do for now: */
775     int i;
776
777     if (nr > ms->mpb->ibuf_alloc)
778     {
779         ms->mpb->ibuf_alloc = nr;
780         srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
781     }
782     MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
783     for (i = 0; i < nr; i++)
784     {
785         r[i] = ms->mpb->ibuf[i];
786     }
787 #endif
788 #endif
789 }
790
791 void gmx_sumli_sim(int nr, gmx_large_int_t r[], const gmx_multisim_t *ms)
792 {
793 #ifndef GMX_MPI
794     gmx_call("gmx_sumli_sim");
795 #else
796 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
797     MPI_Allreduce(MPI_IN_PLACE, r, nr, GMX_MPI_LARGE_INT, MPI_SUM,
798                   ms->mpi_comm_masters);
799 #else
800     /* this is thread-unsafe, but it will do for now: */
801     int i;
802
803     if (nr > ms->mpb->libuf_alloc)
804     {
805         ms->mpb->libuf_alloc = nr;
806         srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
807     }
808     MPI_Allreduce(r, ms->mpb->libuf, nr, GMX_MPI_LARGE_INT, MPI_SUM,
809                   ms->mpi_comm_masters);
810     for (i = 0; i < nr; i++)
811     {
812         r[i] = ms->mpb->libuf[i];
813     }
814 #endif
815 #endif
816 }