Make domdec a proper module
[alexxy/gromacs.git] / src / gromacs / domdec / domdec_network.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2008,2009,2010,2012,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35
36 /*! \libinternal \file
37  *
38  * \brief This file declares functions for (mostly) the domdec module
39  * to use MPI functionality
40  *
41  * \todo Wrap the raw dd_bcast in md.cpp into a higher-level function
42  * in the domdec module, then this file can be module-internal.
43  *
44  * \author Berk Hess <hess@kth.se>
45  * \inlibraryapi
46  * \ingroup module_domdec
47  */
48 #ifndef GMX_DOMDEC_DOMDEC_NETWORK_H
49 #define GMX_DOMDEC_DOMDEC_NETWORK_H
50
51 #include "gromacs/legacyheaders/typedefs.h"
52
53 /* \brief */
54 enum {
55     dddirForward, dddirBackward
56 };
57
58 /*! \brief Move integers in the communication region one cell along
59  * the domain decomposition
60  *
61  * Moves in the dimension indexed by ddimind, either forward
62  * (direction=dddirFoward) or backward (direction=dddirBackward).
63  */
64 void
65 dd_sendrecv_int(const gmx_domdec_t *dd,
66                 int ddimind, int direction,
67                 int *buf_s, int n_s,
68                 int *buf_r, int n_r);
69
70 /*! \brief Move reals in the comm. region one cell along the domain decomposition
71  *
72  * Moves in the dimension indexed by ddimind, either forward
73  * (direction=dddirFoward) or backward (direction=dddirBackward).
74  */
75 void
76 dd_sendrecv_real(const gmx_domdec_t *dd,
77                  int ddimind, int direction,
78                  real *buf_s, int n_s,
79                  real *buf_r, int n_r);
80
81 /*! \brief Move revc's in the comm. region one cell along the domain decomposition
82  *
83  * Moves in dimension indexed by ddimind, either forward
84  * (direction=dddirFoward) or backward (direction=dddirBackward).
85  */
86 void
87 dd_sendrecv_rvec(const gmx_domdec_t *dd,
88                  int ddimind, int direction,
89                  rvec *buf_s, int n_s,
90                  rvec *buf_r, int n_r);
91
92
93 /*! \brief Move revc's in the comm. region one cell along the domain decomposition
94  *
95  * Moves in dimension indexed by ddimind, simultaneously in the forward
96  * and backward directions.
97  */
98 void
99 dd_sendrecv2_rvec(const gmx_domdec_t *dd,
100                   int ddimind,
101                   rvec *buf_s_fw, int n_s_fw,
102                   rvec *buf_r_fw, int n_r_fw,
103                   rvec *buf_s_bw, int n_s_bw,
104                   rvec *buf_r_bw, int n_r_bw);
105
106
107 /* The functions below perform the same operations as the MPI functions
108  * with the same name appendices, but over the domain decomposition
109  * nodes only.
110  * The DD master node is the master for these operations.
111  */
112
113 /*! \brief Broadcasts \p nbytes from \p data on \p DDMASTERRANK to all PP ranks */
114 void
115 dd_bcast(gmx_domdec_t *dd, int nbytes, void *data);
116
117 /*! \brief Copies \p nbytes from \p src to \p dest on \p DDMASTERRANK
118  * and then broadcasts to \p dest on all PP ranks */
119 void
120 dd_bcastc(gmx_domdec_t *dd, int nbytes, void *src, void *dest);
121
122 /*! \brief Scatters \p nbytes from \p src on \p DDMASTERRANK to all PP ranks, received in \p dest */
123 void
124 dd_scatter(gmx_domdec_t *dd, int nbytes, void *src, void *dest);
125
126 /*! \brief Gathers \p nbytes from \p src on all PP ranks, received in \p dest on \p DDMASTERRANK */
127 void
128 dd_gather(gmx_domdec_t *dd, int nbytes, void *src, void *dest);
129
130 /*! \brief Scatters \p scounts bytes from \p src on \p DDMASTERRANK to all PP ranks, receiving \p rcount bytes in \p dest.
131  *
132  * See man MPI_Scatterv for details of how to construct scounts and disps.
133  * If rcount==0, rbuf is allowed to be NULL */
134 void
135 dd_scatterv(gmx_domdec_t *dd,
136             int *scounts, int *disps, void *sbuf,
137             int rcount, void *rbuf);
138
139 /*! \brief Gathers \p rcount bytes from \p src on all PP ranks, received in \p scounts bytes in \p dest on \p DDMASTERRANK.
140  *
141  * See man MPI_Gatherv for details of how to construct scounts and disps.
142  *
143  * If scount==0, sbuf is allowed to be NULL */
144 void
145 dd_gatherv(gmx_domdec_t *dd,
146            int scount, void *sbuf,
147            int *rcounts, int *disps, void *rbuf);
148
149 #endif