439e86b5e7fb0fb608c2061894eb1866dcf50033
[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/math/vectypes.h"
66 #include "gromacs/utility/arrayref.h"
67 #include "gromacs/utility/basedefinitions.h"
68 #include "gromacs/utility/real.h"
69
70 struct cginfo_mb_t;
71 struct gmx_domdec_t;
72 struct gmx_ddbox_t;
73 struct gmx_domdec_zones_t;
74 struct gmx_localtop_t;
75 struct gmx_mtop_t;
76 struct gmx_vsite_t;
77 struct t_block;
78 struct t_blocka;
79 struct t_commrec;
80 struct t_forcerec;
81 struct t_inputrec;
82 struct t_mdatoms;
83 struct t_nrnb;
84 struct gmx_wallcycle;
85 enum class PbcType : int;
86 class t_state;
87
88 namespace gmx
89 {
90 class ForceWithShiftForces;
91 class MDLogger;
92 class RangePartitioning;
93 } // namespace gmx
94
95 /*! \brief Returns the global topology atom number belonging to local atom index i.
96  *
97  * This function is intended for writing ASCII output
98  * and returns atom numbers starting at 1.
99  * When dd=NULL returns i+1.
100  */
101 int ddglatnr(const gmx_domdec_t* dd, int i);
102
103 /*! \brief Returns a list of update group partitioning for each molecule type or empty when update groups are not used */
104 gmx::ArrayRef<const gmx::RangePartitioning> getUpdateGroupingPerMoleculetype(const gmx_domdec_t& dd);
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(struct 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 constraints, not including settles, cross domain boundaries */
150 bool ddHaveSplitConstraints(const gmx_domdec_t& dd);
151
152 /*! \brief Return whether update groups are used */
153 bool ddUsesUpdateGroups(const gmx_domdec_t& dd);
154
155 /*! \brief Return whether the DD has a single dimension with a single pulse
156  *
157  * The GPU halo exchange code requires a 1D single-pulse DD, and its
158  * setup code can use the returned value to understand what it should
159  * do. */
160 bool is1DAnd1PulseDD(const gmx_domdec_t& dd);
161
162 /*! \brief Initialize data structures for bonded interactions */
163 void dd_init_bondeds(FILE*                      fplog,
164                      gmx_domdec_t*              dd,
165                      const gmx_mtop_t&          mtop,
166                      const gmx_vsite_t*         vsite,
167                      const t_inputrec*          ir,
168                      gmx_bool                   bBCheck,
169                      gmx::ArrayRef<cginfo_mb_t> cginfo_mb);
170
171 /*! \brief Returns whether molecules are always whole, i.e. not broken by PBC */
172 bool dd_moleculesAreAlwaysWhole(const gmx_domdec_t& dd);
173
174 /*! \brief Returns if we need to do pbc for calculating bonded interactions */
175 gmx_bool dd_bonded_molpbc(const gmx_domdec_t* dd, PbcType pbcType);
176
177 /*! \brief Change the DD non-bonded communication cut-off.
178  *
179  * This could fail when trying to increase the cut-off,
180  * then FALSE will be returned and the cut-off is not modified.
181  *
182  * \param[in] cr               Communication recrod
183  * \param[in] box              Box matrix, used for computing the dimensions of the system
184  * \param[in] x                Position vector, used for computing the dimensions of the system
185  * \param[in] cutoffRequested  The requested atom to atom cut-off distance, usually the pair-list cutoff distance
186  */
187 gmx_bool change_dd_cutoff(t_commrec* cr, const matrix box, gmx::ArrayRef<const gmx::RVec> x, real cutoffRequested);
188
189 /*! \brief Set up communication for averaging GPU wait times over domains
190  *
191  * When domains (PP MPI ranks) share a GPU, the individual GPU wait times
192  * are meaningless, as it depends on the order in which tasks on the same
193  * GPU finish. Therefore there wait times need to be averaged over the ranks
194  * sharing the same GPU. This function sets up the communication for that.
195  */
196 void dd_setup_dlb_resource_sharing(const t_commrec* cr, int gpu_id);
197
198 /*! \brief Cycle counter indices used internally in the domain decomposition */
199 enum
200 {
201     ddCyclStep,
202     ddCyclPPduringPME,
203     ddCyclF,
204     ddCyclWaitGPU,
205     ddCyclPME,
206     ddCyclNr
207 };
208
209 /*! \brief Add the wallcycle count to the DD counter */
210 void dd_cycles_add(const gmx_domdec_t* dd, float cycles, int ddCycl);
211
212 /*! \brief Communicate the coordinates to the neighboring cells and do pbc. */
213 void dd_move_x(struct gmx_domdec_t* dd, const matrix box, gmx::ArrayRef<gmx::RVec> x, gmx_wallcycle* wcycle);
214
215 /*! \brief Sum the forces over the neighboring cells.
216  *
217  * When fshift!=NULL the shift forces are updated to obtain
218  * the correct virial from the single sum including f.
219  */
220 void dd_move_f(struct gmx_domdec_t* dd, gmx::ForceWithShiftForces* forceWithShiftForces, gmx_wallcycle* wcycle);
221
222 /*! \brief Communicate a real for each atom to the neighboring cells. */
223 void dd_atom_spread_real(struct gmx_domdec_t* dd, real v[]);
224
225 /*! \brief Sum the contributions to a real for each atom over the neighboring cells. */
226 void dd_atom_sum_real(struct gmx_domdec_t* dd, real v[]);
227
228 /*! \brief Reset all the statistics and counters for total run counting */
229 void reset_dd_statistics_counters(struct gmx_domdec_t* dd);
230
231 /* In domdec_con.c */
232
233 /*! \brief Communicates the virtual site forces, reduces the shift forces when \p fshift != NULL */
234 void dd_move_f_vsites(struct gmx_domdec_t* dd, rvec* f, rvec* fshift);
235
236 /*! \brief Clears the forces for virtual sites */
237 void dd_clear_f_vsites(struct gmx_domdec_t* dd, rvec* f);
238
239 /*! \brief Move x0 and also x1 if x1!=NULL. bX1IsCoord tells if to do PBC on x1 */
240 void dd_move_x_constraints(struct gmx_domdec_t*     dd,
241                            const matrix             box,
242                            gmx::ArrayRef<gmx::RVec> x0,
243                            gmx::ArrayRef<gmx::RVec> x1,
244                            gmx_bool                 bX1IsCoord);
245
246 /*! \brief Communicates the coordinates involved in virtual sites */
247 void dd_move_x_vsites(struct gmx_domdec_t* dd, const matrix box, rvec* x);
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                                                 const 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_vsite_t* vsite,
275                          const t_inputrec*  ir,
276                          gmx_bool           bBCheck);
277
278 /*! \brief Generate the local topology and virtual site data */
279 void dd_make_local_top(struct gmx_domdec_t*       dd,
280                        struct gmx_domdec_zones_t* zones,
281                        int                        npbcdim,
282                        matrix                     box,
283                        rvec                       cellsize_min,
284                        const ivec                 npulse,
285                        t_forcerec*                fr,
286                        rvec*                      cgcm_or_x,
287                        const gmx_mtop_t&          top,
288                        gmx_localtop_t*            ltop);
289
290 /*! \brief Sort ltop->ilist when we are doing free energy. */
291 void dd_sort_local_top(gmx_domdec_t* dd, const t_mdatoms* mdatoms, gmx_localtop_t* ltop);
292
293 /*! \brief Initialize local topology
294  *
295  * \param[in] top_global Reference to global topology.
296  * \param[in,out] top Pointer to new local topology
297  */
298 void dd_init_local_top(const gmx_mtop_t& top_global, gmx_localtop_t* top);
299
300 /*! \brief Construct local state */
301 void dd_init_local_state(struct gmx_domdec_t* dd, const t_state* state_global, t_state* local_state);
302
303 /*! \brief Generate a list of links between atoms that are linked by bonded interactions
304  *
305  * Also stores whether atoms are linked in \p cginfo_mb.
306  */
307 t_blocka* makeBondedLinks(const gmx_mtop_t& mtop, gmx::ArrayRef<cginfo_mb_t> cginfo_mb);
308
309 /*! \brief Calculate the maximum distance involved in 2-body and multi-body bonded interactions */
310 void dd_bonded_cg_distance(const gmx::MDLogger& mdlog,
311                            const gmx_mtop_t*    mtop,
312                            const t_inputrec*    ir,
313                            const rvec*          x,
314                            const matrix         box,
315                            gmx_bool             bBCheck,
316                            real*                r_2b,
317                            real*                r_mb);
318
319 #endif