2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2008,2009,2010,2012,2014,2015,2016,2018, 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.
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.
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.
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.
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.
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.
36 /*! \libinternal \file
38 * \brief This file declares functions for (mostly) the domdec module
39 * to use MPI functionality
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.
44 * \author Berk Hess <hess@kth.se>
46 * \ingroup module_domdec
48 #ifndef GMX_DOMDEC_DOMDEC_NETWORK_H
49 #define GMX_DOMDEC_DOMDEC_NETWORK_H
51 #include "gromacs/math/vectypes.h"
52 #include "gromacs/utility/arrayref.h"
58 dddirForward
, dddirBackward
61 /*! \brief Move T values in the communication region one cell along
62 * the domain decomposition
64 * Moves in the dimension indexed by ddDimensionIndex, either forward
65 * (direction=dddirFoward) or backward (direction=dddirBackward).
67 * \todo This function template is deprecated, new calls should be
68 * made to the version taking ArrayRef parameters and this function
69 * template removed when unused.
73 ddSendrecv(const gmx_domdec_t
*dd
,
77 int numElementsToSend
,
79 int numElementsToReceive
);
81 //! Extern declaration for int specialization
84 ddSendrecv
<int>(const gmx_domdec_t
*dd
,
85 int ddDimensionIndex
, int direction
,
89 //! Extern declaration for real specialization
92 ddSendrecv
<real
>(const gmx_domdec_t
*dd
,
93 int ddDimensionIndex
, int direction
,
95 real
*buf_r
, int n_r
);
97 //! Extern declaration for rvec specialization
100 ddSendrecv
<rvec
>(const gmx_domdec_t
*dd
,
101 int ddDimensionIndex
, int direction
,
102 rvec
*buf_s
, int n_s
,
103 rvec
*buf_r
, int n_r
);
105 /*! \brief Move a view of T values in the communication region one
106 * cell along the domain decomposition
108 * Moves in the dimension indexed by ddDimensionIndex, either forward
109 * (direction=dddirFoward) or backward (direction=dddirBackward).
111 template <typename T
>
113 ddSendrecv(const gmx_domdec_t
*dd
,
114 int ddDimensionIndex
,
116 gmx::ArrayRef
<T
> sendBuffer
,
117 gmx::ArrayRef
<T
> receiveBuffer
);
119 //! Extern declaration for int specialization
122 ddSendrecv
<int>(const gmx_domdec_t
*dd
,
123 int ddDimensionIndex
,
125 gmx::ArrayRef
<int> sendBuffer
,
126 gmx::ArrayRef
<int> receiveBuffer
);
128 //! Extern declaration for real specialization
131 ddSendrecv
<real
>(const gmx_domdec_t
*dd
,
132 int ddDimensionIndex
,
134 gmx::ArrayRef
<real
> sendBuffer
,
135 gmx::ArrayRef
<real
> receiveBuffer
);
137 //! Extern declaration for gmx::RVec specialization
140 ddSendrecv
<gmx::RVec
>(const gmx_domdec_t
*dd
,
141 int ddDimensionIndex
,
143 gmx::ArrayRef
<gmx::RVec
> sendBuffer
,
144 gmx::ArrayRef
<gmx::RVec
> receiveBuffer
);
146 /*! \brief Move revc's in the comm. region one cell along the domain decomposition
148 * Moves in dimension indexed by ddimind, simultaneously in the forward
149 * and backward directions.
152 dd_sendrecv2_rvec(const struct gmx_domdec_t
*dd
,
154 rvec
*buf_s_fw
, int n_s_fw
,
155 rvec
*buf_r_fw
, int n_r_fw
,
156 rvec
*buf_s_bw
, int n_s_bw
,
157 rvec
*buf_r_bw
, int n_r_bw
);
160 /* The functions below perform the same operations as the MPI functions
161 * with the same name appendices, but over the domain decomposition
163 * The DD master node is the master for these operations.
166 /*! \brief Broadcasts \p nbytes from \p data on \p DDMASTERRANK to all PP ranks */
168 dd_bcast(const gmx_domdec_t
*dd
, int nbytes
, void *data
);
170 /*! \brief Copies \p nbytes from \p src to \p dest on \p DDMASTERRANK
171 * and then broadcasts to \p dest on all PP ranks */
173 dd_bcastc(const gmx_domdec_t
*dd
, int nbytes
, void *src
, void *dest
);
175 /*! \brief Scatters \p nbytes from \p src on \p DDMASTERRANK to all PP ranks, received in \p dest */
177 dd_scatter(const gmx_domdec_t
*dd
, int nbytes
, const void *src
, void *dest
);
179 /*! \brief Gathers \p nbytes from \p src on all PP ranks, received in \p dest on \p DDMASTERRANK */
181 dd_gather(const gmx_domdec_t
*dd
, int nbytes
, const void *src
, void *dest
);
183 /*! \brief Scatters \p scounts bytes from \p src on \p DDMASTERRANK to all PP ranks, receiving \p rcount bytes in \p dest.
185 * See man MPI_Scatterv for details of how to construct scounts and disps.
186 * If rcount==0, rbuf is allowed to be NULL */
188 dd_scatterv(const gmx_domdec_t
*dd
,
189 int *scounts
, int *disps
, const void *sbuf
,
190 int rcount
, void *rbuf
);
192 /*! \brief Gathers \p rcount bytes from \p src on all PP ranks, received in \p scounts bytes in \p dest on \p DDMASTERRANK.
194 * See man MPI_Gatherv for details of how to construct scounts and disps.
196 * If scount==0, sbuf is allowed to be NULL */
198 dd_gatherv(const gmx_domdec_t
*dd
,
199 int scount
, const void *sbuf
,
200 int *rcounts
, int *disps
, void *rbuf
);