Merge branch 'master' into rotation
[gromacs.git] / src / gmxlib / network.c
blob175f75bb666f5c82ebf103b88c17757f881b904e
1 /*
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.
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
32 * And Hey:
33 * GROningen Mixture of Alchemy and Childrens' Stories
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
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"
48 #ifdef GMX_LIB_MPI
49 #include <mpi.h>
50 #endif
52 #ifdef GMX_THREADS
53 #include "tmpi.h"
54 #endif
56 #include "mpelogging.h"
58 /* The source code in this file should be thread-safe.
59 Please keep it that way. */
61 bool gmx_mpi_initialized(void)
63 int n;
64 #ifndef GMX_MPI
65 return 0;
66 #else
67 MPI_Initialized(&n);
69 return n;
70 #endif
73 int gmx_setup(int *argc,char **argv,int *nnodes)
75 #ifndef GMX_MPI
76 gmx_call("gmx_setup");
77 return 0;
78 #else
79 char buf[256];
80 int resultlen; /* actual length of node name */
81 int i,flag;
82 int mpi_num_nodes;
83 int mpi_my_rank;
84 char mpi_hostname[MPI_MAX_PROCESSOR_NAME];
86 /* Call the MPI routines */
87 #ifdef GMX_LIB_MPI
88 #ifdef GMX_FAHCORE
89 (void) fah_MPI_Init(argc,&argv);
90 #else
91 (void) MPI_Init(argc,&argv);
92 #endif
93 #endif
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 #ifdef USE_MPE
100 /* MPE logging routines. Get event IDs from MPE: */
101 /* General events */
102 ev_timestep1 = MPE_Log_get_event_number( );
103 ev_timestep2 = MPE_Log_get_event_number( );
104 ev_force_start = MPE_Log_get_event_number( );
105 ev_force_finish = MPE_Log_get_event_number( );
106 ev_do_fnbf_start = MPE_Log_get_event_number( );
107 ev_do_fnbf_finish = MPE_Log_get_event_number( );
108 ev_ns_start = MPE_Log_get_event_number( );
109 ev_ns_finish = MPE_Log_get_event_number( );
110 ev_calc_bonds_start = MPE_Log_get_event_number( );
111 ev_calc_bonds_finish = MPE_Log_get_event_number( );
112 ev_global_stat_start = MPE_Log_get_event_number( );
113 ev_global_stat_finish = MPE_Log_get_event_number( );
114 ev_virial_start = MPE_Log_get_event_number( );
115 ev_virial_finish = MPE_Log_get_event_number( );
117 /* Enforced rotation */
118 ev_flexll_start = MPE_Log_get_event_number( );
119 ev_flexll_finish = MPE_Log_get_event_number( );
120 ev_add_rot_forces_start = MPE_Log_get_event_number( );
121 ev_add_rot_forces_finish = MPE_Log_get_event_number( );
122 ev_rotcycles_start = MPE_Log_get_event_number( );
123 ev_rotcycles_finish = MPE_Log_get_event_number( );
124 ev_forcecycles_start = MPE_Log_get_event_number( );
125 ev_forcecycles_finish = MPE_Log_get_event_number( );
127 /* Shift related events */
128 ev_shift_start = MPE_Log_get_event_number( );
129 ev_shift_finish = MPE_Log_get_event_number( );
130 ev_unshift_start = MPE_Log_get_event_number( );
131 ev_unshift_finish = MPE_Log_get_event_number( );
132 ev_mk_mshift_start = MPE_Log_get_event_number( );
133 ev_mk_mshift_finish = MPE_Log_get_event_number( );
135 /* PME related events */
136 ev_pme_start = MPE_Log_get_event_number( );
137 ev_pme_finish = MPE_Log_get_event_number( );
138 ev_spread_on_grid_start = MPE_Log_get_event_number( );
139 ev_spread_on_grid_finish = MPE_Log_get_event_number( );
140 ev_sum_qgrid_start = MPE_Log_get_event_number( );
141 ev_sum_qgrid_finish = MPE_Log_get_event_number( );
142 ev_gmxfft3d_start = MPE_Log_get_event_number( );
143 ev_gmxfft3d_finish = MPE_Log_get_event_number( );
144 ev_solve_pme_start = MPE_Log_get_event_number( );
145 ev_solve_pme_finish = MPE_Log_get_event_number( );
146 ev_gather_f_bsplines_start = MPE_Log_get_event_number( );
147 ev_gather_f_bsplines_finish= MPE_Log_get_event_number( );
148 ev_reduce_start = MPE_Log_get_event_number( );
149 ev_reduce_finish = MPE_Log_get_event_number( );
150 ev_rscatter_start = MPE_Log_get_event_number( );
151 ev_rscatter_finish = MPE_Log_get_event_number( );
152 ev_alltoall_start = MPE_Log_get_event_number( );
153 ev_alltoall_finish = MPE_Log_get_event_number( );
154 ev_pmeredist_start = MPE_Log_get_event_number( );
155 ev_pmeredist_finish = MPE_Log_get_event_number( );
156 ev_init_pme_start = MPE_Log_get_event_number( );
157 ev_init_pme_finish = MPE_Log_get_event_number( );
158 ev_send_coordinates_start = MPE_Log_get_event_number( );
159 ev_send_coordinates_finish = MPE_Log_get_event_number( );
160 ev_update_fr_start = MPE_Log_get_event_number( );
161 ev_update_fr_finish = MPE_Log_get_event_number( );
162 ev_clear_rvecs_start = MPE_Log_get_event_number( );
163 ev_clear_rvecs_finish = MPE_Log_get_event_number( );
164 ev_update_start = MPE_Log_get_event_number( );
165 ev_update_finish = MPE_Log_get_event_number( );
166 ev_output_start = MPE_Log_get_event_number( );
167 ev_output_finish = MPE_Log_get_event_number( );
168 ev_sum_lrforces_start = MPE_Log_get_event_number( );
169 ev_sum_lrforces_finish = MPE_Log_get_event_number( );
170 ev_sort_start = MPE_Log_get_event_number( );
171 ev_sort_finish = MPE_Log_get_event_number( );
172 ev_sum_qgrid_start = MPE_Log_get_event_number( );
173 ev_sum_qgrid_finish = MPE_Log_get_event_number( );
175 /* Essential dynamics related events */
176 ev_edsam_start = MPE_Log_get_event_number( );
177 ev_edsam_finish = MPE_Log_get_event_number( );
178 ev_get_coords_start = MPE_Log_get_event_number( );
179 ev_get_coords_finish = MPE_Log_get_event_number( );
180 ev_ed_apply_cons_start = MPE_Log_get_event_number( );
181 ev_ed_apply_cons_finish = MPE_Log_get_event_number( );
182 ev_fit_to_reference_start = MPE_Log_get_event_number( );
183 ev_fit_to_reference_finish = MPE_Log_get_event_number( );
185 /* describe events: */
186 if ( mpi_my_rank == 0 )
188 /* General events */
189 MPE_Describe_state(ev_timestep1, ev_timestep2, "timestep START", "magenta" );
190 MPE_Describe_state(ev_force_start, ev_force_finish, "force", "cornflower blue" );
191 MPE_Describe_state(ev_do_fnbf_start, ev_do_fnbf_finish, "do_fnbf", "navy" );
192 MPE_Describe_state(ev_ns_start, ev_ns_finish, "neighbor search", "tomato" );
193 MPE_Describe_state(ev_calc_bonds_start, ev_calc_bonds_finish, "bonded forces", "slate blue" );
194 MPE_Describe_state(ev_global_stat_start, ev_global_stat_finish, "global stat", "firebrick3");
195 MPE_Describe_state(ev_update_fr_start, ev_update_fr_finish, "update forcerec", "goldenrod");
196 MPE_Describe_state(ev_clear_rvecs_start, ev_clear_rvecs_finish, "clear rvecs", "bisque");
197 MPE_Describe_state(ev_update_start, ev_update_finish, "update", "cornsilk");
198 MPE_Describe_state(ev_output_start, ev_output_finish, "output", "black");
199 MPE_Describe_state(ev_virial_start, ev_virial_finish, "calc_virial", "thistle4");
201 /* Enforced rotation */
202 MPE_Describe_state(ev_flexll_start, ev_flexll_finish, "flex lowlevel", "navajo white");
203 MPE_Describe_state(ev_add_rot_forces_start, ev_add_rot_forces_finish, "add rot forces", "green");
204 MPE_Describe_state(ev_rotcycles_start, ev_rotcycles_finish, "count rot cyc", "moccasin");
205 MPE_Describe_state(ev_forcecycles_start, ev_forcecycles_finish, "count force cyc", "powder blue");
207 /* PME related events */
208 MPE_Describe_state(ev_pme_start, ev_pme_finish, "doing PME", "grey" );
209 MPE_Describe_state(ev_spread_on_grid_start, ev_spread_on_grid_finish, "spread", "dark orange" );
210 MPE_Describe_state(ev_sum_qgrid_start, ev_sum_qgrid_finish, "sum qgrid", "slate blue");
211 MPE_Describe_state(ev_gmxfft3d_start, ev_gmxfft3d_finish, "fft3d", "snow2" );
212 MPE_Describe_state(ev_solve_pme_start, ev_solve_pme_finish, "solve PME", "indian red" );
213 MPE_Describe_state(ev_gather_f_bsplines_start, ev_gather_f_bsplines_finish, "bsplines", "light sea green" );
214 MPE_Describe_state(ev_reduce_start, ev_reduce_finish, "reduce", "cyan1" );
215 MPE_Describe_state(ev_rscatter_start, ev_rscatter_finish, "rscatter", "cyan3" );
216 MPE_Describe_state(ev_alltoall_start, ev_alltoall_finish, "alltoall", "LightCyan4" );
217 MPE_Describe_state(ev_pmeredist_start, ev_pmeredist_finish, "pmeredist", "thistle" );
218 MPE_Describe_state(ev_init_pme_start, ev_init_pme_finish, "init PME", "snow4");
219 MPE_Describe_state(ev_send_coordinates_start, ev_send_coordinates_finish, "send_coordinates","blue");
220 MPE_Describe_state(ev_sum_lrforces_start, ev_sum_lrforces_finish, "sum_LRforces", "lime green");
221 MPE_Describe_state(ev_sort_start, ev_sort_finish, "sort pme atoms", "brown");
222 MPE_Describe_state(ev_sum_qgrid_start, ev_sum_qgrid_finish, "sum charge grid", "medium orchid");
224 /* Shift related events */
225 MPE_Describe_state(ev_shift_start, ev_shift_finish, "shift", "orange");
226 MPE_Describe_state(ev_unshift_start, ev_unshift_finish, "unshift", "dark orange");
227 MPE_Describe_state(ev_mk_mshift_start, ev_mk_mshift_finish, "mk_mshift", "maroon");
229 /* Essential dynamics related events */
230 MPE_Describe_state(ev_edsam_start, ev_edsam_finish, "EDSAM", "deep sky blue");
231 MPE_Describe_state(ev_get_coords_start, ev_get_coords_finish, "ED get coords", "steel blue");
232 MPE_Describe_state(ev_ed_apply_cons_start, ev_ed_apply_cons_finish, "ED apply constr", "forest green");
233 MPE_Describe_state(ev_fit_to_reference_start, ev_fit_to_reference_finish, "ED fit to ref", "lavender");
236 MPE_Init_log();
237 #endif
239 fprintf(stderr,"NNODES=%d, MYRANK=%d, HOSTNAME=%s\n",
240 mpi_num_nodes,mpi_my_rank,mpi_hostname);
242 *nnodes=mpi_num_nodes;
244 return mpi_my_rank;
245 #endif
248 int gmx_node_num(void)
250 #ifndef GMX_MPI
251 return 1;
252 #else
253 int i;
254 (void) MPI_Comm_size(MPI_COMM_WORLD, &i);
255 return i;
256 #endif
259 int gmx_node_rank(void)
261 #ifndef GMX_MPI
262 return 0;
263 #else
264 int i;
265 (void) MPI_Comm_rank(MPI_COMM_WORLD, &i);
266 return i;
267 #endif
270 void gmx_setup_nodecomm(FILE *fplog,t_commrec *cr)
272 gmx_nodecomm_t *nc;
273 int n,rank,resultlen,hostnum,i,ng,ni;
274 #ifdef GMX_MPI
275 char mpi_hostname[MPI_MAX_PROCESSOR_NAME];
276 #endif
278 /* Many MPI implementations do not optimize MPI_Allreduce
279 * (and probably also other global communication calls)
280 * for multi-core nodes connected by a network.
281 * We can optimize such communication by using one MPI call
282 * within each node and one between the nodes.
283 * For MVAPICH2 and Intel MPI this reduces the time for
284 * the global_stat communication by 25%
285 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
286 * B. Hess, November 2007
289 nc = &cr->nc;
291 nc->bUse = FALSE;
292 #ifndef GMX_THREADS
293 if (getenv("GMX_NO_NODECOMM") == NULL) {
294 #ifdef GMX_MPI
295 MPI_Comm_size(cr->mpi_comm_mygroup,&n);
296 MPI_Comm_rank(cr->mpi_comm_mygroup,&rank);
297 MPI_Get_processor_name(mpi_hostname,&resultlen);
298 /* This procedure can only differentiate nodes with host names
299 * that end on unique numbers.
301 i = resultlen - 1;
302 while(i > 0 && isdigit(mpi_hostname[i-1]))
303 i--;
304 if (isdigit(mpi_hostname[i])) {
305 hostnum = strtol(mpi_hostname+i, NULL, 0);
306 } else {
307 hostnum = 0;
310 if (debug) {
311 fprintf(debug,
312 "In gmx_setup_nodecomm: splitting communicator of size %d\n",
314 fprintf(debug,"In gmx_setup_nodecomm: hostname '%s', hostnum %d\n",
315 mpi_hostname,hostnum);
318 /* The intra-node communicator, split on node number */
319 MPI_Comm_split(cr->mpi_comm_mygroup,hostnum,rank,&nc->comm_intra);
320 MPI_Comm_rank(nc->comm_intra,&nc->rank_intra);
321 if (debug) {
322 fprintf(debug,"In gmx_setup_nodecomm: node rank %d rank_intra %d\n",
323 rank,nc->rank_intra);
325 /* The inter-node communicator, split on rank_intra.
326 * We actually only need the one for rank=0,
327 * but it is easier to create them all.
329 MPI_Comm_split(cr->mpi_comm_mygroup,nc->rank_intra,rank,&nc->comm_inter);
330 /* Check if this really created two step communication */
331 MPI_Comm_size(nc->comm_inter,&ng);
332 MPI_Comm_size(nc->comm_intra,&ni);
333 if ((ng > 1 && ng < n) || (ni > 1 && ni < n) ) {
334 nc->bUse = TRUE;
335 if (fplog)
336 fprintf(fplog,"Using two step summing over %d groups of on average %.1f processes\n\n",ng,(real)n/(real)ng);
337 if (nc->rank_intra > 0)
338 MPI_Comm_free(&nc->comm_inter);
339 } else {
340 /* One group or all processes in a separate group, use normal summing */
341 MPI_Comm_free(&nc->comm_inter);
342 MPI_Comm_free(&nc->comm_intra);
344 #endif
346 #endif
349 void gmx_barrier(const t_commrec *cr)
351 #ifndef GMX_MPI
352 gmx_call("gmx_barrier");
353 #else
354 MPI_Barrier(cr->mpi_comm_mygroup);
355 #endif
358 void gmx_abort(int noderank,int nnodes,int errorno)
360 #ifndef GMX_MPI
361 gmx_call("gmx_abort");
362 #else
363 if (nnodes > 1)
364 fprintf(stderr,"Halting parallel program %s on CPU %d out of %d\n",
365 ShortProgram(),noderank,nnodes);
366 else
367 fprintf(stderr,"Halting program %s\n",ShortProgram());
369 thanx(stderr);
370 MPI_Abort(MPI_COMM_WORLD,errorno);
371 exit(1);
372 #endif
375 void gmx_bcast(int nbytes,void *b,const t_commrec *cr)
377 #ifndef GMX_MPI
378 gmx_call("gmx_bast");
379 #else
380 MPI_Bcast(b,nbytes,MPI_BYTE,MASTERRANK(cr),cr->mpi_comm_mygroup);
381 #endif
384 void gmx_bcast_sim(int nbytes,void *b,const t_commrec *cr)
386 #ifndef GMX_MPI
387 gmx_call("gmx_bast");
388 #else
389 MPI_Bcast(b,nbytes,MPI_BYTE,MASTERRANK(cr),cr->mpi_comm_mysim);
390 #endif
393 void gmx_sumd(int nr,double r[],const t_commrec *cr)
395 #ifndef GMX_MPI
396 gmx_call("gmx_sumd");
397 #else
398 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
399 if (cr->nc.bUse) {
400 if (cr->nc.rank_intra == 0)
402 /* Use two step summing. */
403 MPI_Reduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,0,
404 cr->nc.comm_intra);
405 /* Sum the roots of the internal (intra) buffers. */
406 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,
407 cr->nc.comm_inter);
409 else
411 /* This is here because of the silly MPI specification
412 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
413 MPI_Reduce(r,NULL,nr,MPI_DOUBLE,MPI_SUM,0,cr->nc.comm_intra);
415 MPI_Bcast(r,nr,MPI_DOUBLE,0,cr->nc.comm_intra);
417 else
419 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,
420 cr->mpi_comm_mygroup);
422 #else
423 /* this is thread-unsafe, but it will do for now: */
424 static double *buf=NULL;
425 static int nalloc=0;
426 int i;
428 if (nr > nalloc) {
429 nalloc = nr;
430 srenew(buf,nalloc);
432 if (cr->nc.bUse) {
433 /* Use two step summing */
434 MPI_Allreduce(r,buf,nr,MPI_DOUBLE,MPI_SUM,cr->nc.comm_intra);
435 if (cr->nc.rank_intra == 0) {
436 /* Sum with the buffers reversed */
437 MPI_Allreduce(buf,r,nr,MPI_DOUBLE,MPI_SUM,cr->nc.comm_inter);
439 MPI_Bcast(r,nr,MPI_DOUBLE,0,cr->nc.comm_intra);
440 } else {
441 MPI_Allreduce(r,buf,nr,MPI_DOUBLE,MPI_SUM,cr->mpi_comm_mygroup);
442 for(i=0; i<nr; i++)
443 r[i] = buf[i];
445 #endif
446 #endif
449 void gmx_sumf(int nr,float r[],const t_commrec *cr)
451 #ifndef GMX_MPI
452 gmx_call("gmx_sumf");
453 #else
454 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
455 if (cr->nc.bUse) {
456 /* Use two step summing. */
457 if (cr->nc.rank_intra == 0)
459 MPI_Reduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,0,
460 cr->nc.comm_intra);
461 /* Sum the roots of the internal (intra) buffers */
462 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,
463 cr->nc.comm_inter);
465 else
467 /* This is here because of the silly MPI specification
468 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
469 MPI_Reduce(r,NULL,nr,MPI_FLOAT,MPI_SUM,0,cr->nc.comm_intra);
471 MPI_Bcast(r,nr,MPI_FLOAT,0,cr->nc.comm_intra);
473 else
475 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,cr->mpi_comm_mygroup);
477 #else
478 /* this is thread-unsafe, but it will do for now: */
479 static float *buf=NULL;
480 static int nalloc=0;
481 int i;
483 if (nr > nalloc) {
484 nalloc = nr;
485 srenew(buf,nalloc);
487 if (cr->nc.bUse) {
488 /* Use two step summing */
489 MPI_Allreduce(r,buf,nr,MPI_FLOAT,MPI_SUM,cr->nc.comm_intra);
490 if (cr->nc.rank_intra == 0) {
491 /* Sum with the buffers reversed */
492 MPI_Allreduce(buf,r,nr,MPI_FLOAT,MPI_SUM,cr->nc.comm_inter);
494 MPI_Bcast(r,nr,MPI_FLOAT,0,cr->nc.comm_intra);
495 } else {
496 MPI_Allreduce(r,buf,nr,MPI_FLOAT,MPI_SUM,cr->mpi_comm_mygroup);
497 for(i=0; i<nr; i++)
498 r[i] = buf[i];
500 #endif
501 #endif
504 void gmx_sumi(int nr,int r[],const t_commrec *cr)
506 #ifndef GMX_MPI
507 gmx_call("gmx_sumi");
508 #else
509 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
510 if (cr->nc.bUse) {
511 /* Use two step summing */
512 if (cr->nc.rank_intra == 0)
514 MPI_Reduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,0,cr->nc.comm_intra);
515 /* Sum with the buffers reversed */
516 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,cr->nc.comm_inter);
518 else
520 /* This is here because of the silly MPI specification
521 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
522 MPI_Reduce(r,NULL,nr,MPI_INT,MPI_SUM,0,cr->nc.comm_intra);
524 MPI_Bcast(r,nr,MPI_INT,0,cr->nc.comm_intra);
526 else
528 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,cr->mpi_comm_mygroup);
530 #else
531 /* this is thread-unsafe, but it will do for now: */
532 static int *buf=NULL;
533 static int nalloc=0;
534 int i;
536 if (nr > nalloc) {
537 nalloc = nr;
538 srenew(buf,nalloc);
540 if (cr->nc.bUse) {
541 /* Use two step summing */
542 MPI_Allreduce(r,buf,nr,MPI_INT,MPI_SUM,cr->nc.comm_intra);
543 if (cr->nc.rank_intra == 0) {
544 /* Sum with the buffers reversed */
545 MPI_Allreduce(buf,r,nr,MPI_INT,MPI_SUM,cr->nc.comm_inter);
547 MPI_Bcast(r,nr,MPI_INT,0,cr->nc.comm_intra);
548 } else {
549 MPI_Allreduce(r,buf,nr,MPI_INT,MPI_SUM,cr->mpi_comm_mygroup);
550 for(i=0; i<nr; i++)
551 r[i] = buf[i];
553 #endif
554 #endif
557 #ifdef GMX_MPI
558 void gmx_sumd_comm(int nr,double r[],MPI_Comm mpi_comm)
560 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
561 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,mpi_comm);
562 #else
563 /* this is thread-unsafe, but it will do for now: */
564 static double *buf=NULL;
565 static int nalloc=0;
566 int i;
568 if (nr > nalloc) {
569 nalloc = nr;
570 srenew(buf,nalloc);
572 MPI_Allreduce(r,buf,nr,MPI_DOUBLE,MPI_SUM,mpi_comm);
573 for(i=0; i<nr; i++)
574 r[i] = buf[i];
575 #endif
577 #endif
579 #ifdef GMX_MPI
580 void gmx_sumf_comm(int nr,float r[],MPI_Comm mpi_comm)
582 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
583 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,mpi_comm);
584 #else
585 /* this is thread-unsafe, but it will do for now: */
586 static float *buf=NULL;
587 static int nalloc=0;
588 int i;
590 if (nr > nalloc) {
591 nalloc = nr;
592 srenew(buf,nalloc);
594 MPI_Allreduce(r,buf,nr,MPI_FLOAT,MPI_SUM,mpi_comm);
595 for(i=0; i<nr; i++)
596 r[i] = buf[i];
597 #endif
599 #endif
601 void gmx_sumd_sim(int nr,double r[],const gmx_multisim_t *ms)
603 #ifndef GMX_MPI
604 gmx_call("gmx_sumd");
605 #else
606 gmx_sumd_comm(nr,r,ms->mpi_comm_masters);
607 #endif
610 void gmx_sumf_sim(int nr,float r[],const gmx_multisim_t *ms)
612 #ifndef GMX_MPI
613 gmx_call("gmx_sumf");
614 #else
615 gmx_sumf_comm(nr,r,ms->mpi_comm_masters);
616 #endif
619 void gmx_sumi_sim(int nr,int r[],const gmx_multisim_t *ms)
621 #ifndef GMX_MPI
622 gmx_call("gmx_sumd");
623 #else
624 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
625 MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,ms->mpi_comm_masters);
626 #else
627 /* this is thread-unsafe, but it will do for now: */
628 static int *buf=NULL;
629 static int nalloc=0;
630 int i;
632 if (nr > nalloc) {
633 nalloc = nr;
634 srenew(buf,nalloc);
636 MPI_Allreduce(r,buf,nr,MPI_INT,MPI_SUM,ms->mpi_comm_masters);
637 for(i=0; i<nr; i++)
638 r[i] = buf[i];
639 #endif
640 #endif
643 void gmx_finalize(void)
645 #ifndef GMX_MPI
646 gmx_call("gmx_finalize");
647 #else
648 int ret;
650 /* just as a check; we don't want to finalize twice */
651 int finalized;
652 MPI_Finalized(&finalized);
653 if (finalized)
654 return;
656 /* We sync the processes here to try to avoid problems
657 * with buggy MPI implementations that could cause
658 * unfinished processes to terminate.
660 MPI_Barrier(MPI_COMM_WORLD);
663 if (DOMAINDECOMP(cr)) {
664 if (cr->npmenodes > 0 || cr->dd->bCartesian)
665 MPI_Comm_free(&cr->mpi_comm_mygroup);
666 if (cr->dd->bCartesian)
667 MPI_Comm_free(&cr->mpi_comm_mysim);
671 /* Apparently certain mpich implementations cause problems
672 * with MPI_Finalize. In that case comment out MPI_Finalize.
674 if (debug)
675 fprintf(debug,"Will call MPI_Finalize now\n");
677 ret = MPI_Finalize();
678 if (debug)
679 fprintf(debug,"Return code from MPI_Finalize = %d\n",ret);
680 #endif