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