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,2021, 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/real.h"
72 struct gmx_domdec_zones_t;
73 struct gmx_localtop_t;
83 enum class PbcType : int;
86 class GpuEventSynchronizer;
90 class DeviceStreamManager;
91 class ForceWithShiftForces;
93 class RangePartitioning;
94 class VirtualSitesHandler;
97 enum class DDBondedChecking : bool;
100 /*! \brief Returns the global topology atom number belonging to local atom index i.
102 * This function is intended for writing ASCII output
103 * and returns atom numbers starting at 1.
104 * When dd=NULL returns i+1.
106 int ddglatnr(const gmx_domdec_t* dd, int i);
108 /*! \brief Returns a list of update group partitioning for each molecule type or empty when update groups are not used */
109 gmx::ArrayRef<const gmx::RangePartitioning> getUpdateGroupingsPerMoleculeType(const gmx_domdec_t& dd);
111 /*! \brief Store the global cg indices of the home cgs in state,
113 * This means it can be reset, even after a new DD partitioning.
115 void dd_store_state(const gmx_domdec_t& dd, t_state* state);
117 /*! \brief Returns a pointer to the gmx_domdec_zones_t struct */
118 struct gmx_domdec_zones_t* domdec_zones(struct gmx_domdec_t* dd);
120 /*! \brief Returns the range for atoms in zones*/
121 int dd_numAtomsZones(const gmx_domdec_t& dd);
123 /*! \brief Returns the number of home atoms */
124 int dd_numHomeAtoms(const gmx_domdec_t& dd);
126 /*! \brief Returns the atom range in the local state for atoms that need to be present in mdatoms */
127 int dd_natoms_mdatoms(const gmx_domdec_t& dd);
129 /*! \brief Returns the atom range in the local state for atoms involved in virtual sites */
130 int dd_natoms_vsite(const gmx_domdec_t& dd);
132 /*! \brief Sets the atom range for atom in the local state for atoms received in constraints communication */
133 void dd_get_constraint_range(const gmx_domdec_t& dd, int* at_start, int* at_end);
135 /*! \libinternal \brief Struct for passing around the number of PME domains */
138 int x; //!< The number of PME domains along dimension x
139 int y; //!< The number of PME domains along dimension y
142 /*! \brief Returns the number of PME domains, can be called with dd=NULL */
143 NumPmeDomains getNumPmeDomains(const gmx_domdec_t* dd);
145 /*! \brief Returns the set of DD ranks that communicate with pme node cr->nodeid */
146 std::vector<int> get_pme_ddranks(const t_commrec* cr, int pmenodeid);
148 /*! \brief Returns the maximum shift for coordinate communication in PME, dim x */
149 int dd_pme_maxshift_x(const gmx_domdec_t& dd);
151 /*! \brief Returns the maximum shift for coordinate communication in PME, dim y */
152 int dd_pme_maxshift_y(const gmx_domdec_t& dd);
154 /*! \brief Return whether constraints, not including settles, may cross domain boundaries */
155 bool ddMayHaveSplitConstraints(const gmx_domdec_t& dd);
157 /*! \brief Return whether update groups are used */
158 bool ddUsesUpdateGroups(const gmx_domdec_t& dd);
160 /*! \brief Initialize data structures for bonded interactions */
161 void dd_init_bondeds(FILE* fplog,
163 const gmx_mtop_t& mtop,
164 const gmx::VirtualSitesHandler* vsite,
165 const t_inputrec& inputrec,
166 gmx::DDBondedChecking ddBondedChecking,
167 gmx::ArrayRef<cginfo_mb_t> cginfo_mb);
169 /*! \brief Returns whether molecules are always whole, i.e. not broken by PBC */
170 bool dd_moleculesAreAlwaysWhole(const gmx_domdec_t& dd);
172 /*! \brief Returns if we need to do pbc for calculating bonded interactions */
173 bool dd_bonded_molpbc(const gmx_domdec_t& dd, PbcType pbcType);
175 /*! \brief Change the DD non-bonded communication cut-off.
177 * This could fail when trying to increase the cut-off,
178 * then FALSE will be returned and the cut-off is not modified.
180 * \param[in] cr Communication recrod
181 * \param[in] box Box matrix, used for computing the dimensions of the system
182 * \param[in] x Position vector, used for computing the dimensions of the system
183 * \param[in] cutoffRequested The requested atom to atom cut-off distance, usually the pair-list cutoff distance
185 bool change_dd_cutoff(t_commrec* cr, const matrix box, gmx::ArrayRef<const gmx::RVec> x, real cutoffRequested);
187 /*! \brief Set up communication for averaging GPU wait times over domains
189 * When domains (PP MPI ranks) share a GPU, the individual GPU wait times
190 * are meaningless, as it depends on the order in which tasks on the same
191 * GPU finish. Therefore there wait times need to be averaged over the ranks
192 * sharing the same GPU. This function sets up the communication for that.
194 void dd_setup_dlb_resource_sharing(const t_commrec* cr, int gpu_id);
196 /*! \brief Cycle counter indices used internally in the domain decomposition */
207 /*! \brief Add the wallcycle count to the DD counter */
208 void dd_cycles_add(const gmx_domdec_t* dd, float cycles, int ddCycl);
210 /*! \brief Communicate the coordinates to the neighboring cells and do pbc. */
211 void dd_move_x(struct gmx_domdec_t* dd, const matrix box, gmx::ArrayRef<gmx::RVec> x, gmx_wallcycle* wcycle);
213 /*! \brief Sum the forces over the neighboring cells.
215 * When fshift!=NULL the shift forces are updated to obtain
216 * the correct virial from the single sum including f.
218 void dd_move_f(struct gmx_domdec_t* dd, gmx::ForceWithShiftForces* forceWithShiftForces, gmx_wallcycle* wcycle);
220 /*! \brief Communicate a real for each atom to the neighboring cells. */
221 void dd_atom_spread_real(struct gmx_domdec_t* dd, real v[]);
223 /*! \brief Sum the contributions to a real for each atom over the neighboring cells. */
224 void dd_atom_sum_real(struct gmx_domdec_t* dd, real v[]);
226 /*! \brief Reset all the statistics and counters for total run counting */
227 void reset_dd_statistics_counters(struct gmx_domdec_t* dd);
229 /* In domdec_con.c */
231 /*! \brief Communicates the virtual site forces, reduces the shift forces when \p fshift != NULL */
232 void dd_move_f_vsites(const gmx_domdec_t& dd, gmx::ArrayRef<gmx::RVec> f, gmx::ArrayRef<gmx::RVec> fshift);
234 /*! \brief Clears the forces for virtual sites */
235 void dd_clear_f_vsites(const gmx_domdec_t& dd, gmx::ArrayRef<gmx::RVec> f);
237 /*! \brief Move x0 and also x1 if x1!=NULL. bX1IsCoord tells if to do PBC on x1 */
238 void dd_move_x_constraints(struct gmx_domdec_t* dd,
240 gmx::ArrayRef<gmx::RVec> x0,
241 gmx::ArrayRef<gmx::RVec> x1,
244 /*! \brief Communicates the coordinates involved in virtual sites */
245 void dd_move_x_vsites(const gmx_domdec_t& dd, const matrix box, rvec* x);
246 /*! \brief Communicates the positions and velocities involved in virtual sites */
247 void dd_move_x_and_v_vsites(const gmx_domdec_t& dd, const matrix box, rvec* x, rvec* v);
249 /*! \brief Returns the local atom count array for all constraints
251 * The local atom count for a constraint, possible values 2/1/0, is needed
252 * to avoid not/double-counting contributions linked to the Lagrange
253 * multiplier, such as the virial and free-energy derivatives.
255 * \note When \p dd = nullptr, an empty reference is returned.
257 gmx::ArrayRef<const int> dd_constraints_nlocalatoms(const gmx_domdec_t* dd);
259 /* In domdec_top.c */
261 /*! \brief Print error output when interactions are missing */
262 [[noreturn]] void dd_print_missing_interactions(const gmx::MDLogger& mdlog,
265 const gmx_mtop_t& top_global,
266 const gmx_localtop_t* top_local,
267 gmx::ArrayRef<const gmx::RVec> x,
270 /*! \brief Generate and store the reverse topology */
271 void dd_make_reverse_top(FILE* fplog,
273 const gmx_mtop_t& mtop,
274 const gmx::VirtualSitesHandler* vsite,
275 const t_inputrec& inputrec,
276 gmx::DDBondedChecking ddBondedChecking);
278 /*! \brief Return whether the total bonded interaction count across
279 * domains should be checked this step. */
280 bool shouldCheckNumberOfBondedInteractions(const gmx_domdec_t& dd);
282 //! Return the number of bonded interactions in this domain.
283 int numBondedInteractions(const gmx_domdec_t& dd);
285 /*! \brief Set total bonded interaction count across domains. */
286 void setNumberOfBondedInteractionsOverAllDomains(const gmx_domdec_t& dd, int newValue);
288 /*! \brief Check whether bonded interactions are missing from the reverse topology
289 * produced by domain decomposition.
291 * Must only be called when DD is active.
293 * \param[in] mdlog Logger
294 * \param[in] cr Communication object
295 * \param[in] top_global Global topology for the error message
296 * \param[in] top_local Local topology for the error message
297 * \param[in] x Position vector for the error message
298 * \param[in] box Box matrix for the error message
300 void checkNumberOfBondedInteractions(const gmx::MDLogger& mdlog,
302 const gmx_mtop_t& top_global,
303 const gmx_localtop_t* top_local,
304 gmx::ArrayRef<const gmx::RVec> x,
307 /*! \brief Generate the local topology and virtual site data */
308 void dd_make_local_top(struct gmx_domdec_t* dd,
309 struct gmx_domdec_zones_t* zones,
315 gmx::ArrayRef<const gmx::RVec> coordinates,
316 const gmx_mtop_t& top,
317 gmx_localtop_t* ltop);
319 /*! \brief Sort ltop->ilist when we are doing free energy. */
320 void dd_sort_local_top(const gmx_domdec_t& dd, const t_mdatoms* mdatoms, gmx_localtop_t* ltop);
322 /*! \brief Construct local state */
323 void dd_init_local_state(const gmx_domdec_t& dd, const t_state* state_global, t_state* local_state);
325 /*! \brief Generate a list of links between atoms that are linked by bonded interactions
327 * Also stores whether atoms are linked in \p cginfo_mb.
329 t_blocka* makeBondedLinks(const gmx_mtop_t& mtop, gmx::ArrayRef<cginfo_mb_t> cginfo_mb);
331 /*! \brief Calculate the maximum distance involved in 2-body and multi-body bonded interactions */
332 void dd_bonded_cg_distance(const gmx::MDLogger& mdlog,
333 const gmx_mtop_t& mtop,
334 const t_inputrec& ir,
335 gmx::ArrayRef<const gmx::RVec> x,
337 gmx::DDBondedChecking ddBondedChecking,
341 /*! \brief Construct the GPU halo exchange object(s).
343 * \param[in] mdlog The logger object.
344 * \param[in] cr The commrec object.
345 * \param[in] deviceStreamManager Manager of the GPU context and streams.
346 * \param[in] wcycle The wallclock counter.
348 void constructGpuHaloExchange(const gmx::MDLogger& mdlog,
350 const gmx::DeviceStreamManager& deviceStreamManager,
351 gmx_wallcycle* wcycle);
354 * (Re-) Initialization for GPU halo exchange
355 * \param [in] cr The commrec object
356 * \param [in] d_coordinatesBuffer pointer to coordinates buffer in GPU memory
357 * \param [in] d_forcesBuffer pointer to forces buffer in GPU memory
359 void reinitGpuHaloExchange(const t_commrec& cr,
360 DeviceBuffer<gmx::RVec> d_coordinatesBuffer,
361 DeviceBuffer<gmx::RVec> d_forcesBuffer);
364 /*! \brief GPU halo exchange of coordinates buffer.
365 * \param [in] cr The commrec object
366 * \param [in] box Coordinate box (from which shifts will be constructed)
367 * \param [in] coordinatesReadyOnDeviceEvent event recorded when coordinates have been copied to device
369 void communicateGpuHaloCoordinates(const t_commrec& cr,
371 GpuEventSynchronizer* coordinatesReadyOnDeviceEvent);
374 /*! \brief GPU halo exchange of force buffer.
375 * \param [in] cr The commrec object
376 * \param [in] accumulateForces True if forces should accumulate, otherwise they are set
378 void communicateGpuHaloForces(const t_commrec& cr, bool accumulateForces);