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