Ensure restart with update groups always works
[alexxy/gromacs.git] / src / gromacs / domdec / domdec.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 /*! \defgroup module_domdec Spatial domain decomposition (for parallelization over MPI)
38  * \ingroup group_mdrun
39  *
40  * \brief Manages the decomposition of the simulation volume over MPI
41  * ranks to try to distribute work evenly with minimal communication
42  * overheads.
43  *
44  * \todo Get domdec stuff out of mdtypes/commrec.h
45  *
46  * \author Berk Hess <hess@kth.se>
47  *
48  */
49
50 /*! \libinternal \file
51  *
52  * \brief This file declares functions for mdrun to call to manage the
53  * details of its domain decomposition.
54  *
55  * \author Berk Hess <hess@kth.se>
56  * \inlibraryapi
57  * \ingroup module_domdec
58  */
59
60 #ifndef GMX_DOMDEC_DOMDEC_H
61 #define GMX_DOMDEC_DOMDEC_H
62
63 #include <vector>
64
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"
70
71 struct cginfo_mb_t;
72 struct gmx_domdec_t;
73 struct gmx_ddbox_t;
74 struct gmx_domdec_zones_t;
75 struct gmx_localtop_t;
76 struct gmx_mtop_t;
77 struct t_block;
78 struct t_blocka;
79 struct t_commrec;
80 struct t_forcerec;
81 struct t_inputrec;
82 struct t_mdatoms;
83 struct t_nrnb;
84 struct gmx_wallcycle;
85 enum class PbcType : int;
86 class t_state;
87 class DeviceContext;
88 class GpuEventSynchronizer;
89
90 namespace gmx
91 {
92 class DeviceStreamManager;
93 class ForceWithShiftForces;
94 class MDLogger;
95 class RangePartitioning;
96 class VirtualSitesHandler;
97 } // namespace gmx
98
99 /*! \brief Returns the global topology atom number belonging to local atom index i.
100  *
101  * This function is intended for writing ASCII output
102  * and returns atom numbers starting at 1.
103  * When dd=NULL returns i+1.
104  */
105 int ddglatnr(const gmx_domdec_t* dd, int i);
106
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);
109
110 /*! \brief Store the global cg indices of the home cgs in state,
111  *
112  * This means it can be reset, even after a new DD partitioning.
113  */
114 void dd_store_state(struct gmx_domdec_t* dd, t_state* state);
115
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);
118
119 /*! \brief Returns the range for atoms in zones*/
120 int dd_numAtomsZones(const gmx_domdec_t& dd);
121
122 /*! \brief Returns the number of home atoms */
123 int dd_numHomeAtoms(const gmx_domdec_t& dd);
124
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);
127
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);
130
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);
133
134 /*! \libinternal \brief Struct for passing around the number of PME domains */
135 struct NumPmeDomains
136 {
137     int x; //!< The number of PME domains along dimension x
138     int y; //!< The number of PME domains along dimension y
139 };
140
141 /*! \brief Returns the number of PME domains, can be called with dd=NULL */
142 NumPmeDomains getNumPmeDomains(const gmx_domdec_t* dd);
143
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);
146
147 /*! \brief Returns the maximum shift for coordinate communication in PME, dim x */
148 int dd_pme_maxshift_x(const gmx_domdec_t* dd);
149
150 /*! \brief Returns the maximum shift for coordinate communication in PME, dim y */
151 int dd_pme_maxshift_y(const gmx_domdec_t* dd);
152
153 /*! \brief Return whether constraints, not including settles, cross domain boundaries */
154 bool ddHaveSplitConstraints(const gmx_domdec_t& dd);
155
156 /*! \brief Return whether update groups are used */
157 bool ddUsesUpdateGroups(const gmx_domdec_t& dd);
158
159 /*! \brief Initialize data structures for bonded interactions */
160 void dd_init_bondeds(FILE*                           fplog,
161                      gmx_domdec_t*                   dd,
162                      const gmx_mtop_t&               mtop,
163                      const gmx::VirtualSitesHandler* vsite,
164                      const t_inputrec*               ir,
165                      gmx_bool                        bBCheck,
166                      gmx::ArrayRef<cginfo_mb_t>      cginfo_mb);
167
168 /*! \brief Returns whether molecules are always whole, i.e. not broken by PBC */
169 bool dd_moleculesAreAlwaysWhole(const gmx_domdec_t& dd);
170
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);
173
174 /*! \brief Change the DD non-bonded communication cut-off.
175  *
176  * This could fail when trying to increase the cut-off,
177  * then FALSE will be returned and the cut-off is not modified.
178  *
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
183  */
184 gmx_bool change_dd_cutoff(t_commrec* cr, const matrix box, gmx::ArrayRef<const gmx::RVec> x, real cutoffRequested);
185
186 /*! \brief Set up communication for averaging GPU wait times over domains
187  *
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.
192  */
193 void dd_setup_dlb_resource_sharing(const t_commrec* cr, int gpu_id);
194
195 /*! \brief Cycle counter indices used internally in the domain decomposition */
196 enum
197 {
198     ddCyclStep,
199     ddCyclPPduringPME,
200     ddCyclF,
201     ddCyclWaitGPU,
202     ddCyclPME,
203     ddCyclNr
204 };
205
206 /*! \brief Add the wallcycle count to the DD counter */
207 void dd_cycles_add(const gmx_domdec_t* dd, float cycles, int ddCycl);
208
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);
211
212 /*! \brief Sum the forces over the neighboring cells.
213  *
214  * When fshift!=NULL the shift forces are updated to obtain
215  * the correct virial from the single sum including f.
216  */
217 void dd_move_f(struct gmx_domdec_t* dd, gmx::ForceWithShiftForces* forceWithShiftForces, gmx_wallcycle* wcycle);
218
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[]);
221
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[]);
224
225 /*! \brief Reset all the statistics and counters for total run counting */
226 void reset_dd_statistics_counters(struct gmx_domdec_t* dd);
227
228 /* In domdec_con.c */
229
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);
232
233 /*! \brief Clears the forces for virtual sites */
234 void dd_clear_f_vsites(const gmx_domdec_t& dd, gmx::ArrayRef<gmx::RVec> f);
235
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,
238                            const matrix             box,
239                            gmx::ArrayRef<gmx::RVec> x0,
240                            gmx::ArrayRef<gmx::RVec> x1,
241                            gmx_bool                 bX1IsCoord);
242
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
246 /*! \brief Returns the local atom count array for all constraints
247  *
248  * The local atom count for a constraint, possible values 2/1/0, is needed
249  * to avoid not/double-counting contributions linked to the Lagrange
250  * multiplier, such as the virial and free-energy derivatives.
251  *
252  * \note When \p dd = nullptr, an empty reference is returned.
253  */
254 gmx::ArrayRef<const int> dd_constraints_nlocalatoms(const gmx_domdec_t* dd);
255
256 /* In domdec_top.c */
257
258 /*! \brief Print error output when interactions are missing */
259 [[noreturn]] void dd_print_missing_interactions(const gmx::MDLogger&           mdlog,
260                                                 t_commrec*                     cr,
261                                                 int                            local_count,
262                                                 const gmx_mtop_t*              top_global,
263                                                 const gmx_localtop_t*          top_local,
264                                                 gmx::ArrayRef<const gmx::RVec> x,
265                                                 const matrix                   box);
266
267 /*! \brief Generate and store the reverse topology */
268 void dd_make_reverse_top(FILE*                           fplog,
269                          gmx_domdec_t*                   dd,
270                          const gmx_mtop_t*               mtop,
271                          const gmx::VirtualSitesHandler* vsite,
272                          const t_inputrec*               ir,
273                          gmx_bool                        bBCheck);
274
275 /*! \brief Generate the local topology and virtual site data */
276 void dd_make_local_top(struct gmx_domdec_t*       dd,
277                        struct gmx_domdec_zones_t* zones,
278                        int                        npbcdim,
279                        matrix                     box,
280                        rvec                       cellsize_min,
281                        const ivec                 npulse,
282                        t_forcerec*                fr,
283                        rvec*                      cgcm_or_x,
284                        const gmx_mtop_t&          top,
285                        gmx_localtop_t*            ltop);
286
287 /*! \brief Sort ltop->ilist when we are doing free energy. */
288 void dd_sort_local_top(gmx_domdec_t* dd, const t_mdatoms* mdatoms, gmx_localtop_t* ltop);
289
290 /*! \brief Construct local state */
291 void dd_init_local_state(struct gmx_domdec_t* dd, const t_state* state_global, t_state* local_state);
292
293 /*! \brief Generate a list of links between atoms that are linked by bonded interactions
294  *
295  * Also stores whether atoms are linked in \p cginfo_mb.
296  */
297 t_blocka* makeBondedLinks(const gmx_mtop_t& mtop, gmx::ArrayRef<cginfo_mb_t> cginfo_mb);
298
299 /*! \brief Calculate the maximum distance involved in 2-body and multi-body bonded interactions */
300 void dd_bonded_cg_distance(const gmx::MDLogger&           mdlog,
301                            const gmx_mtop_t*              mtop,
302                            const t_inputrec*              ir,
303                            gmx::ArrayRef<const gmx::RVec> x,
304                            const matrix                   box,
305                            gmx_bool                       bBCheck,
306                            real*                          r_2b,
307                            real*                          r_mb);
308
309 /*! \brief Construct the GPU halo exchange object(s).
310  *
311  * \param[in] mdlog               The logger object.
312  * \param[in] cr                  The commrec object.
313  * \param[in] deviceStreamManager Manager of the GPU context and streams.
314  * \param[in] wcycle              The wallclock counter.
315  */
316 void constructGpuHaloExchange(const gmx::MDLogger&            mdlog,
317                               const t_commrec&                cr,
318                               const gmx::DeviceStreamManager& deviceStreamManager,
319                               gmx_wallcycle*                  wcycle);
320
321 /*! \brief
322  * (Re-) Initialization for GPU halo exchange
323  * \param [in] cr                   The commrec object
324  * \param [in] d_coordinatesBuffer  pointer to coordinates buffer in GPU memory
325  * \param [in] d_forcesBuffer       pointer to forces buffer in GPU memory
326  */
327 void reinitGpuHaloExchange(const t_commrec&        cr,
328                            DeviceBuffer<gmx::RVec> d_coordinatesBuffer,
329                            DeviceBuffer<gmx::RVec> d_forcesBuffer);
330
331
332 /*! \brief GPU halo exchange of coordinates buffer.
333  * \param [in] cr                             The commrec object
334  * \param [in] box                            Coordinate box (from which shifts will be constructed)
335  * \param [in] coordinatesReadyOnDeviceEvent  event recorded when coordinates have been copied to device
336  */
337 void communicateGpuHaloCoordinates(const t_commrec&      cr,
338                                    const matrix          box,
339                                    GpuEventSynchronizer* coordinatesReadyOnDeviceEvent);
340
341
342 /*! \brief GPU halo exchange of force buffer.
343  * \param [in] cr                The commrec object
344  * \param [in] accumulateForces  True if forces should accumulate, otherwise they are set
345  */
346 void communicateGpuHaloForces(const t_commrec& cr, bool accumulateForces);
347
348 /*! \brief Wraps the \c positions so that atoms from the same
349  * update group share the same periodic image wrt \c box.
350  *
351  * When DD and update groups are in use, the simulation master rank
352  * should call this to ensure that e.g. when restarting a simulation
353  * that did not use update groups that the coordinates satisfy the new
354  * requirements.
355  *
356  * This function can probably be removed when even single-rank
357  * simulations use domain decomposition, because then the choice of
358  * whether update groups are used is probably going to be the same
359  * regardless of the rank count.
360  *
361  * \param[in]    dd         The DD manager
362  * \param[in]    mtop       The system topology
363  * \param[in]    box        The global system box
364  * \param[in]    positions  The global system positions
365  */
366 void putUpdateGroupAtomsInSamePeriodicImage(const gmx_domdec_t&      dd,
367                                             const gmx_mtop_t&        mtop,
368                                             const matrix             box,
369                                             gmx::ArrayRef<gmx::RVec> positions);
370
371 #endif