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