0a7aa3202ef0c6eb57957a3540f95671bb96592a
[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, 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 gmx_vsite_t;
78 struct t_block;
79 struct t_blocka;
80 struct t_commrec;
81 struct t_forcerec;
82 struct t_inputrec;
83 struct t_mdatoms;
84 struct t_nrnb;
85 struct gmx_wallcycle;
86 enum class PbcType : int;
87 class t_state;
88 class DeviceContext;
89 class GpuEventSynchronizer;
90
91 namespace gmx
92 {
93 class ForceWithShiftForces;
94 class MDLogger;
95 class RangePartitioning;
96 } // namespace gmx
97
98 /*! \brief Returns the global topology atom number belonging to local atom index i.
99  *
100  * This function is intended for writing ASCII output
101  * and returns atom numbers starting at 1.
102  * When dd=NULL returns i+1.
103  */
104 int ddglatnr(const gmx_domdec_t* dd, int i);
105
106 /*! \brief Returns a list of update group partitioning for each molecule type or empty when update groups are not used */
107 gmx::ArrayRef<const gmx::RangePartitioning> getUpdateGroupingPerMoleculetype(const gmx_domdec_t& dd);
108
109 /*! \brief Store the global cg indices of the home cgs in state,
110  *
111  * This means it can be reset, even after a new DD partitioning.
112  */
113 void dd_store_state(struct gmx_domdec_t* dd, t_state* state);
114
115 /*! \brief Returns a pointer to the gmx_domdec_zones_t struct */
116 struct gmx_domdec_zones_t* domdec_zones(struct gmx_domdec_t* dd);
117
118 /*! \brief Returns the range for atoms in zones*/
119 int dd_numAtomsZones(const gmx_domdec_t& dd);
120
121 /*! \brief Returns the number of home atoms */
122 int dd_numHomeAtoms(const gmx_domdec_t& dd);
123
124 /*! \brief Returns the atom range in the local state for atoms that need to be present in mdatoms */
125 int dd_natoms_mdatoms(const gmx_domdec_t* dd);
126
127 /*! \brief Returns the atom range in the local state for atoms involved in virtual sites */
128 int dd_natoms_vsite(const gmx_domdec_t* dd);
129
130 /*! \brief Sets the atom range for atom in the local state for atoms received in constraints communication */
131 void dd_get_constraint_range(const gmx_domdec_t* dd, int* at_start, int* at_end);
132
133 /*! \libinternal \brief Struct for passing around the number of PME domains */
134 struct NumPmeDomains
135 {
136     int x; //!< The number of PME domains along dimension x
137     int y; //!< The number of PME domains along dimension y
138 };
139
140 /*! \brief Returns the number of PME domains, can be called with dd=NULL */
141 NumPmeDomains getNumPmeDomains(const gmx_domdec_t* dd);
142
143 /*! \brief Returns the set of DD ranks that communicate with pme node cr->nodeid */
144 std::vector<int> get_pme_ddranks(const t_commrec* cr, int pmenodeid);
145
146 /*! \brief Returns the maximum shift for coordinate communication in PME, dim x */
147 int dd_pme_maxshift_x(const gmx_domdec_t* dd);
148
149 /*! \brief Returns the maximum shift for coordinate communication in PME, dim y */
150 int dd_pme_maxshift_y(const gmx_domdec_t* dd);
151
152 /*! \brief Return whether constraints, not including settles, cross domain boundaries */
153 bool ddHaveSplitConstraints(const gmx_domdec_t& dd);
154
155 /*! \brief Return whether update groups are used */
156 bool ddUsesUpdateGroups(const gmx_domdec_t& dd);
157
158 /*! \brief Return whether the DD has a single dimension
159  *
160  * The GPU halo exchange code requires a 1D DD, and its setup code can
161  * use the returned value to understand what it should do.
162  */
163 bool is1D(const gmx_domdec_t& dd);
164
165 /*! \brief Initialize data structures for bonded interactions */
166 void dd_init_bondeds(FILE*                      fplog,
167                      gmx_domdec_t*              dd,
168                      const gmx_mtop_t&          mtop,
169                      const gmx_vsite_t*         vsite,
170                      const t_inputrec*          ir,
171                      gmx_bool                   bBCheck,
172                      gmx::ArrayRef<cginfo_mb_t> cginfo_mb);
173
174 /*! \brief Returns whether molecules are always whole, i.e. not broken by PBC */
175 bool dd_moleculesAreAlwaysWhole(const gmx_domdec_t& dd);
176
177 /*! \brief Returns if we need to do pbc for calculating bonded interactions */
178 gmx_bool dd_bonded_molpbc(const gmx_domdec_t* dd, PbcType pbcType);
179
180 /*! \brief Change the DD non-bonded communication cut-off.
181  *
182  * This could fail when trying to increase the cut-off,
183  * then FALSE will be returned and the cut-off is not modified.
184  *
185  * \param[in] cr               Communication recrod
186  * \param[in] box              Box matrix, used for computing the dimensions of the system
187  * \param[in] x                Position vector, used for computing the dimensions of the system
188  * \param[in] cutoffRequested  The requested atom to atom cut-off distance, usually the pair-list cutoff distance
189  */
190 gmx_bool change_dd_cutoff(t_commrec* cr, const matrix box, gmx::ArrayRef<const gmx::RVec> x, real cutoffRequested);
191
192 /*! \brief Set up communication for averaging GPU wait times over domains
193  *
194  * When domains (PP MPI ranks) share a GPU, the individual GPU wait times
195  * are meaningless, as it depends on the order in which tasks on the same
196  * GPU finish. Therefore there wait times need to be averaged over the ranks
197  * sharing the same GPU. This function sets up the communication for that.
198  */
199 void dd_setup_dlb_resource_sharing(const t_commrec* cr, int gpu_id);
200
201 /*! \brief Cycle counter indices used internally in the domain decomposition */
202 enum
203 {
204     ddCyclStep,
205     ddCyclPPduringPME,
206     ddCyclF,
207     ddCyclWaitGPU,
208     ddCyclPME,
209     ddCyclNr
210 };
211
212 /*! \brief Add the wallcycle count to the DD counter */
213 void dd_cycles_add(const gmx_domdec_t* dd, float cycles, int ddCycl);
214
215 /*! \brief Communicate the coordinates to the neighboring cells and do pbc. */
216 void dd_move_x(struct gmx_domdec_t* dd, const matrix box, gmx::ArrayRef<gmx::RVec> x, gmx_wallcycle* wcycle);
217
218 /*! \brief Sum the forces over the neighboring cells.
219  *
220  * When fshift!=NULL the shift forces are updated to obtain
221  * the correct virial from the single sum including f.
222  */
223 void dd_move_f(struct gmx_domdec_t* dd, gmx::ForceWithShiftForces* forceWithShiftForces, gmx_wallcycle* wcycle);
224
225 /*! \brief Communicate a real for each atom to the neighboring cells. */
226 void dd_atom_spread_real(struct gmx_domdec_t* dd, real v[]);
227
228 /*! \brief Sum the contributions to a real for each atom over the neighboring cells. */
229 void dd_atom_sum_real(struct gmx_domdec_t* dd, real v[]);
230
231 /*! \brief Reset all the statistics and counters for total run counting */
232 void reset_dd_statistics_counters(struct gmx_domdec_t* dd);
233
234 /* In domdec_con.c */
235
236 /*! \brief Communicates the virtual site forces, reduces the shift forces when \p fshift != NULL */
237 void dd_move_f_vsites(struct gmx_domdec_t* dd, rvec* f, rvec* fshift);
238
239 /*! \brief Clears the forces for virtual sites */
240 void dd_clear_f_vsites(struct gmx_domdec_t* dd, rvec* f);
241
242 /*! \brief Move x0 and also x1 if x1!=NULL. bX1IsCoord tells if to do PBC on x1 */
243 void dd_move_x_constraints(struct gmx_domdec_t*     dd,
244                            const matrix             box,
245                            gmx::ArrayRef<gmx::RVec> x0,
246                            gmx::ArrayRef<gmx::RVec> x1,
247                            gmx_bool                 bX1IsCoord);
248
249 /*! \brief Communicates the coordinates involved in virtual sites */
250 void dd_move_x_vsites(struct gmx_domdec_t* dd, const matrix box, rvec* x);
251
252 /*! \brief Returns the local atom count array for all constraints
253  *
254  * The local atom count for a constraint, possible values 2/1/0, is needed
255  * to avoid not/double-counting contributions linked to the Lagrange
256  * multiplier, such as the virial and free-energy derivatives.
257  *
258  * \note When \p dd = nullptr, an empty reference is returned.
259  */
260 gmx::ArrayRef<const int> dd_constraints_nlocalatoms(const gmx_domdec_t* dd);
261
262 /* In domdec_top.c */
263
264 /*! \brief Print error output when interactions are missing */
265 [[noreturn]] void dd_print_missing_interactions(const gmx::MDLogger&           mdlog,
266                                                 t_commrec*                     cr,
267                                                 int                            local_count,
268                                                 const gmx_mtop_t*              top_global,
269                                                 const gmx_localtop_t*          top_local,
270                                                 gmx::ArrayRef<const gmx::RVec> x,
271                                                 const matrix                   box);
272
273 /*! \brief Generate and store the reverse topology */
274 void dd_make_reverse_top(FILE*              fplog,
275                          gmx_domdec_t*      dd,
276                          const gmx_mtop_t*  mtop,
277                          const gmx_vsite_t* vsite,
278                          const t_inputrec*  ir,
279                          gmx_bool           bBCheck);
280
281 /*! \brief Generate the local topology and virtual site data */
282 void dd_make_local_top(struct gmx_domdec_t*       dd,
283                        struct gmx_domdec_zones_t* zones,
284                        int                        npbcdim,
285                        matrix                     box,
286                        rvec                       cellsize_min,
287                        const ivec                 npulse,
288                        t_forcerec*                fr,
289                        rvec*                      cgcm_or_x,
290                        const gmx_mtop_t&          top,
291                        gmx_localtop_t*            ltop);
292
293 /*! \brief Sort ltop->ilist when we are doing free energy. */
294 void dd_sort_local_top(gmx_domdec_t* dd, const t_mdatoms* mdatoms, gmx_localtop_t* ltop);
295
296 /*! \brief Construct local state */
297 void dd_init_local_state(struct gmx_domdec_t* dd, const t_state* state_global, t_state* local_state);
298
299 /*! \brief Generate a list of links between atoms that are linked by bonded interactions
300  *
301  * Also stores whether atoms are linked in \p cginfo_mb.
302  */
303 t_blocka* makeBondedLinks(const gmx_mtop_t& mtop, gmx::ArrayRef<cginfo_mb_t> cginfo_mb);
304
305 /*! \brief Calculate the maximum distance involved in 2-body and multi-body bonded interactions */
306 void dd_bonded_cg_distance(const gmx::MDLogger& mdlog,
307                            const gmx_mtop_t*    mtop,
308                            const t_inputrec*    ir,
309                            const rvec*          x,
310                            const matrix         box,
311                            gmx_bool             bBCheck,
312                            real*                r_2b,
313                            real*                r_mb);
314
315 /*! \brief Construct the GPU halo exchange object(s)
316  * \param[in] mdlog          The logger object
317  * \param[in] cr             The commrec object
318  * \param[in] deviceContext  GPU device context
319  * \param[in] streamLocal    The local GPU stream
320  * \param[in] streamNonLocal The non-local GPU stream
321  */
322 void constructGpuHaloExchange(const gmx::MDLogger& mdlog,
323                               const t_commrec&     cr,
324                               const DeviceContext& deviceContext,
325                               void*                streamLocal,
326                               void*                streamNonLocal);
327
328 /*! \brief
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
333  */
334 void reinitGpuHaloExchange(const t_commrec&        cr,
335                            DeviceBuffer<gmx::RVec> d_coordinatesBuffer,
336                            DeviceBuffer<gmx::RVec> d_forcesBuffer);
337
338
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
343  */
344 void communicateGpuHaloCoordinates(const t_commrec&      cr,
345                                    const matrix          box,
346                                    GpuEventSynchronizer* coordinatesReadyOnDeviceEvent);
347
348
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
352  */
353 void communicateGpuHaloForces(const t_commrec& cr, bool accumulateForces);
354
355 #endif