2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
43 #include "gmx_fatal.h"
60 #include "mpelogging.h"
62 /* The source code in this file should be thread-safe.
63 Please keep it that way. */
65 gmx_bool gmx_mpi_initialized(void)
77 int gmx_setup(int *argc,char **argv,int *nnodes)
80 gmx_call("gmx_setup");
84 int resultlen; /* actual length of node name */
88 char mpi_hostname[MPI_MAX_PROCESSOR_NAME];
90 /* Call the MPI routines */
93 (void) fah_MPI_Init(argc,&argv);
95 (void) MPI_Init(argc,&argv);
98 (void) MPI_Comm_size( MPI_COMM_WORLD, &mpi_num_nodes );
99 (void) MPI_Comm_rank( MPI_COMM_WORLD, &mpi_my_rank );
100 (void) MPI_Get_processor_name( mpi_hostname, &resultlen );
104 /* MPE logging routines. Get event IDs from MPE: */
106 ev_timestep1 = MPE_Log_get_event_number( );
107 ev_timestep2 = MPE_Log_get_event_number( );
108 ev_force_start = MPE_Log_get_event_number( );
109 ev_force_finish = MPE_Log_get_event_number( );
110 ev_do_fnbf_start = MPE_Log_get_event_number( );
111 ev_do_fnbf_finish = MPE_Log_get_event_number( );
112 ev_ns_start = MPE_Log_get_event_number( );
113 ev_ns_finish = MPE_Log_get_event_number( );
114 ev_calc_bonds_start = MPE_Log_get_event_number( );
115 ev_calc_bonds_finish = MPE_Log_get_event_number( );
116 ev_global_stat_start = MPE_Log_get_event_number( );
117 ev_global_stat_finish = MPE_Log_get_event_number( );
118 ev_virial_start = MPE_Log_get_event_number( );
119 ev_virial_finish = MPE_Log_get_event_number( );
121 /* Shift related events */
122 ev_shift_start = MPE_Log_get_event_number( );
123 ev_shift_finish = MPE_Log_get_event_number( );
124 ev_unshift_start = MPE_Log_get_event_number( );
125 ev_unshift_finish = MPE_Log_get_event_number( );
126 ev_mk_mshift_start = MPE_Log_get_event_number( );
127 ev_mk_mshift_finish = MPE_Log_get_event_number( );
129 /* PME related events */
130 ev_pme_start = MPE_Log_get_event_number( );
131 ev_pme_finish = MPE_Log_get_event_number( );
132 ev_spread_on_grid_start = MPE_Log_get_event_number( );
133 ev_spread_on_grid_finish = MPE_Log_get_event_number( );
134 ev_sum_qgrid_start = MPE_Log_get_event_number( );
135 ev_sum_qgrid_finish = MPE_Log_get_event_number( );
136 ev_gmxfft3d_start = MPE_Log_get_event_number( );
137 ev_gmxfft3d_finish = MPE_Log_get_event_number( );
138 ev_solve_pme_start = MPE_Log_get_event_number( );
139 ev_solve_pme_finish = MPE_Log_get_event_number( );
140 ev_gather_f_bsplines_start = MPE_Log_get_event_number( );
141 ev_gather_f_bsplines_finish= MPE_Log_get_event_number( );
142 ev_reduce_start = MPE_Log_get_event_number( );
143 ev_reduce_finish = MPE_Log_get_event_number( );
144 ev_rscatter_start = MPE_Log_get_event_number( );
145 ev_rscatter_finish = MPE_Log_get_event_number( );
146 ev_alltoall_start = MPE_Log_get_event_number( );
147 ev_alltoall_finish = MPE_Log_get_event_number( );
148 ev_pmeredist_start = MPE_Log_get_event_number( );
149 ev_pmeredist_finish = MPE_Log_get_event_number( );
150 ev_init_pme_start = MPE_Log_get_event_number( );
151 ev_init_pme_finish = MPE_Log_get_event_number( );
152 ev_send_coordinates_start = MPE_Log_get_event_number( );
153 ev_send_coordinates_finish = MPE_Log_get_event_number( );
154 ev_update_fr_start = MPE_Log_get_event_number( );
155 ev_update_fr_finish = MPE_Log_get_event_number( );
156 ev_clear_rvecs_start = MPE_Log_get_event_number( );
157 ev_clear_rvecs_finish = MPE_Log_get_event_number( );
158 ev_update_start = MPE_Log_get_event_number( );
159 ev_update_finish = MPE_Log_get_event_number( );
160 ev_output_start = MPE_Log_get_event_number( );
161 ev_output_finish = MPE_Log_get_event_number( );
162 ev_sum_lrforces_start = MPE_Log_get_event_number( );
163 ev_sum_lrforces_finish = MPE_Log_get_event_number( );
164 ev_sort_start = MPE_Log_get_event_number( );
165 ev_sort_finish = MPE_Log_get_event_number( );
166 ev_sum_qgrid_start = MPE_Log_get_event_number( );
167 ev_sum_qgrid_finish = MPE_Log_get_event_number( );
169 /* Essential dynamics related events */
170 ev_edsam_start = MPE_Log_get_event_number( );
171 ev_edsam_finish = MPE_Log_get_event_number( );
172 ev_get_coords_start = MPE_Log_get_event_number( );
173 ev_get_coords_finish = MPE_Log_get_event_number( );
174 ev_ed_apply_cons_start = MPE_Log_get_event_number( );
175 ev_ed_apply_cons_finish = MPE_Log_get_event_number( );
176 ev_fit_to_reference_start = MPE_Log_get_event_number( );
177 ev_fit_to_reference_finish = MPE_Log_get_event_number( );
179 /* describe events: */
180 if ( mpi_my_rank == 0 )
183 MPE_Describe_state(ev_timestep1, ev_timestep2, "timestep START", "magenta" );
184 MPE_Describe_state(ev_force_start, ev_force_finish, "force", "cornflower blue" );
185 MPE_Describe_state(ev_do_fnbf_start, ev_do_fnbf_finish, "do_fnbf", "navy" );
186 MPE_Describe_state(ev_ns_start, ev_ns_finish, "neighbor search", "tomato" );
187 MPE_Describe_state(ev_calc_bonds_start, ev_calc_bonds_finish, "bonded forces", "slate blue" );
188 MPE_Describe_state(ev_global_stat_start, ev_global_stat_finish, "global stat", "firebrick3");
189 MPE_Describe_state(ev_update_fr_start, ev_update_fr_finish, "update forcerec", "goldenrod");
190 MPE_Describe_state(ev_clear_rvecs_start, ev_clear_rvecs_finish, "clear rvecs", "bisque");
191 MPE_Describe_state(ev_update_start, ev_update_finish, "update", "cornsilk");
192 MPE_Describe_state(ev_output_start, ev_output_finish, "output", "black");
193 MPE_Describe_state(ev_virial_start, ev_virial_finish, "calc_virial", "thistle4");
195 /* PME related events */
196 MPE_Describe_state(ev_pme_start, ev_pme_finish, "doing PME", "grey" );
197 MPE_Describe_state(ev_spread_on_grid_start, ev_spread_on_grid_finish, "spread", "dark orange" );
198 MPE_Describe_state(ev_sum_qgrid_start, ev_sum_qgrid_finish, "sum qgrid", "slate blue");
199 MPE_Describe_state(ev_gmxfft3d_start, ev_gmxfft3d_finish, "fft3d", "snow2" );
200 MPE_Describe_state(ev_solve_pme_start, ev_solve_pme_finish, "solve PME", "indian red" );
201 MPE_Describe_state(ev_gather_f_bsplines_start, ev_gather_f_bsplines_finish, "bsplines", "light sea green" );
202 MPE_Describe_state(ev_reduce_start, ev_reduce_finish, "reduce", "cyan1" );
203 MPE_Describe_state(ev_rscatter_start, ev_rscatter_finish, "rscatter", "cyan3" );
204 MPE_Describe_state(ev_alltoall_start, ev_alltoall_finish, "alltoall", "LightCyan4" );
205 MPE_Describe_state(ev_pmeredist_start, ev_pmeredist_finish, "pmeredist", "thistle" );
206 MPE_Describe_state(ev_init_pme_start, ev_init_pme_finish, "init PME", "snow4");
207 MPE_Describe_state(ev_send_coordinates_start, ev_send_coordinates_finish, "send_coordinates","blue");
208 MPE_Describe_state(ev_sum_lrforces_start, ev_sum_lrforces_finish, "sum_LRforces", "lime green");
209 MPE_Describe_state(ev_sort_start, ev_sort_finish, "sort pme atoms", "brown");
210 MPE_Describe_state(ev_sum_qgrid_start, ev_sum_qgrid_finish, "sum charge grid", "medium orchid");
212 /* Shift related events */
213 MPE_Describe_state(ev_shift_start, ev_shift_finish, "shift", "orange");
214 MPE_Describe_state(ev_unshift_start, ev_unshift_finish, "unshift", "dark orange");
215 MPE_Describe_state(ev_mk_mshift_start, ev_mk_mshift_finish, "mk_mshift", "maroon");
217 /* Essential dynamics related events */
218 MPE_Describe_state(ev_edsam_start, ev_edsam_finish, "EDSAM", "deep sky blue");
219 MPE_Describe_state(ev_get_coords_start, ev_get_coords_finish, "ED get coords", "steel blue");
220 MPE_Describe_state(ev_ed_apply_cons_start, ev_ed_apply_cons_finish, "ED apply constr", "forest green");
221 MPE_Describe_state(ev_fit_to_reference_start, ev_fit_to_reference_finish, "ED fit to ref", "lavender");
230 fprintf(debug,"NNODES=%d, MYRANK=%d, HOSTNAME=%s\n",
231 mpi_num_nodes,mpi_my_rank,mpi_hostname);
235 *nnodes=mpi_num_nodes;
241 int gmx_node_num(void)
247 (void) MPI_Comm_size(MPI_COMM_WORLD, &i);
252 int gmx_node_rank(void)
258 (void) MPI_Comm_rank(MPI_COMM_WORLD, &i);
264 int gmx_hostname_num()
269 #ifdef GMX_THREAD_MPI
270 /* thread-MPI currently puts the thread number in the process name,
271 * we might want to change this, as this is inconsistent with what
272 * most MPI implementations would do when running on a single node.
276 int resultlen,hostnum,i,j;
277 char mpi_hostname[MPI_MAX_PROCESSOR_NAME],hostnum_str[MPI_MAX_PROCESSOR_NAME];
279 MPI_Get_processor_name(mpi_hostname,&resultlen);
280 /* This procedure can only differentiate nodes with host names
281 * that end on unique numbers.
285 /* Only parse the host name up to the first dot */
286 while(i < resultlen && mpi_hostname[i] != '.') {
287 if (isdigit(mpi_hostname[i])) {
288 hostnum_str[j++] = mpi_hostname[i];
292 hostnum_str[j] = '\0';
296 /* Use only the last 9 decimals, so we don't overflow an int */
297 hostnum = strtol(hostnum_str + max(0,j-9), NULL, 10);
301 fprintf(debug,"In gmx_setup_nodecomm: hostname '%s', hostnum %d\n",
302 mpi_hostname,hostnum);
309 void gmx_setup_nodecomm(FILE *fplog,t_commrec *cr)
312 int n,rank,hostnum,ng,ni;
314 /* Many MPI implementations do not optimize MPI_Allreduce
315 * (and probably also other global communication calls)
316 * for multi-core nodes connected by a network.
317 * We can optimize such communication by using one MPI call
318 * within each node and one between the nodes.
319 * For MVAPICH2 and Intel MPI this reduces the time for
320 * the global_stat communication by 25%
321 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
322 * B. Hess, November 2007
328 #ifndef GMX_THREAD_MPI
330 MPI_Comm_size(cr->mpi_comm_mygroup,&n);
331 MPI_Comm_rank(cr->mpi_comm_mygroup,&rank);
333 hostnum = gmx_hostname_num();
337 fprintf(debug,"In gmx_setup_nodecomm: splitting communicator of size %d\n",n);
341 /* The intra-node communicator, split on node number */
342 MPI_Comm_split(cr->mpi_comm_mygroup,hostnum,rank,&nc->comm_intra);
343 MPI_Comm_rank(nc->comm_intra,&nc->rank_intra);
346 fprintf(debug,"In gmx_setup_nodecomm: node rank %d rank_intra %d\n",
347 rank,nc->rank_intra);
349 /* The inter-node communicator, split on rank_intra.
350 * We actually only need the one for rank=0,
351 * but it is easier to create them all.
353 MPI_Comm_split(cr->mpi_comm_mygroup,nc->rank_intra,rank,&nc->comm_inter);
354 /* Check if this really created two step communication */
355 MPI_Comm_size(nc->comm_inter,&ng);
356 MPI_Comm_size(nc->comm_intra,&ni);
359 fprintf(debug,"In gmx_setup_nodecomm: groups %d, my group size %d\n",
363 if (getenv("GMX_NO_NODECOMM") == NULL &&
364 ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
369 fprintf(fplog,"Using two step summing over %d groups of on average %.1f processes\n\n",
370 ng,(real)n/(real)ng);
372 if (nc->rank_intra > 0)
374 MPI_Comm_free(&nc->comm_inter);
379 /* One group or all processes in a separate group, use normal summing */
380 MPI_Comm_free(&nc->comm_inter);
381 MPI_Comm_free(&nc->comm_intra);
384 fprintf(debug,"In gmx_setup_nodecomm: not unsing separate inter- and intra-node communicators.\n");
389 /* tMPI runs only on a single node so just use the nodeid */
390 nc->rank_intra = cr->nodeid;
394 void gmx_init_intranode_counters(t_commrec *cr)
396 /* counters for PP+PME and PP-only processes on my physical node */
397 int nrank_intranode, rank_intranode;
398 int nrank_pp_intranode, rank_pp_intranode;
399 /* thread-MPI is not initialized when not running in parallel */
400 #if defined GMX_MPI && !defined GMX_THREAD_MPI
401 int nrank_world, rank_world;
402 int i, mynum, *num, *num_s, *num_pp, *num_pp_s;
404 MPI_Comm_size(MPI_COMM_WORLD,&nrank_world);
405 MPI_Comm_rank(MPI_COMM_WORLD,&rank_world);
407 /* Get the node number from the hostname to identify the nodes */
408 mynum = gmx_hostname_num();
410 /* We can't rely on MPI_IN_PLACE, so we need send and receive buffers */
411 snew(num, nrank_world);
412 snew(num_s, nrank_world);
413 snew(num_pp, nrank_world);
414 snew(num_pp_s, nrank_world);
416 num_s[rank_world] = mynum;
417 num_pp_s[rank_world] = (cr->duty & DUTY_PP) ? mynum : -1;
419 MPI_Allreduce(num_s, num, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
420 MPI_Allreduce(num_pp_s, num_pp, nrank_world, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
424 nrank_pp_intranode = 0;
425 rank_pp_intranode = 0;
426 for(i=0; i<nrank_world; i++)
436 if ((cr->duty & DUTY_PP) && num_pp[i] == mynum)
438 nrank_pp_intranode++;
450 /* Serial or thread-MPI code: we run within a single physical node */
451 nrank_intranode = cr->nnodes;
452 rank_intranode = cr->sim_nodeid;
453 nrank_pp_intranode = cr->nnodes - cr->npmenodes;
454 rank_pp_intranode = cr->nodeid;
460 if (cr->duty & DUTY_PP && cr->duty & DUTY_PME)
462 sprintf(sbuf, "PP+PME");
466 sprintf(sbuf, "%s", cr->duty & DUTY_PP ? "PP" : "PME");
468 fprintf(debug, "On %3s node %d: nrank_intranode=%d, rank_intranode=%d, "
469 "nrank_pp_intranode=%d, rank_pp_intranode=%d\n",
470 sbuf, cr->sim_nodeid,
471 nrank_intranode, rank_intranode,
472 nrank_pp_intranode, rank_pp_intranode);
475 cr->nrank_intranode = nrank_intranode;
476 cr->rank_intranode = rank_intranode;
477 cr->nrank_pp_intranode = nrank_pp_intranode;
478 cr->rank_pp_intranode = rank_pp_intranode;
482 void gmx_barrier(const t_commrec *cr)
485 gmx_call("gmx_barrier");
487 MPI_Barrier(cr->mpi_comm_mygroup);
491 void gmx_abort(int noderank,int nnodes,int errorno)
494 gmx_call("gmx_abort");
496 #ifdef GMX_THREAD_MPI
497 fprintf(stderr,"Halting program %s\n",ShortProgram());
503 fprintf(stderr,"Halting parallel program %s on CPU %d out of %d\n",
504 ShortProgram(),noderank,nnodes);
508 fprintf(stderr,"Halting program %s\n",ShortProgram());
512 MPI_Abort(MPI_COMM_WORLD,errorno);
518 void gmx_bcast(int nbytes,void *b,const t_commrec *cr)
521 gmx_call("gmx_bast");
523 MPI_Bcast(b,nbytes,MPI_BYTE,MASTERRANK(cr),cr->mpi_comm_mygroup);
527 void gmx_bcast_sim(int nbytes,void *b,const t_commrec *cr)
530 gmx_call("gmx_bast");
532 MPI_Bcast(b,nbytes,MPI_BYTE,MASTERRANK(cr),cr->mpi_comm_mysim);
536 void gmx_sumd(int nr,double r[],const t_commrec *cr)
539 gmx_call("gmx_sumd");
541 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
543 if (cr->nc.rank_intra == 0)
545 /* Use two step summing. */
546 MPI_Reduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,0,
548 /* Sum the roots of the internal (intra) buffers. */
549 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,
554 /* This is here because of the silly MPI specification
555 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
556 MPI_Reduce(r,NULL,nr,MPI_DOUBLE,MPI_SUM,0,cr->nc.comm_intra);
558 MPI_Bcast(r,nr,MPI_DOUBLE,0,cr->nc.comm_intra);
562 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,
563 cr->mpi_comm_mygroup);
568 if (nr > cr->mpb->dbuf_alloc) {
569 cr->mpb->dbuf_alloc = nr;
570 srenew(cr->mpb->dbuf,cr->mpb->dbuf_alloc);
573 /* Use two step summing */
574 MPI_Allreduce(r,cr->mpb->dbuf,nr,MPI_DOUBLE,MPI_SUM,cr->nc.comm_intra);
575 if (cr->nc.rank_intra == 0) {
576 /* Sum with the buffers reversed */
577 MPI_Allreduce(cr->mpb->dbuf,r,nr,MPI_DOUBLE,MPI_SUM,
580 MPI_Bcast(r,nr,MPI_DOUBLE,0,cr->nc.comm_intra);
582 MPI_Allreduce(r,cr->mpb->dbuf,nr,MPI_DOUBLE,MPI_SUM,
583 cr->mpi_comm_mygroup);
585 r[i] = cr->mpb->dbuf[i];
591 void gmx_sumf(int nr,float r[],const t_commrec *cr)
594 gmx_call("gmx_sumf");
596 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
598 /* Use two step summing. */
599 if (cr->nc.rank_intra == 0)
601 MPI_Reduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,0,
603 /* Sum the roots of the internal (intra) buffers */
604 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,
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,MPI_FLOAT,MPI_SUM,0,cr->nc.comm_intra);
613 MPI_Bcast(r,nr,MPI_FLOAT,0,cr->nc.comm_intra);
617 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,cr->mpi_comm_mygroup);
622 if (nr > cr->mpb->fbuf_alloc) {
623 cr->mpb->fbuf_alloc = nr;
624 srenew(cr->mpb->fbuf,cr->mpb->fbuf_alloc);
627 /* Use two step summing */
628 MPI_Allreduce(r,cr->mpb->fbuf,nr,MPI_FLOAT,MPI_SUM,cr->nc.comm_intra);
629 if (cr->nc.rank_intra == 0) {
630 /* Sum with the buffers reversed */
631 MPI_Allreduce(cr->mpb->fbuf,r,nr,MPI_FLOAT,MPI_SUM,
634 MPI_Bcast(r,nr,MPI_FLOAT,0,cr->nc.comm_intra);
636 MPI_Allreduce(r,cr->mpb->fbuf,nr,MPI_FLOAT,MPI_SUM,
637 cr->mpi_comm_mygroup);
639 r[i] = cr->mpb->fbuf[i];
645 void gmx_sumi(int nr,int r[],const t_commrec *cr)
648 gmx_call("gmx_sumi");
650 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
652 /* Use two step summing */
653 if (cr->nc.rank_intra == 0)
655 MPI_Reduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,0,cr->nc.comm_intra);
656 /* Sum with the buffers reversed */
657 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,cr->nc.comm_inter);
661 /* This is here because of the silly MPI specification
662 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
663 MPI_Reduce(r,NULL,nr,MPI_INT,MPI_SUM,0,cr->nc.comm_intra);
665 MPI_Bcast(r,nr,MPI_INT,0,cr->nc.comm_intra);
669 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,cr->mpi_comm_mygroup);
674 if (nr > cr->mpb->ibuf_alloc) {
675 cr->mpb->ibuf_alloc = nr;
676 srenew(cr->mpb->ibuf,cr->mpb->ibuf_alloc);
679 /* Use two step summing */
680 MPI_Allreduce(r,cr->mpb->ibuf,nr,MPI_INT,MPI_SUM,cr->nc.comm_intra);
681 if (cr->nc.rank_intra == 0) {
682 /* Sum with the buffers reversed */
683 MPI_Allreduce(cr->mpb->ibuf,r,nr,MPI_INT,MPI_SUM,cr->nc.comm_inter);
685 MPI_Bcast(r,nr,MPI_INT,0,cr->nc.comm_intra);
687 MPI_Allreduce(r,cr->mpb->ibuf,nr,MPI_INT,MPI_SUM,cr->mpi_comm_mygroup);
689 r[i] = cr->mpb->ibuf[i];
695 void gmx_sumli(int nr,gmx_large_int_t r[],const t_commrec *cr)
698 gmx_call("gmx_sumli");
700 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
702 /* Use two step summing */
703 if (cr->nc.rank_intra == 0)
705 MPI_Reduce(MPI_IN_PLACE,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,0,
707 /* Sum with the buffers reversed */
708 MPI_Allreduce(MPI_IN_PLACE,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,
713 /* This is here because of the silly MPI specification
714 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
715 MPI_Reduce(r,NULL,nr,GMX_MPI_LARGE_INT,MPI_SUM,0,cr->nc.comm_intra);
717 MPI_Bcast(r,nr,GMX_MPI_LARGE_INT,0,cr->nc.comm_intra);
721 MPI_Allreduce(MPI_IN_PLACE,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,cr->mpi_comm_mygroup);
726 if (nr > cr->mpb->libuf_alloc) {
727 cr->mpb->libuf_alloc = nr;
728 srenew(cr->mpb->libuf,cr->mpb->libuf_alloc);
731 /* Use two step summing */
732 MPI_Allreduce(r,cr->mpb->libuf,nr,GMX_MPI_LARGE_INT,MPI_SUM,
734 if (cr->nc.rank_intra == 0) {
735 /* Sum with the buffers reversed */
736 MPI_Allreduce(cr->mpb->libuf,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,
739 MPI_Bcast(r,nr,GMX_MPI_LARGE_INT,0,cr->nc.comm_intra);
741 MPI_Allreduce(r,cr->mpb->libuf,nr,GMX_MPI_LARGE_INT,MPI_SUM,
742 cr->mpi_comm_mygroup);
744 r[i] = cr->mpb->libuf[i];
753 void gmx_sumd_comm(int nr,double r[],MPI_Comm mpi_comm)
755 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
756 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,mpi_comm);
758 /* this function is only used in code that is not performance critical,
759 (during setup, when comm_rec is not the appropriate communication
760 structure), so this isn't as bad as it looks. */
765 MPI_Allreduce(r,buf,nr,MPI_DOUBLE,MPI_SUM,mpi_comm);
774 void gmx_sumf_comm(int nr,float r[],MPI_Comm mpi_comm)
776 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
777 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,mpi_comm);
779 /* this function is only used in code that is not performance critical,
780 (during setup, when comm_rec is not the appropriate communication
781 structure), so this isn't as bad as it looks. */
786 MPI_Allreduce(r,buf,nr,MPI_FLOAT,MPI_SUM,mpi_comm);
794 void gmx_sumd_sim(int nr,double r[],const gmx_multisim_t *ms)
797 gmx_call("gmx_sumd_sim");
799 gmx_sumd_comm(nr,r,ms->mpi_comm_masters);
803 void gmx_sumf_sim(int nr,float r[],const gmx_multisim_t *ms)
806 gmx_call("gmx_sumf_sim");
808 gmx_sumf_comm(nr,r,ms->mpi_comm_masters);
812 void gmx_sumi_sim(int nr,int r[], const gmx_multisim_t *ms)
815 gmx_call("gmx_sumi_sim");
817 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
818 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,ms->mpi_comm_masters);
820 /* this is thread-unsafe, but it will do for now: */
823 if (nr > ms->mpb->ibuf_alloc) {
824 ms->mpb->ibuf_alloc = nr;
825 srenew(ms->mpb->ibuf,ms->mpb->ibuf_alloc);
827 MPI_Allreduce(r,ms->mpb->ibuf,nr,MPI_INT,MPI_SUM,ms->mpi_comm_masters);
829 r[i] = ms->mpb->ibuf[i];
834 void gmx_sumli_sim(int nr,gmx_large_int_t r[], const gmx_multisim_t *ms)
837 gmx_call("gmx_sumli_sim");
839 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREAD_MPI)
840 MPI_Allreduce(MPI_IN_PLACE,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,
841 ms->mpi_comm_masters);
843 /* this is thread-unsafe, but it will do for now: */
846 if (nr > ms->mpb->libuf_alloc) {
847 ms->mpb->libuf_alloc = nr;
848 srenew(ms->mpb->libuf,ms->mpb->libuf_alloc);
850 MPI_Allreduce(r,ms->mpb->libuf,nr,GMX_MPI_LARGE_INT,MPI_SUM,
851 ms->mpi_comm_masters);
853 r[i] = ms->mpb->libuf[i];
859 void gmx_finalize_par(void)
862 /* Compiled without MPI, no MPI finalizing needed */
865 int initialized,finalized;
868 MPI_Initialized(&initialized);
873 /* just as a check; we don't want to finalize twice */
874 MPI_Finalized(&finalized);
880 /* We sync the processes here to try to avoid problems
881 * with buggy MPI implementations that could cause
882 * unfinished processes to terminate.
884 MPI_Barrier(MPI_COMM_WORLD);
887 if (DOMAINDECOMP(cr)) {
888 if (cr->npmenodes > 0 || cr->dd->bCartesian)
889 MPI_Comm_free(&cr->mpi_comm_mygroup);
890 if (cr->dd->bCartesian)
891 MPI_Comm_free(&cr->mpi_comm_mysim);
895 /* Apparently certain mpich implementations cause problems
896 * with MPI_Finalize. In that case comment out MPI_Finalize.
899 fprintf(debug,"Will call MPI_Finalize now\n");
901 ret = MPI_Finalize();
903 fprintf(debug,"Return code from MPI_Finalize = %d\n",ret);