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
;
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 Return whether the DD has a single dimension
161 * The GPU halo exchange code requires a 1D DD, and its setup code can
162 * use the returned value to understand what it should do.
164 bool is1D(const gmx_domdec_t
& dd
);
166 /*! \brief Initialize data structures for bonded interactions */
167 void dd_init_bondeds(FILE* fplog
,
169 const gmx_mtop_t
& mtop
,
170 const gmx::VirtualSitesHandler
* vsite
,
171 const t_inputrec
* ir
,
173 gmx::ArrayRef
<cginfo_mb_t
> cginfo_mb
);
175 /*! \brief Returns whether molecules are always whole, i.e. not broken by PBC */
176 bool dd_moleculesAreAlwaysWhole(const gmx_domdec_t
& dd
);
178 /*! \brief Returns if we need to do pbc for calculating bonded interactions */
179 gmx_bool
dd_bonded_molpbc(const gmx_domdec_t
* dd
, PbcType pbcType
);
181 /*! \brief Change the DD non-bonded communication cut-off.
183 * This could fail when trying to increase the cut-off,
184 * then FALSE will be returned and the cut-off is not modified.
186 * \param[in] cr Communication recrod
187 * \param[in] box Box matrix, used for computing the dimensions of the system
188 * \param[in] x Position vector, used for computing the dimensions of the system
189 * \param[in] cutoffRequested The requested atom to atom cut-off distance, usually the pair-list cutoff distance
191 gmx_bool
change_dd_cutoff(t_commrec
* cr
, const matrix box
, gmx::ArrayRef
<const gmx::RVec
> x
, real cutoffRequested
);
193 /*! \brief Set up communication for averaging GPU wait times over domains
195 * When domains (PP MPI ranks) share a GPU, the individual GPU wait times
196 * are meaningless, as it depends on the order in which tasks on the same
197 * GPU finish. Therefore there wait times need to be averaged over the ranks
198 * sharing the same GPU. This function sets up the communication for that.
200 void dd_setup_dlb_resource_sharing(const t_commrec
* cr
, int gpu_id
);
202 /*! \brief Cycle counter indices used internally in the domain decomposition */
213 /*! \brief Add the wallcycle count to the DD counter */
214 void dd_cycles_add(const gmx_domdec_t
* dd
, float cycles
, int ddCycl
);
216 /*! \brief Communicate the coordinates to the neighboring cells and do pbc. */
217 void dd_move_x(struct gmx_domdec_t
* dd
, const matrix box
, gmx::ArrayRef
<gmx::RVec
> x
, gmx_wallcycle
* wcycle
);
219 /*! \brief Sum the forces over the neighboring cells.
221 * When fshift!=NULL the shift forces are updated to obtain
222 * the correct virial from the single sum including f.
224 void dd_move_f(struct gmx_domdec_t
* dd
, gmx::ForceWithShiftForces
* forceWithShiftForces
, gmx_wallcycle
* wcycle
);
226 /*! \brief Communicate a real for each atom to the neighboring cells. */
227 void dd_atom_spread_real(struct gmx_domdec_t
* dd
, real v
[]);
229 /*! \brief Sum the contributions to a real for each atom over the neighboring cells. */
230 void dd_atom_sum_real(struct gmx_domdec_t
* dd
, real v
[]);
232 /*! \brief Reset all the statistics and counters for total run counting */
233 void reset_dd_statistics_counters(struct gmx_domdec_t
* dd
);
235 /* In domdec_con.c */
237 /*! \brief Communicates the virtual site forces, reduces the shift forces when \p fshift != NULL */
238 void dd_move_f_vsites(const gmx_domdec_t
& dd
, gmx::ArrayRef
<gmx::RVec
> f
, gmx::ArrayRef
<gmx::RVec
> fshift
);
240 /*! \brief Clears the forces for virtual sites */
241 void dd_clear_f_vsites(const gmx_domdec_t
& dd
, gmx::ArrayRef
<gmx::RVec
> f
);
243 /*! \brief Move x0 and also x1 if x1!=NULL. bX1IsCoord tells if to do PBC on x1 */
244 void dd_move_x_constraints(struct gmx_domdec_t
* dd
,
246 gmx::ArrayRef
<gmx::RVec
> x0
,
247 gmx::ArrayRef
<gmx::RVec
> x1
,
248 gmx_bool bX1IsCoord
);
250 /*! \brief Communicates the coordinates involved in virtual sites */
251 void dd_move_x_vsites(const gmx_domdec_t
& dd
, const matrix box
, rvec
* x
);
253 /*! \brief Returns the local atom count array for all constraints
255 * The local atom count for a constraint, possible values 2/1/0, is needed
256 * to avoid not/double-counting contributions linked to the Lagrange
257 * multiplier, such as the virial and free-energy derivatives.
259 * \note When \p dd = nullptr, an empty reference is returned.
261 gmx::ArrayRef
<const int> dd_constraints_nlocalatoms(const gmx_domdec_t
* dd
);
263 /* In domdec_top.c */
265 /*! \brief Print error output when interactions are missing */
266 [[noreturn
]] void dd_print_missing_interactions(const gmx::MDLogger
& mdlog
,
269 const gmx_mtop_t
* top_global
,
270 const gmx_localtop_t
* top_local
,
271 gmx::ArrayRef
<const gmx::RVec
> x
,
274 /*! \brief Generate and store the reverse topology */
275 void dd_make_reverse_top(FILE* fplog
,
277 const gmx_mtop_t
* mtop
,
278 const gmx::VirtualSitesHandler
* vsite
,
279 const t_inputrec
* ir
,
282 /*! \brief Generate the local topology and virtual site data */
283 void dd_make_local_top(struct gmx_domdec_t
* dd
,
284 struct gmx_domdec_zones_t
* zones
,
291 const gmx_mtop_t
& top
,
292 gmx_localtop_t
* ltop
);
294 /*! \brief Sort ltop->ilist when we are doing free energy. */
295 void dd_sort_local_top(gmx_domdec_t
* dd
, const t_mdatoms
* mdatoms
, gmx_localtop_t
* ltop
);
297 /*! \brief Construct local state */
298 void dd_init_local_state(struct gmx_domdec_t
* dd
, const t_state
* state_global
, t_state
* local_state
);
300 /*! \brief Generate a list of links between atoms that are linked by bonded interactions
302 * Also stores whether atoms are linked in \p cginfo_mb.
304 t_blocka
* makeBondedLinks(const gmx_mtop_t
& mtop
, gmx::ArrayRef
<cginfo_mb_t
> cginfo_mb
);
306 /*! \brief Calculate the maximum distance involved in 2-body and multi-body bonded interactions */
307 void dd_bonded_cg_distance(const gmx::MDLogger
& mdlog
,
308 const gmx_mtop_t
* mtop
,
309 const t_inputrec
* ir
,
310 gmx::ArrayRef
<const gmx::RVec
> x
,
316 /*! \brief Construct the GPU halo exchange object(s).
318 * \param[in] mdlog The logger object.
319 * \param[in] cr The commrec object.
320 * \param[in] deviceStreamManager Manager of the GPU context and streams.
321 * \param[in] wcycle The wallclock counter.
323 void constructGpuHaloExchange(const gmx::MDLogger
& mdlog
,
325 const gmx::DeviceStreamManager
& deviceStreamManager
,
326 gmx_wallcycle
* wcycle
);
329 * (Re-) Initialization for GPU halo exchange
330 * \param [in] cr The commrec object
331 * \param [in] d_coordinatesBuffer pointer to coordinates buffer in GPU memory
332 * \param [in] d_forcesBuffer pointer to forces buffer in GPU memory
334 void reinitGpuHaloExchange(const t_commrec
& cr
,
335 DeviceBuffer
<gmx::RVec
> d_coordinatesBuffer
,
336 DeviceBuffer
<gmx::RVec
> d_forcesBuffer
);
339 /*! \brief GPU halo exchange of coordinates buffer.
340 * \param [in] cr The commrec object
341 * \param [in] box Coordinate box (from which shifts will be constructed)
342 * \param [in] coordinatesReadyOnDeviceEvent event recorded when coordinates have been copied to device
344 void communicateGpuHaloCoordinates(const t_commrec
& cr
,
346 GpuEventSynchronizer
* coordinatesReadyOnDeviceEvent
);
349 /*! \brief GPU halo exchange of force buffer.
350 * \param [in] cr The commrec object
351 * \param [in] accumulateForces True if forces should accumulate, otherwise they are set
353 void communicateGpuHaloForces(const t_commrec
& cr
, bool accumulateForces
);