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