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/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;
85 enum class PbcType : int;
88 class GpuEventSynchronizer;
92 class DeviceStreamManager;
93 class ForceWithShiftForces;
95 class RangePartitioning;
96 class VirtualSitesHandler;
99 /*! \brief Returns the global topology atom number belonging to local atom index i.
101 * This function is intended for writing ASCII output
102 * and returns atom numbers starting at 1.
103 * When dd=NULL returns i+1.
105 int ddglatnr(const gmx_domdec_t* dd, int i);
107 /*! \brief Returns a list of update group partitioning for each molecule type or empty when update groups are not used */
108 gmx::ArrayRef<const gmx::RangePartitioning> getUpdateGroupingPerMoleculetype(const gmx_domdec_t& dd);
110 /*! \brief Store the global cg indices of the home cgs in state,
112 * This means it can be reset, even after a new DD partitioning.
114 void dd_store_state(struct gmx_domdec_t* dd, t_state* state);
116 /*! \brief Returns a pointer to the gmx_domdec_zones_t struct */
117 struct gmx_domdec_zones_t* domdec_zones(struct gmx_domdec_t* dd);
119 /*! \brief Returns the range for atoms in zones*/
120 int dd_numAtomsZones(const gmx_domdec_t& dd);
122 /*! \brief Returns the number of home atoms */
123 int dd_numHomeAtoms(const gmx_domdec_t& dd);
125 /*! \brief Returns the atom range in the local state for atoms that need to be present in mdatoms */
126 int dd_natoms_mdatoms(const gmx_domdec_t* dd);
128 /*! \brief Returns the atom range in the local state for atoms involved in virtual sites */
129 int dd_natoms_vsite(const gmx_domdec_t* dd);
131 /*! \brief Sets the atom range for atom in the local state for atoms received in constraints communication */
132 void dd_get_constraint_range(const gmx_domdec_t* dd, int* at_start, int* at_end);
134 /*! \libinternal \brief Struct for passing around the number of PME domains */
137 int x; //!< The number of PME domains along dimension x
138 int y; //!< The number of PME domains along dimension y
141 /*! \brief Returns the number of PME domains, can be called with dd=NULL */
142 NumPmeDomains getNumPmeDomains(const gmx_domdec_t* dd);
144 /*! \brief Returns the set of DD ranks that communicate with pme node cr->nodeid */
145 std::vector<int> get_pme_ddranks(const t_commrec* cr, int pmenodeid);
147 /*! \brief Returns the maximum shift for coordinate communication in PME, dim x */
148 int dd_pme_maxshift_x(const gmx_domdec_t* dd);
150 /*! \brief Returns the maximum shift for coordinate communication in PME, dim y */
151 int dd_pme_maxshift_y(const gmx_domdec_t* dd);
153 /*! \brief Return whether constraints, not including settles, cross domain boundaries */
154 bool ddHaveSplitConstraints(const gmx_domdec_t& dd);
156 /*! \brief Return whether update groups are used */
157 bool ddUsesUpdateGroups(const gmx_domdec_t& dd);
159 /*! \brief Initialize data structures for bonded interactions */
160 void dd_init_bondeds(FILE* fplog,
162 const gmx_mtop_t& mtop,
163 const gmx::VirtualSitesHandler* vsite,
164 const t_inputrec* ir,
166 gmx::ArrayRef<cginfo_mb_t> cginfo_mb);
168 /*! \brief Returns whether molecules are always whole, i.e. not broken by PBC */
169 bool dd_moleculesAreAlwaysWhole(const gmx_domdec_t& dd);
171 /*! \brief Returns if we need to do pbc for calculating bonded interactions */
172 gmx_bool dd_bonded_molpbc(const gmx_domdec_t* dd, PbcType pbcType);
174 /*! \brief Change the DD non-bonded communication cut-off.
176 * This could fail when trying to increase the cut-off,
177 * then FALSE will be returned and the cut-off is not modified.
179 * \param[in] cr Communication recrod
180 * \param[in] box Box matrix, used for computing the dimensions of the system
181 * \param[in] x Position vector, used for computing the dimensions of the system
182 * \param[in] cutoffRequested The requested atom to atom cut-off distance, usually the pair-list cutoff distance
184 gmx_bool change_dd_cutoff(t_commrec* cr, const matrix box, gmx::ArrayRef<const gmx::RVec> x, real cutoffRequested);
186 /*! \brief Set up communication for averaging GPU wait times over domains
188 * When domains (PP MPI ranks) share a GPU, the individual GPU wait times
189 * are meaningless, as it depends on the order in which tasks on the same
190 * GPU finish. Therefore there wait times need to be averaged over the ranks
191 * sharing the same GPU. This function sets up the communication for that.
193 void dd_setup_dlb_resource_sharing(const t_commrec* cr, int gpu_id);
195 /*! \brief Cycle counter indices used internally in the domain decomposition */
206 /*! \brief Add the wallcycle count to the DD counter */
207 void dd_cycles_add(const gmx_domdec_t* dd, float cycles, int ddCycl);
209 /*! \brief Communicate the coordinates to the neighboring cells and do pbc. */
210 void dd_move_x(struct gmx_domdec_t* dd, const matrix box, gmx::ArrayRef<gmx::RVec> x, gmx_wallcycle* wcycle);
212 /*! \brief Sum the forces over the neighboring cells.
214 * When fshift!=NULL the shift forces are updated to obtain
215 * the correct virial from the single sum including f.
217 void dd_move_f(struct gmx_domdec_t* dd, gmx::ForceWithShiftForces* forceWithShiftForces, gmx_wallcycle* wcycle);
219 /*! \brief Communicate a real for each atom to the neighboring cells. */
220 void dd_atom_spread_real(struct gmx_domdec_t* dd, real v[]);
222 /*! \brief Sum the contributions to a real for each atom over the neighboring cells. */
223 void dd_atom_sum_real(struct gmx_domdec_t* dd, real v[]);
225 /*! \brief Reset all the statistics and counters for total run counting */
226 void reset_dd_statistics_counters(struct gmx_domdec_t* dd);
228 /* In domdec_con.c */
230 /*! \brief Communicates the virtual site forces, reduces the shift forces when \p fshift != NULL */
231 void dd_move_f_vsites(const gmx_domdec_t& dd, gmx::ArrayRef<gmx::RVec> f, gmx::ArrayRef<gmx::RVec> fshift);
233 /*! \brief Clears the forces for virtual sites */
234 void dd_clear_f_vsites(const gmx_domdec_t& dd, gmx::ArrayRef<gmx::RVec> f);
236 /*! \brief Move x0 and also x1 if x1!=NULL. bX1IsCoord tells if to do PBC on x1 */
237 void dd_move_x_constraints(struct gmx_domdec_t* dd,
239 gmx::ArrayRef<gmx::RVec> x0,
240 gmx::ArrayRef<gmx::RVec> x1,
241 gmx_bool bX1IsCoord);
243 /*! \brief Communicates the coordinates involved in virtual sites */
244 void dd_move_x_vsites(const gmx_domdec_t& dd, const matrix box, rvec* x);
245 /*! \brief Communicates the positions and velocities involved in virtual sites */
246 void dd_move_x_and_v_vsites(const gmx_domdec_t& dd, const matrix box, rvec* x, rvec* v);
248 /*! \brief Returns the local atom count array for all constraints
250 * The local atom count for a constraint, possible values 2/1/0, is needed
251 * to avoid not/double-counting contributions linked to the Lagrange
252 * multiplier, such as the virial and free-energy derivatives.
254 * \note When \p dd = nullptr, an empty reference is returned.
256 gmx::ArrayRef<const int> dd_constraints_nlocalatoms(const gmx_domdec_t* dd);
258 /* In domdec_top.c */
260 /*! \brief Print error output when interactions are missing */
261 [[noreturn]] void dd_print_missing_interactions(const gmx::MDLogger& mdlog,
264 const gmx_mtop_t* top_global,
265 const gmx_localtop_t* top_local,
266 gmx::ArrayRef<const gmx::RVec> x,
269 /*! \brief Generate and store the reverse topology */
270 void dd_make_reverse_top(FILE* fplog,
272 const gmx_mtop_t* mtop,
273 const gmx::VirtualSitesHandler* vsite,
274 const t_inputrec* ir,
277 /*! \brief Generate the local topology and virtual site data */
278 void dd_make_local_top(struct gmx_domdec_t* dd,
279 struct gmx_domdec_zones_t* zones,
286 const gmx_mtop_t& top,
287 gmx_localtop_t* ltop);
289 /*! \brief Sort ltop->ilist when we are doing free energy. */
290 void dd_sort_local_top(gmx_domdec_t* dd, const t_mdatoms* mdatoms, gmx_localtop_t* ltop);
292 /*! \brief Construct local state */
293 void dd_init_local_state(struct gmx_domdec_t* dd, const t_state* state_global, t_state* local_state);
295 /*! \brief Generate a list of links between atoms that are linked by bonded interactions
297 * Also stores whether atoms are linked in \p cginfo_mb.
299 t_blocka* makeBondedLinks(const gmx_mtop_t& mtop, gmx::ArrayRef<cginfo_mb_t> cginfo_mb);
301 /*! \brief Calculate the maximum distance involved in 2-body and multi-body bonded interactions */
302 void dd_bonded_cg_distance(const gmx::MDLogger& mdlog,
303 const gmx_mtop_t* mtop,
304 const t_inputrec* ir,
305 gmx::ArrayRef<const gmx::RVec> x,
311 /*! \brief Construct the GPU halo exchange object(s).
313 * \param[in] mdlog The logger object.
314 * \param[in] cr The commrec object.
315 * \param[in] deviceStreamManager Manager of the GPU context and streams.
316 * \param[in] wcycle The wallclock counter.
318 void constructGpuHaloExchange(const gmx::MDLogger& mdlog,
320 const gmx::DeviceStreamManager& deviceStreamManager,
321 gmx_wallcycle* wcycle);
324 * (Re-) Initialization for GPU halo exchange
325 * \param [in] cr The commrec object
326 * \param [in] d_coordinatesBuffer pointer to coordinates buffer in GPU memory
327 * \param [in] d_forcesBuffer pointer to forces buffer in GPU memory
329 void reinitGpuHaloExchange(const t_commrec& cr,
330 DeviceBuffer<gmx::RVec> d_coordinatesBuffer,
331 DeviceBuffer<gmx::RVec> d_forcesBuffer);
334 /*! \brief GPU halo exchange of coordinates buffer.
335 * \param [in] cr The commrec object
336 * \param [in] box Coordinate box (from which shifts will be constructed)
337 * \param [in] coordinatesReadyOnDeviceEvent event recorded when coordinates have been copied to device
339 void communicateGpuHaloCoordinates(const t_commrec& cr,
341 GpuEventSynchronizer* coordinatesReadyOnDeviceEvent);
344 /*! \brief GPU halo exchange of force buffer.
345 * \param [in] cr The commrec object
346 * \param [in] accumulateForces True if forces should accumulate, otherwise they are set
348 void communicateGpuHaloForces(const t_commrec& cr, bool accumulateForces);