Use ObservablesReducer for check of DD bonded interaction count.
[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 gmx_domdec_t;
70 struct gmx_ddbox_t;
71 struct gmx_domdec_zones_t;
72 struct gmx_localtop_t;
73 struct gmx_mtop_t;
74 struct t_commrec;
75 struct t_forcerec;
76 struct t_inputrec;
77 struct t_mdatoms;
78 struct t_nrnb;
79 struct gmx_wallcycle;
80 enum class PbcType : int;
81 class t_state;
82 class DeviceContext;
83 class GpuEventSynchronizer;
84
85 namespace gmx
86 {
87 struct AtomInfoWithinMoleculeBlock;
88 class DeviceStreamManager;
89 class ForceWithShiftForces;
90 class MDLogger;
91 class RangePartitioning;
92 class VirtualSitesHandler;
93 template<typename>
94 class ArrayRef;
95 } // namespace gmx
96
97 /*! \brief Returns the global topology atom number belonging to local atom index i.
98  *
99  * This function is intended for writing ASCII output
100  * and returns atom numbers starting at 1.
101  * When dd=NULL returns i+1.
102  */
103 int ddglatnr(const gmx_domdec_t* dd, int i);
104
105 /*! \brief Store the global cg indices of the home cgs in state,
106  *
107  * This means it can be reset, even after a new DD partitioning.
108  */
109 void dd_store_state(const gmx_domdec_t& dd, t_state* state);
110
111 /*! \brief Returns a pointer to the gmx_domdec_zones_t struct */
112 struct gmx_domdec_zones_t* domdec_zones(struct gmx_domdec_t* dd);
113
114 /*! \brief Returns the range for atoms in zones*/
115 int dd_numAtomsZones(const gmx_domdec_t& dd);
116
117 /*! \brief Returns the number of home atoms */
118 int dd_numHomeAtoms(const gmx_domdec_t& dd);
119
120 /*! \brief Returns the atom range in the local state for atoms that need to be present in mdatoms */
121 int dd_natoms_mdatoms(const gmx_domdec_t& dd);
122
123 /*! \brief Returns the atom range in the local state for atoms involved in virtual sites */
124 int dd_natoms_vsite(const gmx_domdec_t& dd);
125
126 /*! \brief Sets the atom range for atom in the local state for atoms received in constraints communication */
127 void dd_get_constraint_range(const gmx_domdec_t& dd, int* at_start, int* at_end);
128
129 /*! \libinternal \brief Struct for passing around the number of PME domains */
130 struct NumPmeDomains
131 {
132     int x; //!< The number of PME domains along dimension x
133     int y; //!< The number of PME domains along dimension y
134 };
135
136 /*! \brief Returns the number of PME domains, can be called with dd=NULL */
137 NumPmeDomains getNumPmeDomains(const gmx_domdec_t* dd);
138
139 /*! \brief Returns the set of DD ranks that communicate with pme node cr->nodeid */
140 std::vector<int> get_pme_ddranks(const t_commrec* cr, int pmenodeid);
141
142 /*! \brief Returns the maximum shift for coordinate communication in PME, dim x */
143 int dd_pme_maxshift_x(const gmx_domdec_t& dd);
144
145 /*! \brief Returns the maximum shift for coordinate communication in PME, dim y */
146 int dd_pme_maxshift_y(const gmx_domdec_t& dd);
147
148 /*! \brief Return whether update groups are used */
149 bool ddUsesUpdateGroups(const gmx_domdec_t& dd);
150
151 /*! \brief Returns whether molecules are always whole, i.e. not broken by PBC */
152 bool dd_moleculesAreAlwaysWhole(const gmx_domdec_t& dd);
153
154 /*! \brief Returns if we need to do pbc for calculating bonded interactions */
155 bool dd_bonded_molpbc(const gmx_domdec_t& dd, PbcType pbcType);
156
157 /*! \brief Change the DD non-bonded communication cut-off.
158  *
159  * This could fail when trying to increase the cut-off,
160  * then FALSE will be returned and the cut-off is not modified.
161  *
162  * \param[in] cr               Communication recrod
163  * \param[in] box              Box matrix, used for computing the dimensions of the system
164  * \param[in] x                Position vector, used for computing the dimensions of the system
165  * \param[in] cutoffRequested  The requested atom to atom cut-off distance, usually the pair-list cutoff distance
166  */
167 bool change_dd_cutoff(t_commrec* cr, const matrix box, gmx::ArrayRef<const gmx::RVec> x, real cutoffRequested);
168
169 /*! \brief Set up communication for averaging GPU wait times over domains
170  *
171  * When domains (PP MPI ranks) share a GPU, the individual GPU wait times
172  * are meaningless, as it depends on the order in which tasks on the same
173  * GPU finish. Therefore there wait times need to be averaged over the ranks
174  * sharing the same GPU. This function sets up the communication for that.
175  */
176 void dd_setup_dlb_resource_sharing(const t_commrec* cr, int gpu_id);
177
178 /*! \brief Cycle counter indices used internally in the domain decomposition */
179 enum
180 {
181     ddCyclStep,
182     ddCyclPPduringPME,
183     ddCyclF,
184     ddCyclWaitGPU,
185     ddCyclPME,
186     ddCyclNr
187 };
188
189 /*! \brief Add the wallcycle count to the DD counter */
190 void dd_cycles_add(const gmx_domdec_t* dd, float cycles, int ddCycl);
191
192 /*! \brief Communicate the coordinates to the neighboring cells and do pbc. */
193 void dd_move_x(struct gmx_domdec_t* dd, const matrix box, gmx::ArrayRef<gmx::RVec> x, gmx_wallcycle* wcycle);
194
195 /*! \brief Sum the forces over the neighboring cells.
196  *
197  * When fshift!=NULL the shift forces are updated to obtain
198  * the correct virial from the single sum including f.
199  */
200 void dd_move_f(struct gmx_domdec_t* dd, gmx::ForceWithShiftForces* forceWithShiftForces, gmx_wallcycle* wcycle);
201
202 /*! \brief Communicate a real for each atom to the neighboring cells. */
203 void dd_atom_spread_real(struct gmx_domdec_t* dd, real v[]);
204
205 /*! \brief Sum the contributions to a real for each atom over the neighboring cells. */
206 void dd_atom_sum_real(struct gmx_domdec_t* dd, real v[]);
207
208 /*! \brief Reset all the statistics and counters for total run counting */
209 void reset_dd_statistics_counters(struct gmx_domdec_t* dd);
210
211 /* In domdec_con.c */
212
213 /*! \brief Communicates the virtual site forces, reduces the shift forces when \p fshift != NULL */
214 void dd_move_f_vsites(const gmx_domdec_t& dd, gmx::ArrayRef<gmx::RVec> f, gmx::ArrayRef<gmx::RVec> fshift);
215
216 /*! \brief Clears the forces for virtual sites */
217 void dd_clear_f_vsites(const gmx_domdec_t& dd, gmx::ArrayRef<gmx::RVec> f);
218
219 /*! \brief Move x0 and also x1 if x1!=NULL. bX1IsCoord tells if to do PBC on x1 */
220 void dd_move_x_constraints(struct gmx_domdec_t*     dd,
221                            const matrix             box,
222                            gmx::ArrayRef<gmx::RVec> x0,
223                            gmx::ArrayRef<gmx::RVec> x1,
224                            bool                     bX1IsCoord);
225
226 /*! \brief Communicates the coordinates involved in virtual sites */
227 void dd_move_x_vsites(const gmx_domdec_t& dd, const matrix box, rvec* x);
228 /*! \brief Communicates the positions and velocities involved in virtual sites */
229 void dd_move_x_and_v_vsites(const gmx_domdec_t& dd, const matrix box, rvec* x, rvec* v);
230
231 /*! \brief Returns the local atom count array for all constraints
232  *
233  * The local atom count for a constraint, possible values 2/1/0, is needed
234  * to avoid not/double-counting contributions linked to the Lagrange
235  * multiplier, such as the virial and free-energy derivatives.
236  *
237  * \note When \p dd = nullptr, an empty reference is returned.
238  */
239 gmx::ArrayRef<const int> dd_constraints_nlocalatoms(const gmx_domdec_t* dd);
240
241 /*! \brief Construct local state */
242 void dd_init_local_state(const gmx_domdec_t& dd, const t_state* state_global, t_state* local_state);
243
244 /*! \brief Construct the GPU halo exchange object(s).
245  *
246  * \param[in] mdlog               The logger object.
247  * \param[in] cr                  The commrec object.
248  * \param[in] deviceStreamManager Manager of the GPU context and streams.
249  * \param[in] wcycle              The wallclock counter.
250  */
251 void constructGpuHaloExchange(const gmx::MDLogger&            mdlog,
252                               const t_commrec&                cr,
253                               const gmx::DeviceStreamManager& deviceStreamManager,
254                               gmx_wallcycle*                  wcycle);
255
256 /*! \brief
257  * (Re-) Initialization for GPU halo exchange
258  * \param [in] cr                   The commrec object
259  * \param [in] d_coordinatesBuffer  pointer to coordinates buffer in GPU memory
260  * \param [in] d_forcesBuffer       pointer to forces buffer in GPU memory
261  */
262 void reinitGpuHaloExchange(const t_commrec&        cr,
263                            DeviceBuffer<gmx::RVec> d_coordinatesBuffer,
264                            DeviceBuffer<gmx::RVec> d_forcesBuffer);
265
266
267 /*! \brief GPU halo exchange of coordinates buffer.
268  * \param [in] cr                             The commrec object
269  * \param [in] box                            Coordinate box (from which shifts will be constructed)
270  * \param [in] coordinatesReadyOnDeviceEvent  event recorded when coordinates have been copied to device
271  */
272 void communicateGpuHaloCoordinates(const t_commrec&      cr,
273                                    const matrix          box,
274                                    GpuEventSynchronizer* coordinatesReadyOnDeviceEvent);
275
276
277 /*! \brief GPU halo exchange of force buffer.
278  * \param [in] cr                The commrec object
279  * \param [in] accumulateForces  True if forces should accumulate, otherwise they are set
280  */
281 void communicateGpuHaloForces(const t_commrec& cr, bool accumulateForces);
282
283 /*! \brief Wraps the \c positions so that atoms from the same
284  * update group share the same periodic image wrt \c box.
285  *
286  * When DD and update groups are in use, the simulation master rank
287  * should call this to ensure that e.g. when restarting a simulation
288  * that did not use update groups that the coordinates satisfy the new
289  * requirements.
290  *
291  * This function can probably be removed when even single-rank
292  * simulations use domain decomposition, because then the choice of
293  * whether update groups are used is probably going to be the same
294  * regardless of the rank count.
295  *
296  * \param[in]    dd         The DD manager
297  * \param[in]    mtop       The system topology
298  * \param[in]    box        The global system box
299  * \param[in]    positions  The global system positions
300  */
301 void putUpdateGroupAtomsInSamePeriodicImage(const gmx_domdec_t&      dd,
302                                             const gmx_mtop_t&        mtop,
303                                             const matrix             box,
304                                             gmx::ArrayRef<gmx::RVec> positions);
305
306 #endif