2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2008-2018, The GROMACS development team.
5 * Copyright (c) 2019,2020, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
37 /*! \libinternal \file
39 * \brief This file declares functions for (mostly) the domdec module
40 * to use MPI functionality
42 * \todo Wrap the raw dd_bcast in md.cpp into a higher-level function
43 * in the domdec module, then this file can be module-internal.
45 * \author Berk Hess <hess@kth.se>
47 * \ingroup module_domdec
49 #ifndef GMX_DOMDEC_DOMDEC_NETWORK_H
50 #define GMX_DOMDEC_DOMDEC_NETWORK_H
52 #include "gromacs/math/vectypes.h"
68 /*! \brief Move T values in the communication region one cell along
69 * the domain decomposition
71 * Moves in the dimension indexed by ddDimensionIndex, either forward
72 * (direction=dddirFoward) or backward (direction=dddirBackward).
74 * \todo This function template is deprecated, new calls should be
75 * made to the version taking ArrayRef parameters and this function
76 * template removed when unused.
79 void ddSendrecv(const gmx_domdec_t* dd,
83 int numElementsToSend,
85 int numElementsToReceive);
87 //! Extern declaration for int specialization
88 extern template void ddSendrecv<int>(const gmx_domdec_t* dd,
96 //! Extern declaration for real specialization
97 extern template void ddSendrecv<real>(const gmx_domdec_t* dd,
105 //! Extern declaration for rvec specialization
106 extern template void ddSendrecv<rvec>(const gmx_domdec_t* dd,
107 int ddDimensionIndex,
114 /*! \brief Move a view of T values in the communication region one
115 * cell along the domain decomposition
117 * Moves in the dimension indexed by ddDimensionIndex, either forward
118 * (direction=dddirFoward) or backward (direction=dddirBackward).
121 void ddSendrecv(const gmx_domdec_t* dd,
122 int ddDimensionIndex,
124 gmx::ArrayRef<T> sendBuffer,
125 gmx::ArrayRef<T> receiveBuffer);
127 //! Extern declaration for int specialization
128 extern template void ddSendrecv<int>(const gmx_domdec_t* dd,
129 int ddDimensionIndex,
131 gmx::ArrayRef<int> sendBuffer,
132 gmx::ArrayRef<int> receiveBuffer);
134 //! Extern declaration for real specialization
135 extern template void ddSendrecv<real>(const gmx_domdec_t* dd,
136 int ddDimensionIndex,
138 gmx::ArrayRef<real> sendBuffer,
139 gmx::ArrayRef<real> receiveBuffer);
141 //! Extern declaration for gmx::RVec specialization
142 extern template void ddSendrecv<gmx::RVec>(const gmx_domdec_t* dd,
143 int ddDimensionIndex,
145 gmx::ArrayRef<gmx::RVec> sendBuffer,
146 gmx::ArrayRef<gmx::RVec> receiveBuffer);
148 /*! \brief Move revc's in the comm. region one cell along the domain decomposition
150 * Moves in dimension indexed by ddimind, simultaneously in the forward
151 * and backward directions.
153 void dd_sendrecv2_rvec(const struct gmx_domdec_t* dd,
165 /* The functions below perform the same operations as the MPI functions
166 * with the same name appendices, but over the domain decomposition
168 * The DD master node is the master for these operations.
171 /*! \brief Broadcasts \p nbytes from \p data on \p DDMASTERRANK to all PP ranks */
172 void dd_bcast(const gmx_domdec_t* dd, int nbytes, void* data);
174 /*! \brief Copies \p nbytes from \p src to \p dest on \p DDMASTERRANK
175 * and then broadcasts to \p dest on all PP ranks */
176 void dd_bcastc(const gmx_domdec_t* dd, int nbytes, void* src, void* dest);
178 /*! \brief Scatters \p nbytes from \p src on \p DDMASTERRANK to all PP ranks, received in \p dest */
179 void dd_scatter(const gmx_domdec_t* dd, int nbytes, const void* src, void* dest);
181 /*! \brief Gathers \p nbytes from \p src on all PP ranks, received in \p dest on \p DDMASTERRANK */
182 void dd_gather(const gmx_domdec_t* dd, int nbytes, const void* src, void* dest);
184 /*! \brief Scatters \p scounts bytes from \p src on \p DDMASTERRANK to all PP ranks, receiving \p rcount bytes in \p dest.
186 * See man MPI_Scatterv for details of how to construct scounts and disps.
187 * If rcount==0, rbuf is allowed to be NULL */
188 void dd_scatterv(const gmx_domdec_t* dd, int* scounts, int* disps, const void* sbuf, int rcount, void* rbuf);
190 /*! \brief Gathers \p rcount bytes from \p src on all PP ranks, received in \p scounts bytes in \p dest on \p DDMASTERRANK.
192 * See man MPI_Gatherv for details of how to construct scounts and disps.
194 * If scount==0, sbuf is allowed to be NULL */
195 void dd_gatherv(const gmx_domdec_t* dd, int scount, const void* sbuf, int* rcounts, int* disps, void* rbuf);