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_intra_counters(t_commrec *cr)
264 /* counters for PP+PME and PP-only processes on my node */
265 int nnodes, nnodes_pp, id_mynode=-1, id_mynode_group=-1, nproc_mynode, nproc_mynode_pp;
266 #if defined GMX_MPI && !defined GMX_THREAD_MPI
267 int i, mynum, *num, *num_s, *num_pp, *num_pp_s;
271 nnodes_pp = nnodes - cr->npmenodes;
273 #if defined GMX_MPI && !defined GMX_THREAD_MPI
274 /* We have MPI and can expect to have different compute nodes */
275 mynum = gmx_hostname_num();
277 /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
280 snew(num_pp, nnodes_pp);
281 snew(num_pp_s, nnodes_pp);
283 num_s[cr->sim_nodeid] = mynum;
284 if (cr->duty & DUTY_PP)
286 num_pp_s[cr->nodeid] = mynum;
289 MPI_Allreduce(num_s, num, nnodes, MPI_INT, MPI_SUM, cr->mpi_comm_mysim);
290 MPI_Allreduce(num_pp_s, num_pp, nnodes_pp, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
296 for(i=0; i<nnodes; i++)
301 if (i < cr->sim_nodeid)
311 for(i=0; i<nnodes_pp; i++)
313 if (num_pp[i] == mynum)
323 /* Serial or thread-MPI code, we are running within a node */
324 id_mynode = cr->sim_nodeid;
325 id_mynode_group = cr->nodeid;
326 nproc_mynode = cr->nnodes;
327 nproc_mynode_pp = cr->nnodes - cr->npmenodes;
333 if (cr->duty & DUTY_PP && cr->duty & DUTY_PME)
335 sprintf(sbuf, "PP+PME");
339 sprintf(sbuf, "%s", cr->duty & DUTY_PP ? "PP" : "PME");
341 fprintf(debug, "On %3s node %d: nodeid_intra=%d, nodeid_group_intra=%d, "
342 "nnodes_intra=%d, nnodes_pp_intra=%d\n", sbuf, cr->sim_nodeid,
343 id_mynode, id_mynode_group, nproc_mynode, nproc_mynode_pp);
346 cr->nodeid_intra = id_mynode;
347 cr->nodeid_group_intra = id_mynode_group;
348 cr->nnodes_intra = nproc_mynode;
349 cr->nnodes_pp_intra = nproc_mynode_pp;
353 void gmx_barrier(const t_commrec *cr)
356 gmx_call("gmx_barrier");
358 MPI_Barrier(cr->mpi_comm_mygroup);
362 void gmx_abort(int noderank,int nnodes,int errorno)
365 gmx_call("gmx_abort");
367 #ifdef GMX_THREAD_MPI
368 fprintf(stderr,"Halting program %s\n",ShortProgram());
374 fprintf(stderr,"Halting parallel program %s on CPU %d out of %d\n",
375 ShortProgram(),noderank,nnodes);
379 fprintf(stderr,"Halting program %s\n",ShortProgram());
383 MPI_Abort(MPI_COMM_WORLD,errorno);
389 void gmx_bcast(int nbytes,void *b,const t_commrec *cr)
392 gmx_call("gmx_bast");
394 MPI_Bcast(b,nbytes,MPI_BYTE,MASTERRANK(cr),cr->mpi_comm_mygroup);
398 void gmx_bcast_sim(int nbytes,void *b,const t_commrec *cr)
401 gmx_call("gmx_bast");
403 MPI_Bcast(b,nbytes,MPI_BYTE,MASTERRANK(cr),cr->mpi_comm_mysim);
407 void gmx_sumd(int nr,double r[],const t_commrec *cr)
410 gmx_call("gmx_sumd");
412 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
414 if (cr->nc.rank_intra == 0)
416 /* Use two step summing. */
417 MPI_Reduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,0,
419 /* Sum the roots of the internal (intra) buffers. */
420 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,
425 /* This is here because of the silly MPI specification
426 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
427 MPI_Reduce(r,NULL,nr,MPI_DOUBLE,MPI_SUM,0,cr->nc.comm_intra);
429 MPI_Bcast(r,nr,MPI_DOUBLE,0,cr->nc.comm_intra);
433 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,
434 cr->mpi_comm_mygroup);
439 if (nr > cr->mpb->dbuf_alloc) {
440 cr->mpb->dbuf_alloc = nr;
441 srenew(cr->mpb->dbuf,cr->mpb->dbuf_alloc);
444 /* Use two step summing */
445 MPI_Allreduce(r,cr->mpb->dbuf,nr,MPI_DOUBLE,MPI_SUM,cr->nc.comm_intra);
446 if (cr->nc.rank_intra == 0) {
447 /* Sum with the buffers reversed */
448 MPI_Allreduce(cr->mpb->dbuf,r,nr,MPI_DOUBLE,MPI_SUM,
451 MPI_Bcast(r,nr,MPI_DOUBLE,0,cr->nc.comm_intra);
453 MPI_Allreduce(r,cr->mpb->dbuf,nr,MPI_DOUBLE,MPI_SUM,
454 cr->mpi_comm_mygroup);
456 r[i] = cr->mpb->dbuf[i];
462 void gmx_sumf(int nr,float r[],const t_commrec *cr)
465 gmx_call("gmx_sumf");
467 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
469 /* Use two step summing. */
470 if (cr->nc.rank_intra == 0)
472 MPI_Reduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,0,
474 /* Sum the roots of the internal (intra) buffers */
475 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,
480 /* This is here because of the silly MPI specification
481 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
482 MPI_Reduce(r,NULL,nr,MPI_FLOAT,MPI_SUM,0,cr->nc.comm_intra);
484 MPI_Bcast(r,nr,MPI_FLOAT,0,cr->nc.comm_intra);
488 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,cr->mpi_comm_mygroup);
493 if (nr > cr->mpb->fbuf_alloc) {
494 cr->mpb->fbuf_alloc = nr;
495 srenew(cr->mpb->fbuf,cr->mpb->fbuf_alloc);
498 /* Use two step summing */
499 MPI_Allreduce(r,cr->mpb->fbuf,nr,MPI_FLOAT,MPI_SUM,cr->nc.comm_intra);
500 if (cr->nc.rank_intra == 0) {
501 /* Sum with the buffers reversed */
502 MPI_Allreduce(cr->mpb->fbuf,r,nr,MPI_FLOAT,MPI_SUM,
505 MPI_Bcast(r,nr,MPI_FLOAT,0,cr->nc.comm_intra);
507 MPI_Allreduce(r,cr->mpb->fbuf,nr,MPI_FLOAT,MPI_SUM,
508 cr->mpi_comm_mygroup);
510 r[i] = cr->mpb->fbuf[i];
516 void gmx_sumi(int nr,int r[],const t_commrec *cr)
519 gmx_call("gmx_sumi");
521 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
523 /* Use two step summing */
524 if (cr->nc.rank_intra == 0)
526 MPI_Reduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,0,cr->nc.comm_intra);
527 /* Sum with the buffers reversed */
528 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,cr->nc.comm_inter);
532 /* This is here because of the silly MPI specification
533 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
534 MPI_Reduce(r,NULL,nr,MPI_INT,MPI_SUM,0,cr->nc.comm_intra);
536 MPI_Bcast(r,nr,MPI_INT,0,cr->nc.comm_intra);
540 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,cr->mpi_comm_mygroup);
545 if (nr > cr->mpb->ibuf_alloc) {
546 cr->mpb->ibuf_alloc = nr;
547 srenew(cr->mpb->ibuf,cr->mpb->ibuf_alloc);
550 /* Use two step summing */
551 MPI_Allreduce(r,cr->mpb->ibuf,nr,MPI_INT,MPI_SUM,cr->nc.comm_intra);
552 if (cr->nc.rank_intra == 0) {
553 /* Sum with the buffers reversed */
554 MPI_Allreduce(cr->mpb->ibuf,r,nr,MPI_INT,MPI_SUM,cr->nc.comm_inter);
556 MPI_Bcast(r,nr,MPI_INT,0,cr->nc.comm_intra);
558 MPI_Allreduce(r,cr->mpb->ibuf,nr,MPI_INT,MPI_SUM,cr->mpi_comm_mygroup);
560 r[i] = cr->mpb->ibuf[i];
566 void gmx_sumli(int nr,gmx_large_int_t r[],const t_commrec *cr)
569 gmx_call("gmx_sumli");
571 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
573 /* Use two step summing */
574 if (cr->nc.rank_intra == 0)
576 MPI_Reduce(MPI_IN_PLACE,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,0,
578 /* Sum with the buffers reversed */
579 MPI_Allreduce(MPI_IN_PLACE,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,
584 /* This is here because of the silly MPI specification
585 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
586 MPI_Reduce(r,NULL,nr,GMX_MPI_LARGE_INT,MPI_SUM,0,cr->nc.comm_intra);
588 MPI_Bcast(r,nr,GMX_MPI_LARGE_INT,0,cr->nc.comm_intra);
592 MPI_Allreduce(MPI_IN_PLACE,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,cr->mpi_comm_mygroup);
597 if (nr > cr->mpb->libuf_alloc) {
598 cr->mpb->libuf_alloc = nr;
599 srenew(cr->mpb->libuf,cr->mpb->libuf_alloc);
602 /* Use two step summing */
603 MPI_Allreduce(r,cr->mpb->libuf,nr,GMX_MPI_LARGE_INT,MPI_SUM,
605 if (cr->nc.rank_intra == 0) {
606 /* Sum with the buffers reversed */
607 MPI_Allreduce(cr->mpb->libuf,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,
610 MPI_Bcast(r,nr,GMX_MPI_LARGE_INT,0,cr->nc.comm_intra);
612 MPI_Allreduce(r,cr->mpb->libuf,nr,GMX_MPI_LARGE_INT,MPI_SUM,
613 cr->mpi_comm_mygroup);
615 r[i] = cr->mpb->libuf[i];
624 void gmx_sumd_comm(int nr,double r[],MPI_Comm mpi_comm)
626 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
627 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,mpi_comm);
629 /* this function is only used in code that is not performance critical,
630 (during setup, when comm_rec is not the appropriate communication
631 structure), so this isn't as bad as it looks. */
636 MPI_Allreduce(r,buf,nr,MPI_DOUBLE,MPI_SUM,mpi_comm);
645 void gmx_sumf_comm(int nr,float r[],MPI_Comm mpi_comm)
647 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
648 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,mpi_comm);
650 /* this function is only used in code that is not performance critical,
651 (during setup, when comm_rec is not the appropriate communication
652 structure), so this isn't as bad as it looks. */
657 MPI_Allreduce(r,buf,nr,MPI_FLOAT,MPI_SUM,mpi_comm);
665 void gmx_sumd_sim(int nr,double r[],const gmx_multisim_t *ms)
668 gmx_call("gmx_sumd_sim");
670 gmx_sumd_comm(nr,r,ms->mpi_comm_masters);
674 void gmx_sumf_sim(int nr,float r[],const gmx_multisim_t *ms)
677 gmx_call("gmx_sumf_sim");
679 gmx_sumf_comm(nr,r,ms->mpi_comm_masters);
683 void gmx_sumi_sim(int nr,int r[], const gmx_multisim_t *ms)
686 gmx_call("gmx_sumi_sim");
688 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
689 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,ms->mpi_comm_masters);
691 /* this is thread-unsafe, but it will do for now: */
694 if (nr > ms->mpb->ibuf_alloc) {
695 ms->mpb->ibuf_alloc = nr;
696 srenew(ms->mpb->ibuf,ms->mpb->ibuf_alloc);
698 MPI_Allreduce(r,ms->mpb->ibuf,nr,MPI_INT,MPI_SUM,ms->mpi_comm_masters);
700 r[i] = ms->mpb->ibuf[i];
705 void gmx_sumli_sim(int nr,gmx_large_int_t r[], const gmx_multisim_t *ms)
708 gmx_call("gmx_sumli_sim");
710 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
711 MPI_Allreduce(MPI_IN_PLACE,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,
712 ms->mpi_comm_masters);
714 /* this is thread-unsafe, but it will do for now: */
717 if (nr > ms->mpb->libuf_alloc) {
718 ms->mpb->libuf_alloc = nr;
719 srenew(ms->mpb->libuf,ms->mpb->libuf_alloc);
721 MPI_Allreduce(r,ms->mpb->libuf,nr,GMX_MPI_LARGE_INT,MPI_SUM,
722 ms->mpi_comm_masters);
724 r[i] = ms->mpb->libuf[i];
730 void gmx_finalize_par(void)
733 /* Compiled without MPI, no MPI finalizing needed */
736 int initialized,finalized;
739 MPI_Initialized(&initialized);
744 /* just as a check; we don't want to finalize twice */
745 MPI_Finalized(&finalized);
751 /* We sync the processes here to try to avoid problems
752 * with buggy MPI implementations that could cause
753 * unfinished processes to terminate.
755 MPI_Barrier(MPI_COMM_WORLD);
758 if (DOMAINDECOMP(cr)) {
759 if (cr->npmenodes > 0 || cr->dd->bCartesian)
760 MPI_Comm_free(&cr->mpi_comm_mygroup);
761 if (cr->dd->bCartesian)
762 MPI_Comm_free(&cr->mpi_comm_mysim);
766 /* Apparently certain mpich implementations cause problems
767 * with MPI_Finalize. In that case comment out MPI_Finalize.
770 fprintf(debug,"Will call MPI_Finalize now\n");
772 ret = MPI_Finalize();
774 fprintf(debug,"Return code from MPI_Finalize = %d\n",ret);