Fixes and updates to BlueGene/Q support
[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 <string.h>
42 #include "gmx_fatal.h"
43 #include "main.h"
44 #include "gromacs/utility/smalloc.h"
45 #include "types/commrec.h"
46 #include "network.h"
47 #include "copyrite.h"
48 #include <ctype.h>
49 #include "macros.h"
50 #include "gromacs/utility/cstringutil.h"
51
52 #include "gromacs/utility/gmxmpi.h"
53
54
55 /* The source code in this file should be thread-safe.
56       Please keep it that way. */
57
58 gmx_bool gmx_mpi_initialized(void)
59 {
60     int n;
61 #ifndef GMX_MPI
62     return 0;
63 #else
64     MPI_Initialized(&n);
65
66     return n;
67 #endif
68 }
69
70 void gmx_fill_commrec_from_mpi(t_commrec gmx_unused *cr)
71 {
72 #ifndef GMX_MPI
73     gmx_call("gmx_fill_commrec_from_mpi");
74 #else
75     if (!gmx_mpi_initialized())
76     {
77         gmx_comm("MPI has not been initialized properly");
78     }
79
80     cr->nnodes           = gmx_node_num();
81     cr->nodeid           = gmx_node_rank();
82     cr->sim_nodeid       = cr->nodeid;
83     cr->mpi_comm_mysim   = MPI_COMM_WORLD;
84     cr->mpi_comm_mygroup = MPI_COMM_WORLD;
85
86 #endif
87 }
88
89 t_commrec *init_commrec()
90 {
91     t_commrec    *cr;
92
93     snew(cr, 1);
94
95 #ifdef GMX_LIB_MPI
96     gmx_fill_commrec_from_mpi(cr);
97 #else
98     cr->mpi_comm_mysim   = NULL;
99     cr->mpi_comm_mygroup = NULL;
100     cr->nnodes           = 1;
101     cr->sim_nodeid       = 0;
102     cr->nodeid           = cr->sim_nodeid;
103 #endif
104
105     // TODO cr->duty should not be initialized here
106     cr->duty = (DUTY_PP | DUTY_PME);
107
108 #if defined GMX_MPI && !defined MPI_IN_PLACE_EXISTS
109     /* initialize the MPI_IN_PLACE replacement buffers */
110     snew(cr->mpb, 1);
111     cr->mpb->ibuf        = NULL;
112     cr->mpb->libuf       = NULL;
113     cr->mpb->fbuf        = NULL;
114     cr->mpb->dbuf        = NULL;
115     cr->mpb->ibuf_alloc  = 0;
116     cr->mpb->libuf_alloc = 0;
117     cr->mpb->fbuf_alloc  = 0;
118     cr->mpb->dbuf_alloc  = 0;
119 #endif
120
121     return cr;
122 }
123
124 t_commrec *reinitialize_commrec_for_this_thread(const t_commrec gmx_unused *cro)
125 {
126 #ifdef GMX_THREAD_MPI
127     t_commrec *cr;
128
129     /* make a thread-specific commrec */
130     snew(cr, 1);
131     /* now copy the whole thing, so settings like the number of PME nodes
132        get propagated. */
133     *cr = *cro;
134
135     /* and we start setting our own thread-specific values for things */
136     gmx_fill_commrec_from_mpi(cr);
137
138     // TODO cr->duty should not be initialized here
139     cr->duty             = (DUTY_PP | DUTY_PME);
140
141     return cr;
142 #else
143     return NULL;
144 #endif
145 }
146
147 int  gmx_node_num(void)
148 {
149 #ifndef GMX_MPI
150     return 1;
151 #else
152     int i;
153     (void) MPI_Comm_size(MPI_COMM_WORLD, &i);
154     return i;
155 #endif
156 }
157
158 int gmx_node_rank(void)
159 {
160 #ifndef GMX_MPI
161     return 0;
162 #else
163     int i;
164     (void) MPI_Comm_rank(MPI_COMM_WORLD, &i);
165     return i;
166 #endif
167 }
168
169 static int mpi_hostname_hash(void)
170 {
171     int hash_int;
172
173 #ifndef GMX_LIB_MPI
174     /* We have a single physical node */
175     hash_int = 0;
176 #else
177     int  resultlen;
178     char mpi_hostname[MPI_MAX_PROCESSOR_NAME];
179
180     /* This procedure can only differentiate nodes with different names.
181      * Architectures where different physical nodes have identical names,
182      * such as IBM Blue Gene, should use an architecture specific solution.
183      */
184     MPI_Get_processor_name(mpi_hostname, &resultlen);
185
186     /* The string hash function returns an unsigned int. We cast to an int.
187      * Negative numbers are converted to positive by setting the sign bit to 0.
188      * This makes the hash one bit smaller.
189      * A 63-bit hash (with 64-bit int) should be enough for unique node hashes,
190      * even on a million node machine. 31 bits might not be enough though!
191      */
192     hash_int =
193         (int)gmx_string_fullhash_func(mpi_hostname, gmx_string_hash_init);
194     if (hash_int < 0)
195     {
196         hash_int -= INT_MIN;
197     }
198 #endif
199
200     return hash_int;
201 }
202
203 #if defined GMX_LIB_MPI && defined GMX_TARGET_BGQ
204
205 #ifdef __clang__
206 /* IBM's declaration of this function in
207  * /bgsys/drivers/V1R2M2/ppc64/spi/include/kernel/process.h
208  * erroneously fails to specify __INLINE__, despite
209  * /bgsys/drivers/V1R2M2/ppc64/spi/include/kernel/cnk/process_impl.h
210  * specifiying __INLINE__, so bgclang thinks they are different enough
211  * to complain about. */
212 static uint64_t Kernel_GetJobID();
213 #endif
214
215 #include <spi/include/kernel/location.h>
216
217 static int bgq_nodenum(void)
218 {
219     int           hostnum;
220     Personality_t personality;
221     Kernel_GetPersonality(&personality, sizeof(personality));
222     /* Each MPI rank has a unique coordinate in a 6-dimensional space
223        (A,B,C,D,E,T), with dimensions A-E corresponding to different
224        physical nodes, and T within each node. Each node has sixteen
225        physical cores, each of which can have up to four hardware
226        threads, so 0 <= T <= 63 (but the maximum value of T depends on
227        the confituration of ranks and OpenMP threads per
228        node). However, T is irrelevant for computing a suitable return
229        value for gmx_hostname_num().
230      */
231     hostnum  = personality.Network_Config.Acoord;
232     hostnum *= personality.Network_Config.Bnodes;
233     hostnum += personality.Network_Config.Bcoord;
234     hostnum *= personality.Network_Config.Cnodes;
235     hostnum += personality.Network_Config.Ccoord;
236     hostnum *= personality.Network_Config.Dnodes;
237     hostnum += personality.Network_Config.Dcoord;
238     hostnum *= personality.Network_Config.Enodes;
239     hostnum += personality.Network_Config.Ecoord;
240
241     if (debug)
242     {
243         fprintf(debug,
244                 "Torus ID A: %d / %d B: %d / %d C: %d / %d D: %d / %d E: %d / %d\nNode ID T: %d / %d core: %d / %d hardware thread: %d / %d\n",
245                 personality.Network_Config.Acoord,
246                 personality.Network_Config.Anodes,
247                 personality.Network_Config.Bcoord,
248                 personality.Network_Config.Bnodes,
249                 personality.Network_Config.Ccoord,
250                 personality.Network_Config.Cnodes,
251                 personality.Network_Config.Dcoord,
252                 personality.Network_Config.Dnodes,
253                 personality.Network_Config.Ecoord,
254                 personality.Network_Config.Enodes,
255                 Kernel_ProcessorCoreID(),
256                 16,
257                 Kernel_ProcessorID(),
258                 64,
259                 Kernel_ProcessorThreadID(),
260                 4);
261     }
262     return hostnum;
263 }
264 #endif
265
266 int gmx_physicalnode_id_hash(void)
267 {
268     int hash;
269
270 #ifndef GMX_MPI
271     hash = 0;
272 #else
273 #ifdef GMX_THREAD_MPI
274     /* thread-MPI currently puts the thread number in the process name,
275      * we might want to change this, as this is inconsistent with what
276      * most MPI implementations would do when running on a single node.
277      */
278     hash = 0;
279 #else
280 #ifdef GMX_TARGET_BGQ
281     hash = bgq_nodenum();
282 #else
283     hash = mpi_hostname_hash();
284 #endif
285 #endif
286 #endif
287
288     if (debug)
289     {
290         fprintf(debug, "In gmx_physicalnode_id_hash: hash %d\n", hash);
291     }
292
293     return hash;
294 }
295
296 void gmx_setup_nodecomm(FILE gmx_unused *fplog, t_commrec *cr)
297 {
298     gmx_nodecomm_t *nc;
299     int             n, rank, nodehash, ng, ni;
300
301     /* Many MPI implementations do not optimize MPI_Allreduce
302      * (and probably also other global communication calls)
303      * for multi-core nodes connected by a network.
304      * We can optimize such communication by using one MPI call
305      * within each node and one between the nodes.
306      * For MVAPICH2 and Intel MPI this reduces the time for
307      * the global_stat communication by 25%
308      * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
309      * B. Hess, November 2007
310      */
311
312     nc = &cr->nc;
313
314     nc->bUse = FALSE;
315 #ifndef GMX_THREAD_MPI
316 #ifdef GMX_MPI
317     MPI_Comm_size(cr->mpi_comm_mygroup, &n);
318     MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
319
320     nodehash = gmx_physicalnode_id_hash();
321
322     if (debug)
323     {
324         fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
325     }
326
327
328     /* The intra-node communicator, split on node number */
329     MPI_Comm_split(cr->mpi_comm_mygroup, nodehash, rank, &nc->comm_intra);
330     MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
331     if (debug)
332     {
333         fprintf(debug, "In gmx_setup_nodecomm: node ID %d rank within node %d\n",
334                 rank, nc->rank_intra);
335     }
336     /* The inter-node communicator, split on rank_intra.
337      * We actually only need the one for rank=0,
338      * but it is easier to create them all.
339      */
340     MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
341     /* Check if this really created two step communication */
342     MPI_Comm_size(nc->comm_inter, &ng);
343     MPI_Comm_size(nc->comm_intra, &ni);
344     if (debug)
345     {
346         fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n",
347                 ng, ni);
348     }
349
350     if (getenv("GMX_NO_NODECOMM") == NULL &&
351         ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
352     {
353         nc->bUse = TRUE;
354         if (fplog)
355         {
356             fprintf(fplog, "Using two step summing over %d groups of on average %.1f ranks\n\n",
357                     ng, (real)n/(real)ng);
358         }
359         if (nc->rank_intra > 0)
360         {
361             MPI_Comm_free(&nc->comm_inter);
362         }
363     }
364     else
365     {
366         /* One group or all processes in a separate group, use normal summing */
367         MPI_Comm_free(&nc->comm_inter);
368         MPI_Comm_free(&nc->comm_intra);
369         if (debug)
370         {
371             fprintf(debug, "In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
372         }
373     }
374 #endif
375 #else
376     /* tMPI runs only on a single node so just use the nodeid */
377     nc->rank_intra = cr->nodeid;
378 #endif
379 }
380
381 void gmx_init_intranode_counters(t_commrec *cr)
382 {
383     /* counters for PP+PME and PP-only processes on my physical node */
384     int nrank_intranode, rank_intranode;
385     int nrank_pp_intranode, rank_pp_intranode;
386     /* thread-MPI is not initialized when not running in parallel */
387 #if defined GMX_MPI && !defined GMX_THREAD_MPI
388     int nrank_world, rank_world;
389     int i, myhash, *hash, *hash_s, *hash_pp, *hash_pp_s;
390
391     MPI_Comm_size(MPI_COMM_WORLD, &nrank_world);
392     MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
393
394     /* Get a (hopefully unique) hash that identifies our physical node */
395     myhash = gmx_physicalnode_id_hash();
396
397     /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
398     snew(hash,   nrank_world);
399     snew(hash_s, nrank_world);
400     snew(hash_pp,   nrank_world);
401     snew(hash_pp_s, nrank_world);
402
403     hash_s[rank_world]    = myhash;
404     hash_pp_s[rank_world] = (cr->duty & DUTY_PP) ? myhash : -1;
405
406     MPI_Allreduce(hash_s,    hash,    nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
407     MPI_Allreduce(hash_pp_s, hash_pp, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
408
409     nrank_intranode    = 0;
410     rank_intranode     = 0;
411     nrank_pp_intranode = 0;
412     rank_pp_intranode  = 0;
413     for (i = 0; i < nrank_world; i++)
414     {
415         if (hash[i] == myhash)
416         {
417             nrank_intranode++;
418             if (i < rank_world)
419             {
420                 rank_intranode++;
421             }
422         }
423         if (hash_pp[i] == myhash)
424         {
425             nrank_pp_intranode++;
426             if ((cr->duty & DUTY_PP) && i < rank_world)
427             {
428                 rank_pp_intranode++;
429             }
430         }
431     }
432     sfree(hash);
433     sfree(hash_s);
434     sfree(hash_pp);
435     sfree(hash_pp_s);
436 #else
437     /* Serial or thread-MPI code: we run within a single physical node */
438     nrank_intranode    = cr->nnodes;
439     rank_intranode     = cr->sim_nodeid;
440     nrank_pp_intranode = cr->nnodes - cr->npmenodes;
441     rank_pp_intranode  = cr->nodeid;
442 #endif
443
444     if (debug)
445     {
446         char sbuf[STRLEN];
447         if (cr->duty & DUTY_PP && cr->duty & DUTY_PME)
448         {
449             sprintf(sbuf, "PP+PME");
450         }
451         else
452         {
453             sprintf(sbuf, "%s", cr->duty & DUTY_PP ? "PP" : "PME");
454         }
455         fprintf(debug, "On %3s rank %d: nrank_intranode=%d, rank_intranode=%d, "
456                 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
457                 sbuf, cr->sim_nodeid,
458                 nrank_intranode, rank_intranode,
459                 nrank_pp_intranode, rank_pp_intranode);
460     }
461
462     cr->nrank_intranode    = nrank_intranode;
463     cr->rank_intranode     = rank_intranode;
464     cr->nrank_pp_intranode = nrank_pp_intranode;
465     cr->rank_pp_intranode  = rank_pp_intranode;
466 }
467
468
469 void gmx_barrier(const t_commrec gmx_unused *cr)
470 {
471 #ifndef GMX_MPI
472     gmx_call("gmx_barrier");
473 #else
474     MPI_Barrier(cr->mpi_comm_mygroup);
475 #endif
476 }
477
478 void gmx_abort(int gmx_unused noderank, int gmx_unused nnodes, int gmx_unused errorno)
479 {
480 #ifndef GMX_MPI
481     gmx_call("gmx_abort");
482 #else
483 #ifdef GMX_THREAD_MPI
484     fprintf(stderr, "Halting program %s\n", ShortProgram());
485     gmx_thanx(stderr);
486     exit(1);
487 #else
488     if (nnodes > 1)
489     {
490         fprintf(stderr, "Halting parallel program %s on CPU %d out of %d\n",
491                 ShortProgram(), noderank, nnodes);
492     }
493     else
494     {
495         fprintf(stderr, "Halting program %s\n", ShortProgram());
496     }
497
498     gmx_thanx(stderr);
499     MPI_Abort(MPI_COMM_WORLD, errorno);
500     exit(1);
501 #endif
502 #endif
503 }
504
505 void gmx_bcast(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
506 {
507 #ifndef GMX_MPI
508     gmx_call("gmx_bast");
509 #else
510     MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
511 #endif
512 }
513
514 void gmx_bcast_sim(int gmx_unused nbytes, void gmx_unused *b, const t_commrec gmx_unused *cr)
515 {
516 #ifndef GMX_MPI
517     gmx_call("gmx_bast");
518 #else
519     MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
520 #endif
521 }
522
523 void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unused *cr)
524 {
525 #ifndef GMX_MPI
526     gmx_call("gmx_sumd");
527 #else
528 #if defined(MPI_IN_PLACE_EXISTS)
529     if (cr->nc.bUse)
530     {
531         if (cr->nc.rank_intra == 0)
532         {
533             /* Use two step summing. */
534             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0,
535                        cr->nc.comm_intra);
536             /* Sum the roots of the internal (intra) buffers. */
537             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
538                           cr->nc.comm_inter);
539         }
540         else
541         {
542             /* This is here because of the silly MPI specification
543                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
544             MPI_Reduce(r, NULL, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
545         }
546         MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
547     }
548     else
549     {
550         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM,
551                       cr->mpi_comm_mygroup);
552     }
553 #else
554     int i;
555
556     if (nr > cr->mpb->dbuf_alloc)
557     {
558         cr->mpb->dbuf_alloc = nr;
559         srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
560     }
561     if (cr->nc.bUse)
562     {
563         /* Use two step summing */
564         MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
565         if (cr->nc.rank_intra == 0)
566         {
567             /* Sum with the buffers reversed */
568             MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM,
569                           cr->nc.comm_inter);
570         }
571         MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
572     }
573     else
574     {
575         MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM,
576                       cr->mpi_comm_mygroup);
577         for (i = 0; i < nr; i++)
578         {
579             r[i] = cr->mpb->dbuf[i];
580         }
581     }
582 #endif
583 #endif
584 }
585
586 void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unused *cr)
587 {
588 #ifndef GMX_MPI
589     gmx_call("gmx_sumf");
590 #else
591 #if defined(MPI_IN_PLACE_EXISTS)
592     if (cr->nc.bUse)
593     {
594         /* Use two step summing.  */
595         if (cr->nc.rank_intra == 0)
596         {
597             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0,
598                        cr->nc.comm_intra);
599             /* Sum the roots of the internal (intra) buffers */
600             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM,
601                           cr->nc.comm_inter);
602         }
603         else
604         {
605             /* This is here because of the silly MPI specification
606                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
607             MPI_Reduce(r, NULL, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
608         }
609         MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
610     }
611     else
612     {
613         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
614     }
615 #else
616     int i;
617
618     if (nr > cr->mpb->fbuf_alloc)
619     {
620         cr->mpb->fbuf_alloc = nr;
621         srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
622     }
623     if (cr->nc.bUse)
624     {
625         /* Use two step summing */
626         MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
627         if (cr->nc.rank_intra == 0)
628         {
629             /* Sum with the buffers reversed */
630             MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM,
631                           cr->nc.comm_inter);
632         }
633         MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
634     }
635     else
636     {
637         MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM,
638                       cr->mpi_comm_mygroup);
639         for (i = 0; i < nr; i++)
640         {
641             r[i] = cr->mpb->fbuf[i];
642         }
643     }
644 #endif
645 #endif
646 }
647
648 void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused *cr)
649 {
650 #ifndef GMX_MPI
651     gmx_call("gmx_sumi");
652 #else
653 #if defined(MPI_IN_PLACE_EXISTS)
654     if (cr->nc.bUse)
655     {
656         /* Use two step summing */
657         if (cr->nc.rank_intra == 0)
658         {
659             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
660             /* Sum with the buffers reversed */
661             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
662         }
663         else
664         {
665             /* This is here because of the silly MPI specification
666                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
667             MPI_Reduce(r, NULL, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
668         }
669         MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
670     }
671     else
672     {
673         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
674     }
675 #else
676     int i;
677
678     if (nr > cr->mpb->ibuf_alloc)
679     {
680         cr->mpb->ibuf_alloc = nr;
681         srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
682     }
683     if (cr->nc.bUse)
684     {
685         /* Use two step summing */
686         MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
687         if (cr->nc.rank_intra == 0)
688         {
689             /* Sum with the buffers reversed */
690             MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
691         }
692         MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
693     }
694     else
695     {
696         MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
697         for (i = 0; i < nr; i++)
698         {
699             r[i] = cr->mpb->ibuf[i];
700         }
701     }
702 #endif
703 #endif
704 }
705
706 void gmx_sumli(int gmx_unused nr, gmx_int64_t gmx_unused r[], const t_commrec gmx_unused *cr)
707 {
708 #ifndef GMX_MPI
709     gmx_call("gmx_sumli");
710 #else
711 #if defined(MPI_IN_PLACE_EXISTS)
712     if (cr->nc.bUse)
713     {
714         /* Use two step summing */
715         if (cr->nc.rank_intra == 0)
716         {
717             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, 0,
718                        cr->nc.comm_intra);
719             /* Sum with the buffers reversed */
720             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
721                           cr->nc.comm_inter);
722         }
723         else
724         {
725             /* This is here because of the silly MPI specification
726                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
727             MPI_Reduce(r, NULL, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
728         }
729         MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
730     }
731     else
732     {
733         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
734     }
735 #else
736     int i;
737
738     if (nr > cr->mpb->libuf_alloc)
739     {
740         cr->mpb->libuf_alloc = nr;
741         srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
742     }
743     if (cr->nc.bUse)
744     {
745         /* Use two step summing */
746         MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
747                       cr->nc.comm_intra);
748         if (cr->nc.rank_intra == 0)
749         {
750             /* Sum with the buffers reversed */
751             MPI_Allreduce(cr->mpb->libuf, r, nr, MPI_INT64_T, MPI_SUM,
752                           cr->nc.comm_inter);
753         }
754         MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
755     }
756     else
757     {
758         MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
759                       cr->mpi_comm_mygroup);
760         for (i = 0; i < nr; i++)
761         {
762             r[i] = cr->mpb->libuf[i];
763         }
764     }
765 #endif
766 #endif
767 }
768
769
770
771 #ifdef GMX_MPI
772 static void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
773 {
774 #if defined(MPI_IN_PLACE_EXISTS)
775     MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
776 #else
777     /* this function is only used in code that is not performance critical,
778        (during setup, when comm_rec is not the appropriate communication
779        structure), so this isn't as bad as it looks. */
780     double *buf;
781     int     i;
782
783     snew(buf, nr);
784     MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
785     for (i = 0; i < nr; i++)
786     {
787         r[i] = buf[i];
788     }
789     sfree(buf);
790 #endif
791 }
792 #endif
793
794 #ifdef GMX_MPI
795 static void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
796 {
797 #if defined(MPI_IN_PLACE_EXISTS)
798     MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
799 #else
800     /* this function is only used in code that is not performance critical,
801        (during setup, when comm_rec is not the appropriate communication
802        structure), so this isn't as bad as it looks. */
803     float *buf;
804     int    i;
805
806     snew(buf, nr);
807     MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
808     for (i = 0; i < nr; i++)
809     {
810         r[i] = buf[i];
811     }
812     sfree(buf);
813 #endif
814 }
815 #endif
816
817 void gmx_sumd_sim(int gmx_unused nr, double gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
818 {
819 #ifndef GMX_MPI
820     gmx_call("gmx_sumd_sim");
821 #else
822     gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
823 #endif
824 }
825
826 void gmx_sumf_sim(int gmx_unused nr, float gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
827 {
828 #ifndef GMX_MPI
829     gmx_call("gmx_sumf_sim");
830 #else
831     gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
832 #endif
833 }
834
835 void gmx_sumi_sim(int gmx_unused nr, int gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
836 {
837 #ifndef GMX_MPI
838     gmx_call("gmx_sumi_sim");
839 #else
840 #if defined(MPI_IN_PLACE_EXISTS)
841     MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
842 #else
843     /* this is thread-unsafe, but it will do for now: */
844     int i;
845
846     if (nr > ms->mpb->ibuf_alloc)
847     {
848         ms->mpb->ibuf_alloc = nr;
849         srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
850     }
851     MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
852     for (i = 0; i < nr; i++)
853     {
854         r[i] = ms->mpb->ibuf[i];
855     }
856 #endif
857 #endif
858 }
859
860 void gmx_sumli_sim(int gmx_unused nr, gmx_int64_t gmx_unused r[], const gmx_multisim_t gmx_unused *ms)
861 {
862 #ifndef GMX_MPI
863     gmx_call("gmx_sumli_sim");
864 #else
865 #if defined(MPI_IN_PLACE_EXISTS)
866     MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM,
867                   ms->mpi_comm_masters);
868 #else
869     /* this is thread-unsafe, but it will do for now: */
870     int i;
871
872     if (nr > ms->mpb->libuf_alloc)
873     {
874         ms->mpb->libuf_alloc = nr;
875         srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
876     }
877     MPI_Allreduce(r, ms->mpb->libuf, nr, MPI_INT64_T, MPI_SUM,
878                   ms->mpi_comm_masters);
879     for (i = 0; i < nr; i++)
880     {
881         r[i] = ms->mpb->libuf[i];
882     }
883 #endif
884 #endif
885 }