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