1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
30 * For more info, check our website at http://www.gromacs.org
33 * GROningen Mixture of Alchemy and Childrens' Stories
40 #include "gmx_fatal.h"
58 /* The source code in this file should be thread-safe.
59 Please keep it that way. */
61 gmx_bool gmx_mpi_initialized(void)
73 int gmx_setup(int *argc,char **argv,int *nnodes)
76 gmx_call("gmx_setup");
80 int resultlen; /* actual length of node name */
84 char mpi_hostname[MPI_MAX_PROCESSOR_NAME];
86 /* Call the MPI routines */
89 (void) fah_MPI_Init(argc,&argv);
91 (void) MPI_Init(argc,&argv);
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 );
99 fprintf(stderr,"NNODES=%d, MYRANK=%d, HOSTNAME=%s\n",
100 mpi_num_nodes,mpi_my_rank,mpi_hostname);
103 *nnodes=mpi_num_nodes;
109 int gmx_node_num(void)
115 (void) MPI_Comm_size(MPI_COMM_WORLD, &i);
120 int gmx_node_rank(void)
126 (void) MPI_Comm_rank(MPI_COMM_WORLD, &i);
132 int gmx_hostname_num()
137 #ifdef GMX_THREAD_MPI
138 /* thread-MPI currently puts the thread number in the process name,
139 * we might want to change this, as this is inconsistent with what
140 * most MPI implementations would do when running on a single node.
144 int resultlen,hostnum,i,j;
145 char mpi_hostname[MPI_MAX_PROCESSOR_NAME],hostnum_str[MPI_MAX_PROCESSOR_NAME];
147 MPI_Get_processor_name(mpi_hostname,&resultlen);
148 /* This procedure can only differentiate nodes with host names
149 * that end on unique numbers.
153 /* Only parse the host name up to the first dot */
154 while(i < resultlen && mpi_hostname[i] != '.') {
155 if (isdigit(mpi_hostname[i])) {
156 hostnum_str[j++] = mpi_hostname[i];
160 hostnum_str[j] = '\0';
164 /* Use only the last 9 decimals, so we don't overflow an int */
165 hostnum = strtol(hostnum_str + max(0,j-9), NULL, 10);
169 fprintf(debug,"In gmx_setup_nodecomm: hostname '%s', hostnum %d\n",
170 mpi_hostname,hostnum);
177 void gmx_setup_nodecomm(FILE *fplog,t_commrec *cr)
180 int n,rank,hostnum,ng,ni;
182 /* Many MPI implementations do not optimize MPI_Allreduce
183 * (and probably also other global communication calls)
184 * for multi-core nodes connected by a network.
185 * We can optimize such communication by using one MPI call
186 * within each node and one between the nodes.
187 * For MVAPICH2 and Intel MPI this reduces the time for
188 * the global_stat communication by 25%
189 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
190 * B. Hess, November 2007
196 #ifndef GMX_THREAD_MPI
198 MPI_Comm_size(cr->mpi_comm_mygroup,&n);
199 MPI_Comm_rank(cr->mpi_comm_mygroup,&rank);
201 hostnum = gmx_hostname_num();
205 fprintf(debug,"In gmx_setup_nodecomm: splitting communicator of size %d\n",n);
209 /* The intra-node communicator, split on node number */
210 MPI_Comm_split(cr->mpi_comm_mygroup,hostnum,rank,&nc->comm_intra);
211 MPI_Comm_rank(nc->comm_intra,&nc->rank_intra);
214 fprintf(debug,"In gmx_setup_nodecomm: node rank %d rank_intra %d\n",
215 rank,nc->rank_intra);
217 /* The inter-node communicator, split on rank_intra.
218 * We actually only need the one for rank=0,
219 * but it is easier to create them all.
221 MPI_Comm_split(cr->mpi_comm_mygroup,nc->rank_intra,rank,&nc->comm_inter);
222 /* Check if this really created two step communication */
223 MPI_Comm_size(nc->comm_inter,&ng);
224 MPI_Comm_size(nc->comm_intra,&ni);
227 fprintf(debug,"In gmx_setup_nodecomm: groups %d, my group size %d\n",
231 if (getenv("GMX_NO_NODECOMM") == NULL &&
232 ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
237 fprintf(fplog,"Using two step summing over %d groups of on average %.1f processes\n\n",
238 ng,(real)n/(real)ng);
240 if (nc->rank_intra > 0)
242 MPI_Comm_free(&nc->comm_inter);
247 /* One group or all processes in a separate group, use normal summing */
248 MPI_Comm_free(&nc->comm_inter);
249 MPI_Comm_free(&nc->comm_intra);
252 fprintf(debug,"In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
257 /* tMPI runs only on a single node so just use the nodeid */
258 nc->rank_intra = cr->nodeid;
262 void gmx_init_intranode_counters(t_commrec *cr)
264 /* counters for PP+PME and PP-only processes on my physical node */
265 int nrank_intranode, rank_intranode;
266 int nrank_pp_intranode, rank_pp_intranode;
267 /* thread-MPI is not initialized when not running in parallel */
268 #if defined GMX_MPI && !defined GMX_THREAD_MPI
269 int nrank_world, rank_world;
270 int i, mynum, *num, *num_s, *num_pp, *num_pp_s;
272 MPI_Comm_size(MPI_COMM_WORLD,&nrank_world);
273 MPI_Comm_rank(MPI_COMM_WORLD,&rank_world);
275 /* Get the node number from the hostname to identify the nodes */
276 mynum = gmx_hostname_num();
278 /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
279 snew(num, nrank_world);
280 snew(num_s, nrank_world);
281 snew(num_pp, nrank_world);
282 snew(num_pp_s, nrank_world);
284 num_s[rank_world] = mynum;
285 num_pp_s[rank_world] = (cr->duty & DUTY_PP) ? mynum : -1;
287 MPI_Allreduce(num_s, num, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
288 MPI_Allreduce(num_pp_s, num_pp, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
292 nrank_pp_intranode = 0;
293 rank_pp_intranode = 0;
294 for(i=0; i<nrank_world; i++)
304 if ((cr->duty & DUTY_PP) && num_pp[i] == mynum)
306 nrank_pp_intranode++;
318 /* Serial or thread-MPI code: we run within a single physical node */
319 nrank_intranode = cr->nnodes;
320 rank_intranode = cr->sim_nodeid;
321 nrank_pp_intranode = cr->nnodes - cr->npmenodes;
322 rank_pp_intranode = cr->nodeid;
328 if (cr->duty & DUTY_PP && cr->duty & DUTY_PME)
330 sprintf(sbuf, "PP+PME");
334 sprintf(sbuf, "%s", cr->duty & DUTY_PP ? "PP" : "PME");
336 fprintf(debug, "On %3s node %d: nrank_intranode=%d, rank_intranode=%d, "
337 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
338 sbuf, cr->sim_nodeid,
339 nrank_intranode, rank_intranode,
340 nrank_pp_intranode, rank_pp_intranode);
343 cr->nrank_intranode = nrank_intranode;
344 cr->rank_intranode = rank_intranode;
345 cr->nrank_pp_intranode = nrank_pp_intranode;
346 cr->rank_pp_intranode = rank_pp_intranode;
350 void gmx_barrier(const t_commrec *cr)
353 gmx_call("gmx_barrier");
355 MPI_Barrier(cr->mpi_comm_mygroup);
359 void gmx_abort(int noderank,int nnodes,int errorno)
362 gmx_call("gmx_abort");
364 #ifdef GMX_THREAD_MPI
365 fprintf(stderr,"Halting program %s\n",ShortProgram());
371 fprintf(stderr,"Halting parallel program %s on CPU %d out of %d\n",
372 ShortProgram(),noderank,nnodes);
376 fprintf(stderr,"Halting program %s\n",ShortProgram());
380 MPI_Abort(MPI_COMM_WORLD,errorno);
386 void gmx_bcast(int nbytes,void *b,const t_commrec *cr)
389 gmx_call("gmx_bast");
391 MPI_Bcast(b,nbytes,MPI_BYTE,MASTERRANK(cr),cr->mpi_comm_mygroup);
395 void gmx_bcast_sim(int nbytes,void *b,const t_commrec *cr)
398 gmx_call("gmx_bast");
400 MPI_Bcast(b,nbytes,MPI_BYTE,MASTERRANK(cr),cr->mpi_comm_mysim);
404 void gmx_sumd(int nr,double r[],const t_commrec *cr)
407 gmx_call("gmx_sumd");
409 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
411 if (cr->nc.rank_intra == 0)
413 /* Use two step summing. */
414 MPI_Reduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,0,
416 /* Sum the roots of the internal (intra) buffers. */
417 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,
422 /* This is here because of the silly MPI specification
423 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
424 MPI_Reduce(r,NULL,nr,MPI_DOUBLE,MPI_SUM,0,cr->nc.comm_intra);
426 MPI_Bcast(r,nr,MPI_DOUBLE,0,cr->nc.comm_intra);
430 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,
431 cr->mpi_comm_mygroup);
436 if (nr > cr->mpb->dbuf_alloc) {
437 cr->mpb->dbuf_alloc = nr;
438 srenew(cr->mpb->dbuf,cr->mpb->dbuf_alloc);
441 /* Use two step summing */
442 MPI_Allreduce(r,cr->mpb->dbuf,nr,MPI_DOUBLE,MPI_SUM,cr->nc.comm_intra);
443 if (cr->nc.rank_intra == 0) {
444 /* Sum with the buffers reversed */
445 MPI_Allreduce(cr->mpb->dbuf,r,nr,MPI_DOUBLE,MPI_SUM,
448 MPI_Bcast(r,nr,MPI_DOUBLE,0,cr->nc.comm_intra);
450 MPI_Allreduce(r,cr->mpb->dbuf,nr,MPI_DOUBLE,MPI_SUM,
451 cr->mpi_comm_mygroup);
453 r[i] = cr->mpb->dbuf[i];
459 void gmx_sumf(int nr,float r[],const t_commrec *cr)
462 gmx_call("gmx_sumf");
464 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
466 /* Use two step summing. */
467 if (cr->nc.rank_intra == 0)
469 MPI_Reduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,0,
471 /* Sum the roots of the internal (intra) buffers */
472 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,
477 /* This is here because of the silly MPI specification
478 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
479 MPI_Reduce(r,NULL,nr,MPI_FLOAT,MPI_SUM,0,cr->nc.comm_intra);
481 MPI_Bcast(r,nr,MPI_FLOAT,0,cr->nc.comm_intra);
485 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,cr->mpi_comm_mygroup);
490 if (nr > cr->mpb->fbuf_alloc) {
491 cr->mpb->fbuf_alloc = nr;
492 srenew(cr->mpb->fbuf,cr->mpb->fbuf_alloc);
495 /* Use two step summing */
496 MPI_Allreduce(r,cr->mpb->fbuf,nr,MPI_FLOAT,MPI_SUM,cr->nc.comm_intra);
497 if (cr->nc.rank_intra == 0) {
498 /* Sum with the buffers reversed */
499 MPI_Allreduce(cr->mpb->fbuf,r,nr,MPI_FLOAT,MPI_SUM,
502 MPI_Bcast(r,nr,MPI_FLOAT,0,cr->nc.comm_intra);
504 MPI_Allreduce(r,cr->mpb->fbuf,nr,MPI_FLOAT,MPI_SUM,
505 cr->mpi_comm_mygroup);
507 r[i] = cr->mpb->fbuf[i];
513 void gmx_sumi(int nr,int r[],const t_commrec *cr)
516 gmx_call("gmx_sumi");
518 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
520 /* Use two step summing */
521 if (cr->nc.rank_intra == 0)
523 MPI_Reduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,0,cr->nc.comm_intra);
524 /* Sum with the buffers reversed */
525 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,cr->nc.comm_inter);
529 /* This is here because of the silly MPI specification
530 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
531 MPI_Reduce(r,NULL,nr,MPI_INT,MPI_SUM,0,cr->nc.comm_intra);
533 MPI_Bcast(r,nr,MPI_INT,0,cr->nc.comm_intra);
537 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,cr->mpi_comm_mygroup);
542 if (nr > cr->mpb->ibuf_alloc) {
543 cr->mpb->ibuf_alloc = nr;
544 srenew(cr->mpb->ibuf,cr->mpb->ibuf_alloc);
547 /* Use two step summing */
548 MPI_Allreduce(r,cr->mpb->ibuf,nr,MPI_INT,MPI_SUM,cr->nc.comm_intra);
549 if (cr->nc.rank_intra == 0) {
550 /* Sum with the buffers reversed */
551 MPI_Allreduce(cr->mpb->ibuf,r,nr,MPI_INT,MPI_SUM,cr->nc.comm_inter);
553 MPI_Bcast(r,nr,MPI_INT,0,cr->nc.comm_intra);
555 MPI_Allreduce(r,cr->mpb->ibuf,nr,MPI_INT,MPI_SUM,cr->mpi_comm_mygroup);
557 r[i] = cr->mpb->ibuf[i];
563 void gmx_sumli(int nr,gmx_large_int_t r[],const t_commrec *cr)
566 gmx_call("gmx_sumli");
568 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
570 /* Use two step summing */
571 if (cr->nc.rank_intra == 0)
573 MPI_Reduce(MPI_IN_PLACE,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,0,
575 /* Sum with the buffers reversed */
576 MPI_Allreduce(MPI_IN_PLACE,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,
581 /* This is here because of the silly MPI specification
582 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
583 MPI_Reduce(r,NULL,nr,GMX_MPI_LARGE_INT,MPI_SUM,0,cr->nc.comm_intra);
585 MPI_Bcast(r,nr,GMX_MPI_LARGE_INT,0,cr->nc.comm_intra);
589 MPI_Allreduce(MPI_IN_PLACE,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,cr->mpi_comm_mygroup);
594 if (nr > cr->mpb->libuf_alloc) {
595 cr->mpb->libuf_alloc = nr;
596 srenew(cr->mpb->libuf,cr->mpb->libuf_alloc);
599 /* Use two step summing */
600 MPI_Allreduce(r,cr->mpb->libuf,nr,GMX_MPI_LARGE_INT,MPI_SUM,
602 if (cr->nc.rank_intra == 0) {
603 /* Sum with the buffers reversed */
604 MPI_Allreduce(cr->mpb->libuf,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,
607 MPI_Bcast(r,nr,GMX_MPI_LARGE_INT,0,cr->nc.comm_intra);
609 MPI_Allreduce(r,cr->mpb->libuf,nr,GMX_MPI_LARGE_INT,MPI_SUM,
610 cr->mpi_comm_mygroup);
612 r[i] = cr->mpb->libuf[i];
621 void gmx_sumd_comm(int nr,double r[],MPI_Comm mpi_comm)
623 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
624 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,mpi_comm);
626 /* this function is only used in code that is not performance critical,
627 (during setup, when comm_rec is not the appropriate communication
628 structure), so this isn't as bad as it looks. */
633 MPI_Allreduce(r,buf,nr,MPI_DOUBLE,MPI_SUM,mpi_comm);
642 void gmx_sumf_comm(int nr,float r[],MPI_Comm mpi_comm)
644 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
645 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,mpi_comm);
647 /* this function is only used in code that is not performance critical,
648 (during setup, when comm_rec is not the appropriate communication
649 structure), so this isn't as bad as it looks. */
654 MPI_Allreduce(r,buf,nr,MPI_FLOAT,MPI_SUM,mpi_comm);
662 void gmx_sumd_sim(int nr,double r[],const gmx_multisim_t *ms)
665 gmx_call("gmx_sumd_sim");
667 gmx_sumd_comm(nr,r,ms->mpi_comm_masters);
671 void gmx_sumf_sim(int nr,float r[],const gmx_multisim_t *ms)
674 gmx_call("gmx_sumf_sim");
676 gmx_sumf_comm(nr,r,ms->mpi_comm_masters);
680 void gmx_sumi_sim(int nr,int r[], const gmx_multisim_t *ms)
683 gmx_call("gmx_sumi_sim");
685 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
686 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,ms->mpi_comm_masters);
688 /* this is thread-unsafe, but it will do for now: */
691 if (nr > ms->mpb->ibuf_alloc) {
692 ms->mpb->ibuf_alloc = nr;
693 srenew(ms->mpb->ibuf,ms->mpb->ibuf_alloc);
695 MPI_Allreduce(r,ms->mpb->ibuf,nr,MPI_INT,MPI_SUM,ms->mpi_comm_masters);
697 r[i] = ms->mpb->ibuf[i];
702 void gmx_sumli_sim(int nr,gmx_large_int_t r[], const gmx_multisim_t *ms)
705 gmx_call("gmx_sumli_sim");
707 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
708 MPI_Allreduce(MPI_IN_PLACE,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,
709 ms->mpi_comm_masters);
711 /* this is thread-unsafe, but it will do for now: */
714 if (nr > ms->mpb->libuf_alloc) {
715 ms->mpb->libuf_alloc = nr;
716 srenew(ms->mpb->libuf,ms->mpb->libuf_alloc);
718 MPI_Allreduce(r,ms->mpb->libuf,nr,GMX_MPI_LARGE_INT,MPI_SUM,
719 ms->mpi_comm_masters);
721 r[i] = ms->mpb->libuf[i];
727 void gmx_finalize_par(void)
730 /* Compiled without MPI, no MPI finalizing needed */
733 int initialized,finalized;
736 MPI_Initialized(&initialized);
741 /* just as a check; we don't want to finalize twice */
742 MPI_Finalized(&finalized);
748 /* We sync the processes here to try to avoid problems
749 * with buggy MPI implementations that could cause
750 * unfinished processes to terminate.
752 MPI_Barrier(MPI_COMM_WORLD);
755 if (DOMAINDECOMP(cr)) {
756 if (cr->npmenodes > 0 || cr->dd->bCartesian)
757 MPI_Comm_free(&cr->mpi_comm_mygroup);
758 if (cr->dd->bCartesian)
759 MPI_Comm_free(&cr->mpi_comm_mysim);
763 /* Apparently certain mpich implementations cause problems
764 * with MPI_Finalize. In that case comment out MPI_Finalize.
767 fprintf(debug,"Will call MPI_Finalize now\n");
769 ret = MPI_Finalize();
771 fprintf(debug,"Return code from MPI_Finalize = %d\n",ret);