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