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 * Gromacs Runs On Most of All Computer Systems
41 * This module defines the interface of the actual communication routines.
45 #include "visibility.h"
46 #include "types/simple.h"
47 #include "types/commrec.h"
50 #include "gmx_fatal.h"
56 int gmx_setup(int *argc
,char **argv
,int *nnodes
);
57 /* Initializes the parallel communication, return the ID of the node */
59 int gmx_node_num(void);
60 /* return the number of nodes in the ring */
62 int gmx_node_rank(void);
63 /* return the rank of the node */
66 int gmx_hostname_num(void);
67 /* If the first part of the hostname (up to the first dot) ends with a number, returns this number.
68 If the first part of the hostname does not ends in a number (0-9 characters), returns 0.
72 void gmx_setup_nodecomm(FILE *fplog
,t_commrec
*cr
);
73 /* Sets up fast global communication for clusters with multi-core nodes */
76 void gmx_init_intra_counters(t_commrec
*cr
);
77 /* Initializes intra-node process counts and ID. */
79 gmx_bool
gmx_mpi_initialized(void);
80 /* return TRUE when MPI_Init has been called.
81 * return FALSE when MPI_Init has not been called OR
82 * when GROMACS was compiled without MPI support.
85 void gmx_barrier(const t_commrec
*cr
);
86 /* Wait till all processes in cr->mpi_comm_mygroup have reached the barrier */
89 void gmx_bcast(int nbytes
,void *b
,const t_commrec
*cr
);
90 /* Broadcast nbytes bytes from the master to cr->mpi_comm_mygroup */
93 void gmx_bcast_sim(int nbytes
,void *b
,const t_commrec
*cr
);
94 /* Broadcast nbytes bytes from the sim master to cr->mpi_comm_mysim */
97 void gmx_sumi(int nr
,int r
[],const t_commrec
*cr
);
98 /* Calculate the global sum of an array of ints */
100 void gmx_sumli(int nr
,gmx_large_int_t r
[],const t_commrec
*cr
);
101 /* Calculate the global sum of an array of large ints */
104 void gmx_sumf(int nr
,float r
[],const t_commrec
*cr
);
105 /* Calculate the global sum of an array of floats */
108 void gmx_sumd(int nr
,double r
[],const t_commrec
*cr
);
109 /* Calculate the global sum of an array of doubles */
111 void gmx_sumf_comm(int nr
,float r
[],MPI_Comm mpi_comm
);
112 /* Calculate the global sum of an array of floats */
114 void gmx_sumd_comm(int nr
,double r
[],MPI_Comm mpi_comm
);
115 /* Calculate the global sum of an array of doubles */
118 void gmx_sumi_sim(int nr
,int r
[],const gmx_multisim_t
*ms
);
119 /* Calculate the sum over the simulations of an array of ints */
122 void gmx_sumli_sim(int nr
,gmx_large_int_t r
[],const gmx_multisim_t
*ms
);
123 /* Calculate the sum over the simulations of an array of large ints */
126 void gmx_sumf_sim(int nr
,float r
[],const gmx_multisim_t
*ms
);
127 /* Calculate the sum over the simulations of an array of floats */
129 void gmx_sumd_sim(int nr
,double r
[],const gmx_multisim_t
*ms
);
130 /* Calculate the sum over the simulations of an array of doubles */
132 void gmx_abort(int nodeid
,int nnodes
,int errorno
);
133 /* Abort the parallel run */
136 void gmx_finalize_par(void);
137 /* Finish the parallel run in an ordered manner */
140 #define gmx_sum_comm gmx_sumd_comm
141 #define gmx_sum gmx_sumd
142 #define gmx_sum_sim gmx_sumd_sim
144 #define gmx_sum_comm gmx_sumf_comm
145 #define gmx_sum gmx_sumf
146 #define gmx_sum_sim gmx_sumf_sim
150 #define debug_gmx() do { FILE *fp=debug ? debug : stderr;\
151 if (bDebugMode()) fprintf(fp,"NODEID=%d, %s %d\n",gmx_mpi_initialized() ? gmx_node_rank() : -1,__FILE__,__LINE__); fflush(fp); } while (0)
161 #endif /* _network_h */