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