2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2005 - 2014, The GROMACS development team.
5 * Copyright (c) 2015,2016,2017,2018,2019 by the GROMACS development team.
6 * Copyright (c) 2020, 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.
37 /*! \defgroup module_domdec Spatial domain decomposition (for parallelization over MPI)
38 * \ingroup group_mdrun
40 * \brief Manages the decomposition of the simulation volume over MPI
41 * ranks to try to distribute work evenly with minimal communication
44 * \todo Get domdec stuff out of mdtypes/commrec.h
46 * \author Berk Hess <hess@kth.se>
50 /*! \libinternal \file
52 * \brief This file declares functions for mdrun to call to manage the
53 * details of its domain decomposition.
55 * \author Berk Hess <hess@kth.se>
57 * \ingroup module_domdec
60 #ifndef GMX_DOMDEC_DOMDEC_H
61 #define GMX_DOMDEC_DOMDEC_H
65 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
66 #include "gromacs/math/vectypes.h"
67 #include "gromacs/utility/arrayref.h"
68 #include "gromacs/utility/basedefinitions.h"
69 #include "gromacs/utility/real.h"
74 struct gmx_domdec_zones_t;
75 struct gmx_localtop_t;
86 enum class PbcType : int;
88 class GpuEventSynchronizer;
92 class ForceWithShiftForces;
94 class RangePartitioning;
97 /*! \brief Returns the global topology atom number belonging to local atom index i.
99 * This function is intended for writing ASCII output
100 * and returns atom numbers starting at 1.
101 * When dd=NULL returns i+1.
103 int ddglatnr(const gmx_domdec_t* dd, int i);
105 /*! \brief Returns a list of update group partitioning for each molecule type or empty when update groups are not used */
106 gmx::ArrayRef<const gmx::RangePartitioning> getUpdateGroupingPerMoleculetype(const gmx_domdec_t& dd);
108 /*! \brief Store the global cg indices of the home cgs in state,
110 * This means it can be reset, even after a new DD partitioning.
112 void dd_store_state(struct gmx_domdec_t* dd, t_state* state);
114 /*! \brief Returns a pointer to the gmx_domdec_zones_t struct */
115 struct gmx_domdec_zones_t* domdec_zones(struct gmx_domdec_t* dd);
117 /*! \brief Returns the range for atoms in zones*/
118 int dd_numAtomsZones(const gmx_domdec_t& dd);
120 /*! \brief Returns the number of home atoms */
121 int dd_numHomeAtoms(const gmx_domdec_t& dd);
123 /*! \brief Returns the atom range in the local state for atoms that need to be present in mdatoms */
124 int dd_natoms_mdatoms(const gmx_domdec_t* dd);
126 /*! \brief Returns the atom range in the local state for atoms involved in virtual sites */
127 int dd_natoms_vsite(const gmx_domdec_t* dd);
129 /*! \brief Sets the atom range for atom in the local state for atoms received in constraints communication */
130 void dd_get_constraint_range(const gmx_domdec_t* dd, int* at_start, int* at_end);
132 /*! \libinternal \brief Struct for passing around the number of PME domains */
135 int x; //!< The number of PME domains along dimension x
136 int y; //!< The number of PME domains along dimension y
139 /*! \brief Returns the number of PME domains, can be called with dd=NULL */
140 NumPmeDomains getNumPmeDomains(const gmx_domdec_t* dd);
142 /*! \brief Returns the set of DD ranks that communicate with pme node cr->nodeid */
143 std::vector<int> get_pme_ddranks(const t_commrec* cr, int pmenodeid);
145 /*! \brief Returns the maximum shift for coordinate communication in PME, dim x */
146 int dd_pme_maxshift_x(const gmx_domdec_t* dd);
148 /*! \brief Returns the maximum shift for coordinate communication in PME, dim y */
149 int dd_pme_maxshift_y(const gmx_domdec_t* dd);
151 /*! \brief Return whether constraints, not including settles, cross domain boundaries */
152 bool ddHaveSplitConstraints(const gmx_domdec_t& dd);
154 /*! \brief Return whether update groups are used */
155 bool ddUsesUpdateGroups(const gmx_domdec_t& dd);
157 /*! \brief Return whether the DD has a single dimension
159 * The GPU halo exchange code requires a 1D DD, and its setup code can
160 * use the returned value to understand what it should do.
162 bool is1D(const gmx_domdec_t& dd);
164 /*! \brief Initialize data structures for bonded interactions */
165 void dd_init_bondeds(FILE* fplog,
167 const gmx_mtop_t& mtop,
168 const gmx_vsite_t* vsite,
169 const t_inputrec* ir,
171 gmx::ArrayRef<cginfo_mb_t> cginfo_mb);
173 /*! \brief Returns whether molecules are always whole, i.e. not broken by PBC */
174 bool dd_moleculesAreAlwaysWhole(const gmx_domdec_t& dd);
176 /*! \brief Returns if we need to do pbc for calculating bonded interactions */
177 gmx_bool dd_bonded_molpbc(const gmx_domdec_t* dd, PbcType pbcType);
179 /*! \brief Change the DD non-bonded communication cut-off.
181 * This could fail when trying to increase the cut-off,
182 * then FALSE will be returned and the cut-off is not modified.
184 * \param[in] cr Communication recrod
185 * \param[in] box Box matrix, used for computing the dimensions of the system
186 * \param[in] x Position vector, used for computing the dimensions of the system
187 * \param[in] cutoffRequested The requested atom to atom cut-off distance, usually the pair-list cutoff distance
189 gmx_bool change_dd_cutoff(t_commrec* cr, const matrix box, gmx::ArrayRef<const gmx::RVec> x, real cutoffRequested);
191 /*! \brief Set up communication for averaging GPU wait times over domains
193 * When domains (PP MPI ranks) share a GPU, the individual GPU wait times
194 * are meaningless, as it depends on the order in which tasks on the same
195 * GPU finish. Therefore there wait times need to be averaged over the ranks
196 * sharing the same GPU. This function sets up the communication for that.
198 void dd_setup_dlb_resource_sharing(const t_commrec* cr, int gpu_id);
200 /*! \brief Cycle counter indices used internally in the domain decomposition */
211 /*! \brief Add the wallcycle count to the DD counter */
212 void dd_cycles_add(const gmx_domdec_t* dd, float cycles, int ddCycl);
214 /*! \brief Communicate the coordinates to the neighboring cells and do pbc. */
215 void dd_move_x(struct gmx_domdec_t* dd, const matrix box, gmx::ArrayRef<gmx::RVec> x, gmx_wallcycle* wcycle);
217 /*! \brief Sum the forces over the neighboring cells.
219 * When fshift!=NULL the shift forces are updated to obtain
220 * the correct virial from the single sum including f.
222 void dd_move_f(struct gmx_domdec_t* dd, gmx::ForceWithShiftForces* forceWithShiftForces, gmx_wallcycle* wcycle);
224 /*! \brief Communicate a real for each atom to the neighboring cells. */
225 void dd_atom_spread_real(struct gmx_domdec_t* dd, real v[]);
227 /*! \brief Sum the contributions to a real for each atom over the neighboring cells. */
228 void dd_atom_sum_real(struct gmx_domdec_t* dd, real v[]);
230 /*! \brief Reset all the statistics and counters for total run counting */
231 void reset_dd_statistics_counters(struct gmx_domdec_t* dd);
233 /* In domdec_con.c */
235 /*! \brief Communicates the virtual site forces, reduces the shift forces when \p fshift != NULL */
236 void dd_move_f_vsites(struct gmx_domdec_t* dd, rvec* f, rvec* fshift);
238 /*! \brief Clears the forces for virtual sites */
239 void dd_clear_f_vsites(struct gmx_domdec_t* dd, rvec* f);
241 /*! \brief Move x0 and also x1 if x1!=NULL. bX1IsCoord tells if to do PBC on x1 */
242 void dd_move_x_constraints(struct gmx_domdec_t* dd,
244 gmx::ArrayRef<gmx::RVec> x0,
245 gmx::ArrayRef<gmx::RVec> x1,
246 gmx_bool bX1IsCoord);
248 /*! \brief Communicates the coordinates involved in virtual sites */
249 void dd_move_x_vsites(struct gmx_domdec_t* dd, const matrix box, rvec* x);
251 /*! \brief Returns the local atom count array for all constraints
253 * The local atom count for a constraint, possible values 2/1/0, is needed
254 * to avoid not/double-counting contributions linked to the Lagrange
255 * multiplier, such as the virial and free-energy derivatives.
257 * \note When \p dd = nullptr, an empty reference is returned.
259 gmx::ArrayRef<const int> dd_constraints_nlocalatoms(const gmx_domdec_t* dd);
261 /* In domdec_top.c */
263 /*! \brief Print error output when interactions are missing */
264 [[noreturn]] void dd_print_missing_interactions(const gmx::MDLogger& mdlog,
267 const gmx_mtop_t* top_global,
268 const gmx_localtop_t* top_local,
269 gmx::ArrayRef<const gmx::RVec> x,
272 /*! \brief Generate and store the reverse topology */
273 void dd_make_reverse_top(FILE* fplog,
275 const gmx_mtop_t* mtop,
276 const gmx_vsite_t* vsite,
277 const t_inputrec* ir,
280 /*! \brief Generate the local topology and virtual site data */
281 void dd_make_local_top(struct gmx_domdec_t* dd,
282 struct gmx_domdec_zones_t* zones,
289 const gmx_mtop_t& top,
290 gmx_localtop_t* ltop);
292 /*! \brief Sort ltop->ilist when we are doing free energy. */
293 void dd_sort_local_top(gmx_domdec_t* dd, const t_mdatoms* mdatoms, gmx_localtop_t* ltop);
295 /*! \brief Construct local state */
296 void dd_init_local_state(struct gmx_domdec_t* dd, const t_state* state_global, t_state* local_state);
298 /*! \brief Generate a list of links between atoms that are linked by bonded interactions
300 * Also stores whether atoms are linked in \p cginfo_mb.
302 t_blocka* makeBondedLinks(const gmx_mtop_t& mtop, gmx::ArrayRef<cginfo_mb_t> cginfo_mb);
304 /*! \brief Calculate the maximum distance involved in 2-body and multi-body bonded interactions */
305 void dd_bonded_cg_distance(const gmx::MDLogger& mdlog,
306 const gmx_mtop_t* mtop,
307 const t_inputrec* ir,
314 /*! \brief Construct the GPU halo exchange object(s)
315 * \param[in] mdlog The logger object
316 * \param[in] cr The commrec object
317 * \param[in] streamLocal The local GPU stream
318 * \param[in] streamNonLocal The non-local GPU stream
320 void constructGpuHaloExchange(const gmx::MDLogger& mdlog, const t_commrec& cr, void* streamLocal, void* streamNonLocal);
323 * (Re-) Initialization for GPU halo exchange
324 * \param [in] cr The commrec object
325 * \param [in] d_coordinatesBuffer pointer to coordinates buffer in GPU memory
326 * \param [in] d_forcesBuffer pointer to forces buffer in GPU memory
328 void reinitGpuHaloExchange(const t_commrec& cr,
329 DeviceBuffer<gmx::RVec> d_coordinatesBuffer,
330 DeviceBuffer<gmx::RVec> d_forcesBuffer);
333 /*! \brief GPU halo exchange of coordinates buffer.
334 * \param [in] cr The commrec object
335 * \param [in] box Coordinate box (from which shifts will be constructed)
336 * \param [in] coordinatesReadyOnDeviceEvent event recorded when coordinates have been copied to device
338 void communicateGpuHaloCoordinates(const t_commrec& cr,
340 GpuEventSynchronizer* coordinatesReadyOnDeviceEvent);
343 /*! \brief GPU halo exchange of force buffer.
344 * \param [in] cr The commrec object
345 * \param [in] accumulateForces True if forces should accumulate, otherwise they are set
347 void communicateGpuHaloForces(const t_commrec& cr, bool accumulateForces);