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);
131 void gmx_setup_nodecomm(FILE *fplog,t_commrec *cr)
134 int n,rank,resultlen,hostnum,i,j,ng,ni;
136 char mpi_hostname[MPI_MAX_PROCESSOR_NAME],num[MPI_MAX_PROCESSOR_NAME];
139 /* Many MPI implementations do not optimize MPI_Allreduce
140 * (and probably also other global communication calls)
141 * for multi-core nodes connected by a network.
142 * We can optimize such communication by using one MPI call
143 * within each node and one between the nodes.
144 * For MVAPICH2 and Intel MPI this reduces the time for
145 * the global_stat communication by 25%
146 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
147 * B. Hess, November 2007
154 if (getenv("GMX_NO_NODECOMM") == NULL) {
156 MPI_Comm_size(cr->mpi_comm_mygroup,&n);
157 MPI_Comm_rank(cr->mpi_comm_mygroup,&rank);
158 MPI_Get_processor_name(mpi_hostname,&resultlen);
159 /* This procedure can only differentiate nodes with host names
160 * that end on unique numbers.
164 /* Only parse the host name up to the first dot */
165 while(i < resultlen && mpi_hostname[i] != '.') {
166 if (isdigit(mpi_hostname[i])) {
167 num[j++] = mpi_hostname[i];
175 /* Use only the last 9 decimals, so we don't overflow an int */
176 hostnum = strtol(num + max(0,j-9), NULL, 10);
181 "In gmx_setup_nodecomm: splitting communicator of size %d\n",
183 fprintf(debug,"In gmx_setup_nodecomm: hostname '%s', hostnum %d\n",
184 mpi_hostname,hostnum);
187 /* The intra-node communicator, split on node number */
188 MPI_Comm_split(cr->mpi_comm_mygroup,hostnum,rank,&nc->comm_intra);
189 MPI_Comm_rank(nc->comm_intra,&nc->rank_intra);
191 fprintf(debug,"In gmx_setup_nodecomm: node rank %d rank_intra %d\n",
192 rank,nc->rank_intra);
194 /* The inter-node communicator, split on rank_intra.
195 * We actually only need the one for rank=0,
196 * but it is easier to create them all.
198 MPI_Comm_split(cr->mpi_comm_mygroup,nc->rank_intra,rank,&nc->comm_inter);
199 /* Check if this really created two step communication */
200 MPI_Comm_size(nc->comm_inter,&ng);
201 MPI_Comm_size(nc->comm_intra,&ni);
203 fprintf(debug,"In gmx_setup_nodecomm: groups %d, my group size %d\n",
206 if ((ng > 1 && ng < n) || (ni > 1 && ni < n)) {
209 fprintf(fplog,"Using two step summing over %d groups of on average %.1f processes\n\n",ng,(real)n/(real)ng);
210 if (nc->rank_intra > 0)
211 MPI_Comm_free(&nc->comm_inter);
213 /* One group or all processes in a separate group, use normal summing */
214 MPI_Comm_free(&nc->comm_inter);
215 MPI_Comm_free(&nc->comm_intra);
222 void gmx_barrier(const t_commrec *cr)
225 gmx_call("gmx_barrier");
227 MPI_Barrier(cr->mpi_comm_mygroup);
231 void gmx_abort(int noderank,int nnodes,int errorno)
234 gmx_call("gmx_abort");
237 fprintf(stderr,"Halting program %s\n",ShortProgram());
243 fprintf(stderr,"Halting parallel program %s on CPU %d out of %d\n",
244 ShortProgram(),noderank,nnodes);
248 fprintf(stderr,"Halting program %s\n",ShortProgram());
252 MPI_Abort(MPI_COMM_WORLD,errorno);
258 void gmx_bcast(int nbytes,void *b,const t_commrec *cr)
261 gmx_call("gmx_bast");
263 MPI_Bcast(b,nbytes,MPI_BYTE,MASTERRANK(cr),cr->mpi_comm_mygroup);
267 void gmx_bcast_sim(int nbytes,void *b,const t_commrec *cr)
270 gmx_call("gmx_bast");
272 MPI_Bcast(b,nbytes,MPI_BYTE,MASTERRANK(cr),cr->mpi_comm_mysim);
276 void gmx_sumd(int nr,double r[],const t_commrec *cr)
279 gmx_call("gmx_sumd");
281 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
283 if (cr->nc.rank_intra == 0)
285 /* Use two step summing. */
286 MPI_Reduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,0,
288 /* Sum the roots of the internal (intra) buffers. */
289 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,
294 /* This is here because of the silly MPI specification
295 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
296 MPI_Reduce(r,NULL,nr,MPI_DOUBLE,MPI_SUM,0,cr->nc.comm_intra);
298 MPI_Bcast(r,nr,MPI_DOUBLE,0,cr->nc.comm_intra);
302 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,
303 cr->mpi_comm_mygroup);
308 if (nr > cr->mpb->dbuf_alloc) {
309 cr->mpb->dbuf_alloc = nr;
310 srenew(cr->mpb->dbuf,cr->mpb->dbuf_alloc);
313 /* Use two step summing */
314 MPI_Allreduce(r,cr->mpb->dbuf,nr,MPI_DOUBLE,MPI_SUM,cr->nc.comm_intra);
315 if (cr->nc.rank_intra == 0) {
316 /* Sum with the buffers reversed */
317 MPI_Allreduce(cr->mpb->dbuf,r,nr,MPI_DOUBLE,MPI_SUM,
320 MPI_Bcast(r,nr,MPI_DOUBLE,0,cr->nc.comm_intra);
322 MPI_Allreduce(r,cr->mpb->dbuf,nr,MPI_DOUBLE,MPI_SUM,
323 cr->mpi_comm_mygroup);
325 r[i] = cr->mpb->dbuf[i];
331 void gmx_sumf(int nr,float r[],const t_commrec *cr)
334 gmx_call("gmx_sumf");
336 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
338 /* Use two step summing. */
339 if (cr->nc.rank_intra == 0)
341 MPI_Reduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,0,
343 /* Sum the roots of the internal (intra) buffers */
344 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,
349 /* This is here because of the silly MPI specification
350 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
351 MPI_Reduce(r,NULL,nr,MPI_FLOAT,MPI_SUM,0,cr->nc.comm_intra);
353 MPI_Bcast(r,nr,MPI_FLOAT,0,cr->nc.comm_intra);
357 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,cr->mpi_comm_mygroup);
362 if (nr > cr->mpb->fbuf_alloc) {
363 cr->mpb->fbuf_alloc = nr;
364 srenew(cr->mpb->fbuf,cr->mpb->fbuf_alloc);
367 /* Use two step summing */
368 MPI_Allreduce(r,cr->mpb->fbuf,nr,MPI_FLOAT,MPI_SUM,cr->nc.comm_intra);
369 if (cr->nc.rank_intra == 0) {
370 /* Sum with the buffers reversed */
371 MPI_Allreduce(cr->mpb->fbuf,r,nr,MPI_FLOAT,MPI_SUM,
374 MPI_Bcast(r,nr,MPI_FLOAT,0,cr->nc.comm_intra);
376 MPI_Allreduce(r,cr->mpb->fbuf,nr,MPI_FLOAT,MPI_SUM,
377 cr->mpi_comm_mygroup);
379 r[i] = cr->mpb->fbuf[i];
385 void gmx_sumi(int nr,int r[],const t_commrec *cr)
388 gmx_call("gmx_sumi");
390 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
392 /* Use two step summing */
393 if (cr->nc.rank_intra == 0)
395 MPI_Reduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,0,cr->nc.comm_intra);
396 /* Sum with the buffers reversed */
397 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,cr->nc.comm_inter);
401 /* This is here because of the silly MPI specification
402 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
403 MPI_Reduce(r,NULL,nr,MPI_INT,MPI_SUM,0,cr->nc.comm_intra);
405 MPI_Bcast(r,nr,MPI_INT,0,cr->nc.comm_intra);
409 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,cr->mpi_comm_mygroup);
414 if (nr > cr->mpb->ibuf_alloc) {
415 cr->mpb->ibuf_alloc = nr;
416 srenew(cr->mpb->ibuf,cr->mpb->ibuf_alloc);
419 /* Use two step summing */
420 MPI_Allreduce(r,cr->mpb->ibuf,nr,MPI_INT,MPI_SUM,cr->nc.comm_intra);
421 if (cr->nc.rank_intra == 0) {
422 /* Sum with the buffers reversed */
423 MPI_Allreduce(cr->mpb->ibuf,r,nr,MPI_INT,MPI_SUM,cr->nc.comm_inter);
425 MPI_Bcast(r,nr,MPI_INT,0,cr->nc.comm_intra);
427 MPI_Allreduce(r,cr->mpb->ibuf,nr,MPI_INT,MPI_SUM,cr->mpi_comm_mygroup);
429 r[i] = cr->mpb->ibuf[i];
435 void gmx_sumli(int nr,gmx_large_int_t r[],const t_commrec *cr)
438 gmx_call("gmx_sumli");
440 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
442 /* Use two step summing */
443 if (cr->nc.rank_intra == 0)
445 MPI_Reduce(MPI_IN_PLACE,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,0,
447 /* Sum with the buffers reversed */
448 MPI_Allreduce(MPI_IN_PLACE,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,
453 /* This is here because of the silly MPI specification
454 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
455 MPI_Reduce(r,NULL,nr,GMX_MPI_LARGE_INT,MPI_SUM,0,cr->nc.comm_intra);
457 MPI_Bcast(r,nr,GMX_MPI_LARGE_INT,0,cr->nc.comm_intra);
461 MPI_Allreduce(MPI_IN_PLACE,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,cr->mpi_comm_mygroup);
466 if (nr > cr->mpb->ibuf_alloc) {
467 cr->mpb->ibuf_alloc = nr;
468 srenew(cr->mpb->ibuf,cr->mpb->ibuf_alloc);
471 /* Use two step summing */
472 MPI_Allreduce(r,cr->mpb->ibuf,nr,GMX_MPI_LARGE_INT,MPI_SUM,
474 if (cr->nc.rank_intra == 0) {
475 /* Sum with the buffers reversed */
476 MPI_Allreduce(cr->mpb->ibuf,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,
479 MPI_Bcast(r,nr,GMX_MPI_LARGE_INT,0,cr->nc.comm_intra);
481 MPI_Allreduce(r,cr->mpb->ibuf,nr,GMX_MPI_LARGE_INT,MPI_SUM,
482 cr->mpi_comm_mygroup);
484 r[i] = cr->mpb->ibuf[i];
493 void gmx_sumd_comm(int nr,double r[],MPI_Comm mpi_comm)
495 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
496 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,mpi_comm);
498 /* this function is only used in code that is not performance critical,
499 (during setup, when comm_rec is not the appropriate communication
500 structure), so this isn't as bad as it looks. */
505 MPI_Allreduce(r,buf,nr,MPI_DOUBLE,MPI_SUM,mpi_comm);
514 void gmx_sumf_comm(int nr,float r[],MPI_Comm mpi_comm)
516 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
517 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,mpi_comm);
519 /* this function is only used in code that is not performance critical,
520 (during setup, when comm_rec is not the appropriate communication
521 structure), so this isn't as bad as it looks. */
526 MPI_Allreduce(r,buf,nr,MPI_FLOAT,MPI_SUM,mpi_comm);
534 void gmx_sumd_sim(int nr,double r[],const gmx_multisim_t *ms)
537 gmx_call("gmx_sumd_sim");
539 gmx_sumd_comm(nr,r,ms->mpi_comm_masters);
543 void gmx_sumf_sim(int nr,float r[],const gmx_multisim_t *ms)
546 gmx_call("gmx_sumf_sim");
548 gmx_sumf_comm(nr,r,ms->mpi_comm_masters);
552 void gmx_sumi_sim(int nr,int r[], const gmx_multisim_t *ms)
555 gmx_call("gmx_sumi_sim");
557 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
558 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,ms->mpi_comm_masters);
560 /* this is thread-unsafe, but it will do for now: */
563 if (nr > ms->mpb->ibuf_alloc) {
564 ms->mpb->ibuf_alloc = nr;
565 srenew(ms->mpb->ibuf,ms->mpb->ibuf_alloc);
567 MPI_Allreduce(r,ms->mpb->ibuf,nr,MPI_INT,MPI_SUM,ms->mpi_comm_masters);
569 r[i] = ms->mpb->ibuf[i];
574 void gmx_sumli_sim(int nr,gmx_large_int_t r[], const gmx_multisim_t *ms)
577 gmx_call("gmx_sumli_sim");
579 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
580 MPI_Allreduce(MPI_IN_PLACE,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,
581 ms->mpi_comm_masters);
583 /* this is thread-unsafe, but it will do for now: */
586 if (nr > ms->mpb->ibuf_alloc) {
587 ms->mpb->ibuf_alloc = nr;
588 srenew(ms->mpb->ibuf,ms->mpb->ibuf_alloc);
590 MPI_Allreduce(r,ms->mpb->ibuf,nr,GMX_MPI_LARGE_INT,MPI_SUM,
591 ms->mpi_comm_masters);
593 r[i] = ms->mpb->ibuf[i];
599 void gmx_finalize(void)
602 gmx_call("gmx_finalize");
606 /* just as a check; we don't want to finalize twice */
608 MPI_Finalized(&finalized);
612 /* We sync the processes here to try to avoid problems
613 * with buggy MPI implementations that could cause
614 * unfinished processes to terminate.
616 MPI_Barrier(MPI_COMM_WORLD);
619 if (DOMAINDECOMP(cr)) {
620 if (cr->npmenodes > 0 || cr->dd->bCartesian)
621 MPI_Comm_free(&cr->mpi_comm_mygroup);
622 if (cr->dd->bCartesian)
623 MPI_Comm_free(&cr->mpi_comm_mysim);
627 /* Apparently certain mpich implementations cause problems
628 * with MPI_Finalize. In that case comment out MPI_Finalize.
631 fprintf(debug,"Will call MPI_Finalize now\n");
633 ret = MPI_Finalize();
635 fprintf(debug,"Return code from MPI_Finalize = %d\n",ret);