a62a0c91f244a55562f1de87514af44989a4e8d4
[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/real.h"
68
69 struct cginfo_mb_t;
70 struct gmx_domdec_t;
71 struct gmx_ddbox_t;
72 struct gmx_domdec_zones_t;
73 struct gmx_localtop_t;
74 struct gmx_mtop_t;
75 struct t_block;
76 struct t_blocka;
77 struct t_commrec;
78 struct t_forcerec;
79 struct t_inputrec;
80 struct t_mdatoms;
81 struct t_nrnb;
82 struct gmx_wallcycle;
83 enum class PbcType : int;
84 class t_state;
85 class DeviceContext;
86 class GpuEventSynchronizer;
87
88 namespace gmx
89 {
90 class DeviceStreamManager;
91 class ForceWithShiftForces;
92 class MDLogger;
93 class RangePartitioning;
94 class VirtualSitesHandler;
95 template<typename>
96 class ArrayRef;
97 enum class DDBondedChecking;
98 } // namespace gmx
99
100 /*! \brief Returns the global topology atom number belonging to local atom index i.
101  *
102  * This function is intended for writing ASCII output
103  * and returns atom numbers starting at 1.
104  * When dd=NULL returns i+1.
105  */
106 int ddglatnr(const gmx_domdec_t* dd, int i);
107
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);
110
111 /*! \brief Store the global cg indices of the home cgs in state,
112  *
113  * This means it can be reset, even after a new DD partitioning.
114  */
115 void dd_store_state(const gmx_domdec_t& dd, t_state* state);
116
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);
119
120 /*! \brief Returns the range for atoms in zones*/
121 int dd_numAtomsZones(const gmx_domdec_t& dd);
122
123 /*! \brief Returns the number of home atoms */
124 int dd_numHomeAtoms(const gmx_domdec_t& dd);
125
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);
128
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);
131
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);
134
135 /*! \libinternal \brief Struct for passing around the number of PME domains */
136 struct NumPmeDomains
137 {
138     int x; //!< The number of PME domains along dimension x
139     int y; //!< The number of PME domains along dimension y
140 };
141
142 /*! \brief Returns the number of PME domains, can be called with dd=NULL */
143 NumPmeDomains getNumPmeDomains(const gmx_domdec_t* dd);
144
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);
147
148 /*! \brief Returns the maximum shift for coordinate communication in PME, dim x */
149 int dd_pme_maxshift_x(const gmx_domdec_t& dd);
150
151 /*! \brief Returns the maximum shift for coordinate communication in PME, dim y */
152 int dd_pme_maxshift_y(const gmx_domdec_t& dd);
153
154 /*! \brief Return whether constraints, not including settles, cross domain boundaries */
155 bool ddHaveSplitConstraints(const gmx_domdec_t& dd);
156
157 /*! \brief Return whether update groups are used */
158 bool ddUsesUpdateGroups(const gmx_domdec_t& dd);
159
160 /*! \brief Initialize data structures for bonded interactions */
161 void dd_init_bondeds(FILE*                           fplog,
162                      gmx_domdec_t*                   dd,
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);
168
169 /*! \brief Returns whether molecules are always whole, i.e. not broken by PBC */
170 bool dd_moleculesAreAlwaysWhole(const gmx_domdec_t& dd);
171
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);
174
175 /*! \brief Change the DD non-bonded communication cut-off.
176  *
177  * This could fail when trying to increase the cut-off,
178  * then FALSE will be returned and the cut-off is not modified.
179  *
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
184  */
185 bool change_dd_cutoff(t_commrec* cr, const matrix box, gmx::ArrayRef<const gmx::RVec> x, real cutoffRequested);
186
187 /*! \brief Set up communication for averaging GPU wait times over domains
188  *
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.
193  */
194 void dd_setup_dlb_resource_sharing(const t_commrec* cr, int gpu_id);
195
196 /*! \brief Cycle counter indices used internally in the domain decomposition */
197 enum
198 {
199     ddCyclStep,
200     ddCyclPPduringPME,
201     ddCyclF,
202     ddCyclWaitGPU,
203     ddCyclPME,
204     ddCyclNr
205 };
206
207 /*! \brief Add the wallcycle count to the DD counter */
208 void dd_cycles_add(const gmx_domdec_t* dd, float cycles, int ddCycl);
209
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);
212
213 /*! \brief Sum the forces over the neighboring cells.
214  *
215  * When fshift!=NULL the shift forces are updated to obtain
216  * the correct virial from the single sum including f.
217  */
218 void dd_move_f(struct gmx_domdec_t* dd, gmx::ForceWithShiftForces* forceWithShiftForces, gmx_wallcycle* wcycle);
219
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[]);
222
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[]);
225
226 /*! \brief Reset all the statistics and counters for total run counting */
227 void reset_dd_statistics_counters(struct gmx_domdec_t* dd);
228
229 /* In domdec_con.c */
230
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);
233
234 /*! \brief Clears the forces for virtual sites */
235 void dd_clear_f_vsites(const gmx_domdec_t& dd, gmx::ArrayRef<gmx::RVec> f);
236
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,
239                            const matrix             box,
240                            gmx::ArrayRef<gmx::RVec> x0,
241                            gmx::ArrayRef<gmx::RVec> x1,
242                            bool                     bX1IsCoord);
243
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);
248
249 /*! \brief Returns the local atom count array for all constraints
250  *
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.
254  *
255  * \note When \p dd = nullptr, an empty reference is returned.
256  */
257 gmx::ArrayRef<const int> dd_constraints_nlocalatoms(const gmx_domdec_t* dd);
258
259 /* In domdec_top.c */
260
261 /*! \brief Print error output when interactions are missing */
262 [[noreturn]] void dd_print_missing_interactions(const gmx::MDLogger&           mdlog,
263                                                 t_commrec*                     cr,
264                                                 int                            local_count,
265                                                 const gmx_mtop_t&              top_global,
266                                                 const gmx_localtop_t*          top_local,
267                                                 gmx::ArrayRef<const gmx::RVec> x,
268                                                 const matrix                   box);
269
270 /*! \brief Generate and store the reverse topology */
271 void dd_make_reverse_top(FILE*                           fplog,
272                          gmx_domdec_t*                   dd,
273                          const gmx_mtop_t&               mtop,
274                          const gmx::VirtualSitesHandler* vsite,
275                          const t_inputrec&               inputrec,
276                          gmx::DDBondedChecking           ddBondedChecking);
277
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);
281
282 //! Return the number of bonded interactions in this domain.
283 int numBondedInteractions(const gmx_domdec_t& dd);
284
285 /*! \brief Set total bonded interaction count across domains. */
286 void setNumberOfBondedInteractionsOverAllDomains(const gmx_domdec_t& dd, int newValue);
287
288 /*! \brief Check whether bonded interactions are missing from the reverse topology
289  * produced by domain decomposition.
290  *
291  * Must only be called when DD is active.
292  *
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
299  */
300 void checkNumberOfBondedInteractions(const gmx::MDLogger&           mdlog,
301                                      t_commrec*                     cr,
302                                      const gmx_mtop_t&              top_global,
303                                      const gmx_localtop_t*          top_local,
304                                      gmx::ArrayRef<const gmx::RVec> x,
305                                      const matrix                   box);
306
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,
310                        int                        npbcdim,
311                        matrix                     box,
312                        rvec                       cellsize_min,
313                        const ivec                 npulse,
314                        t_forcerec*                fr,
315                        rvec*                      cgcm_or_x,
316                        const gmx_mtop_t&          top,
317                        gmx_localtop_t*            ltop);
318
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);
321
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);
324
325 /*! \brief Generate a list of links between atoms that are linked by bonded interactions
326  *
327  * Also stores whether atoms are linked in \p cginfo_mb.
328  */
329 t_blocka* makeBondedLinks(const gmx_mtop_t& mtop, gmx::ArrayRef<cginfo_mb_t> cginfo_mb);
330
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,
336                            const matrix                   box,
337                            gmx::DDBondedChecking          ddBondedChecking,
338                            real*                          r_2b,
339                            real*                          r_mb);
340
341 /*! \brief Construct the GPU halo exchange object(s).
342  *
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.
347  */
348 void constructGpuHaloExchange(const gmx::MDLogger&            mdlog,
349                               const t_commrec&                cr,
350                               const gmx::DeviceStreamManager& deviceStreamManager,
351                               gmx_wallcycle*                  wcycle);
352
353 /*! \brief
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
358  */
359 void reinitGpuHaloExchange(const t_commrec&        cr,
360                            DeviceBuffer<gmx::RVec> d_coordinatesBuffer,
361                            DeviceBuffer<gmx::RVec> d_forcesBuffer);
362
363
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
368  */
369 void communicateGpuHaloCoordinates(const t_commrec&      cr,
370                                    const matrix          box,
371                                    GpuEventSynchronizer* coordinatesReadyOnDeviceEvent);
372
373
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
377  */
378 void communicateGpuHaloForces(const t_commrec& cr, bool accumulateForces);
379
380 #endif