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 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
43 * This module defines the interface of the actual communication routines.
48 #include "types/simple.h"
51 #include "gmx_fatal.h"
57 t_commrec *init_commrec(void);
58 /* Allocate, initialize and return the commrec. */
60 t_commrec *reinitialize_commrec_for_this_thread(const t_commrec *cro);
61 /* Initialize communication records for thread-parallel simulations.
62 Must be called on all threads before any communication takes place by
63 the individual threads. Copies the original commrec to
64 thread-local versions (a small memory leak results because we don't
65 deallocate the old shared version). */
67 void gmx_fill_commrec_from_mpi(t_commrec *cr);
68 /* Continues t_commrec construction */
70 int gmx_node_num(void);
71 /* return the number of nodes in the ring */
73 int gmx_node_rank(void);
74 /* return the rank of the node */
76 int gmx_physicalnode_id_hash(void);
77 /* Return a non-negative hash that is, hopefully, unique for each physical node.
78 * This hash is useful for determining hardware locality.
81 void gmx_setup_nodecomm(FILE *fplog, t_commrec *cr);
82 /* Sets up fast global communication for clusters with multi-core nodes */
84 void gmx_init_intranode_counters(t_commrec *cr);
85 /* Initializes intra-physical-node MPI process/thread counts and ID. */
87 gmx_bool gmx_mpi_initialized(void);
88 /* return TRUE when MPI_Init has been called.
89 * return FALSE when MPI_Init has not been called OR
90 * when GROMACS was compiled without MPI support.
93 void gmx_barrier(const t_commrec *cr);
94 /* Wait till all processes in cr->mpi_comm_mygroup have reached the barrier */
96 void gmx_bcast(int nbytes, void *b, const t_commrec *cr);
97 /* Broadcast nbytes bytes from the master to cr->mpi_comm_mygroup */
99 void gmx_bcast_sim(int nbytes, void *b, const t_commrec *cr);
100 /* Broadcast nbytes bytes from the sim master to cr->mpi_comm_mysim */
102 void gmx_sumi(int nr, int r[], const t_commrec *cr);
103 /* Calculate the global sum of an array of ints */
105 void gmx_sumli(int nr, gmx_int64_t r[], const t_commrec *cr);
106 /* Calculate the global sum of an array of large ints */
108 void gmx_sumf(int nr, float r[], const t_commrec *cr);
109 /* Calculate the global sum of an array of floats */
111 void gmx_sumd(int nr, double r[], const t_commrec *cr);
112 /* Calculate the global sum of an array of doubles */
114 void gmx_sumi_sim(int nr, int r[], const gmx_multisim_t *ms);
115 /* Calculate the sum over the simulations of an array of ints */
117 void gmx_sumli_sim(int nr, gmx_int64_t r[], const gmx_multisim_t *ms);
118 /* Calculate the sum over the simulations of an array of large ints */
120 void gmx_sumf_sim(int nr, float r[], const gmx_multisim_t *ms);
121 /* Calculate the sum over the simulations of an array of floats */
123 void gmx_sumd_sim(int nr, double r[], const gmx_multisim_t *ms);
124 /* Calculate the sum over the simulations of an array of doubles */
126 void gmx_abort(int nodeid, int nnodes, int errorno);
127 /* Abort the parallel run */
130 #define gmx_sum gmx_sumd
131 #define gmx_sum_sim gmx_sumd_sim
133 #define gmx_sum gmx_sumf
134 #define gmx_sum_sim gmx_sumf_sim
138 #define debug_gmx() do { FILE *fp = debug ? debug : stderr; \
139 if (bDebugMode()) { fprintf(fp, "NODEID=%d, %s %d\n", gmx_mpi_initialized() ? gmx_node_rank() : -1, __FILE__, __LINE__); } fflush(fp); } while (0)
149 #endif /* _network_h */