206c4a065b1c2f14084310d9599af9fb34cb0370
[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, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \defgroup module_domdec Spatial domain decomposition (for parallelization over MPI)
37  * \ingroup group_mdrun
38  *
39  * \brief Manages the decomposition of the simulation volume over MPI
40  * ranks to try to distribute work evenly with minimal communication
41  * overheads.
42  *
43  * \todo Get domdec stuff out of mdtypes/commrec.h
44  *
45  * \author Berk Hess <hess@kth.se>
46  *
47  */
48
49 /*! \libinternal \file
50  *
51  * \brief This file declares functions for mdrun to call to manage the
52  * details of its domain decomposition.
53  *
54  * \author Berk Hess <hess@kth.se>
55  * \inlibraryapi
56  * \ingroup module_domdec
57  */
58
59 #ifndef GMX_DOMDEC_DOMDEC_H
60 #define GMX_DOMDEC_DOMDEC_H
61
62 #include <vector>
63
64 #include "gromacs/math/vectypes.h"
65 #include "gromacs/utility/arrayref.h"
66 #include "gromacs/utility/basedefinitions.h"
67 #include "gromacs/utility/real.h"
68
69 struct cginfo_mb_t;
70 struct DDSystemInfo;
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 class t_state;
86
87 namespace gmx
88 {
89 class ForceWithShiftForces;
90 class MDLogger;
91 class LocalAtomSetManager;
92 class RangePartitioning;
93 struct DomdecOptions;
94 struct MdrunOptions;
95 } // namespace
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 Sets the j-charge-group range for i-charge-group \p icg */
118 void dd_get_ns_ranges(const gmx_domdec_t *dd, int icg,
119                       int *jcg0, int *jcg1, ivec shift0, ivec shift1);
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,
132                              int *at_start, int *at_end);
133
134 /*! \libinternal \brief Struct for passing around the number of PME domains */
135 struct NumPmeDomains
136 {
137     int x; //!< The number of PME domains along dimension x
138     int y; //!< The number of PME domains along dimension y
139 };
140
141 /*! \brief Returns the number of PME domains, can be called with dd=NULL */
142 NumPmeDomains getNumPmeDomains(const gmx_domdec_t *dd);
143
144 /*! \brief Returns the set of DD ranks that communicate with pme node cr->nodeid */
145 std::vector<int> get_pme_ddranks(const t_commrec *cr, int pmenodeid);
146
147 /*! \brief Returns the maximum shift for coordinate communication in PME, dim x */
148 int dd_pme_maxshift_x(const gmx_domdec_t *dd);
149
150 /*! \brief Returns the maximum shift for coordinate communication in PME, dim y */
151 int dd_pme_maxshift_y(const gmx_domdec_t *dd);
152
153 /*! \brief Return whether constraints, not including settles, cross domain boundaries */
154 bool ddHaveSplitConstraints(const gmx_domdec_t &dd);
155
156 /*! \brief Initialized the domain decomposition, chooses the DD grid and PME ranks, return the DD struct */
157 gmx_domdec_t *
158 init_domain_decomposition(const gmx::MDLogger            &mdlog,
159                           t_commrec                      *cr,
160                           const gmx::DomdecOptions       &options,
161                           const gmx::MdrunOptions        &mdrunOptions,
162                           bool                            prefer1DAnd1Pulse,
163                           const gmx_mtop_t               *mtop,
164                           const t_inputrec               *ir,
165                           const matrix                    box,
166                           gmx::ArrayRef<const gmx::RVec>  xGlobal,
167                           gmx::LocalAtomSetManager       *atomSets);
168
169 /*! \brief Return whether the DD has a single dimension with a single pulse
170  *
171  * The GPU halo exchange code requires a 1D single-pulse DD, and its
172  * setup code can use the returned value to understand what it should
173  * do. */
174 bool is1DAnd1PulseDD(const gmx_domdec_t &dd);
175
176 /*! \brief Initialize data structures for bonded interactions */
177 void dd_init_bondeds(FILE              *fplog,
178                      gmx_domdec_t      *dd,
179                      const gmx_mtop_t  *mtop,
180                      const gmx_vsite_t *vsite,
181                      const t_inputrec  *ir,
182                      gmx_bool           bBCheck,
183                      cginfo_mb_t       *cginfo_mb);
184
185 /*! \brief Returns whether molecules are always whole, i.e. not broken by PBC */
186 bool dd_moleculesAreAlwaysWhole(const gmx_domdec_t &dd);
187
188 /*! \brief Returns if we need to do pbc for calculating bonded interactions */
189 gmx_bool dd_bonded_molpbc(const gmx_domdec_t *dd, int ePBC);
190
191 /*! \brief Change the DD non-bonded communication cut-off.
192  *
193  * This could fail when trying to increase the cut-off,
194  * then FALSE will be returned and the cut-off is not modified.
195  *
196  * \param[in] cr               Communication recrod
197  * \param[in] box              Box matrix, used for computing the dimensions of the system
198  * \param[in] x                Position vector, used for computing the dimensions of the system
199  * \param[in] cutoffRequested  The requested atom to atom cut-off distance, usually the pair-list cutoff distance
200  */
201 gmx_bool change_dd_cutoff(t_commrec                     *cr,
202                           const matrix                   box,
203                           gmx::ArrayRef<const gmx::RVec> x,
204                           real                           cutoffRequested);
205
206 /*! \brief Set up communication for averaging GPU wait times over domains
207  *
208  * When domains (PP MPI ranks) share a GPU, the individual GPU wait times
209  * are meaningless, as it depends on the order in which tasks on the same
210  * GPU finish. Therefore there wait times need to be averaged over the ranks
211  * sharing the same GPU. This function sets up the communication for that.
212  */
213 void dd_setup_dlb_resource_sharing(const t_commrec     *cr,
214                                    int                  gpu_id);
215
216 /*! \brief Cycle counter indices used internally in the domain decomposition */
217 enum {
218     ddCyclStep, ddCyclPPduringPME, ddCyclF, ddCyclWaitGPU, ddCyclPME, ddCyclNr
219 };
220
221 /*! \brief Add the wallcycle count to the DD counter */
222 void dd_cycles_add(const gmx_domdec_t *dd, float cycles, int ddCycl);
223
224 /*! \brief Communicate the coordinates to the neighboring cells and do pbc. */
225 void dd_move_x(struct gmx_domdec_t      *dd,
226                const matrix              box,
227                gmx::ArrayRef<gmx::RVec>  x,
228                gmx_wallcycle            *wcycle);
229
230 /*! \brief Sum the forces over the neighboring cells.
231  *
232  * When fshift!=NULL the shift forces are updated to obtain
233  * the correct virial from the single sum including f.
234  */
235 void dd_move_f(struct gmx_domdec_t       *dd,
236                gmx::ForceWithShiftForces *forceWithShiftForces,
237                gmx_wallcycle             *wcycle);
238
239 /*! \brief Communicate a real for each atom to the neighboring cells. */
240 void dd_atom_spread_real(struct gmx_domdec_t *dd, real v[]);
241
242 /*! \brief Sum the contributions to a real for each atom over the neighboring cells. */
243 void dd_atom_sum_real(struct gmx_domdec_t *dd, real v[]);
244
245 /*! \brief Reset all the statistics and counters for total run counting */
246 void reset_dd_statistics_counters(struct gmx_domdec_t *dd);
247
248 /* In domdec_con.c */
249
250 /*! \brief Communicates the virtual site forces, reduces the shift forces when \p fshift != NULL */
251 void dd_move_f_vsites(struct gmx_domdec_t *dd, rvec *f, rvec *fshift);
252
253 /*! \brief Clears the forces for virtual sites */
254 void dd_clear_f_vsites(struct gmx_domdec_t *dd, rvec *f);
255
256 /*! \brief Move x0 and also x1 if x1!=NULL. bX1IsCoord tells if to do PBC on x1 */
257 void dd_move_x_constraints(struct gmx_domdec_t *dd, const matrix box,
258                            rvec *x0, rvec *x1, gmx_bool bX1IsCoord);
259
260 /*! \brief Communicates the coordinates involved in virtual sites */
261 void dd_move_x_vsites(struct gmx_domdec_t *dd, const matrix box, rvec *x);
262
263 /*! \brief Returns the local atom count array for all constraints
264  *
265  * The local atom count for a constraint, possible values 2/1/0, is needed
266  * to avoid not/double-counting contributions linked to the Lagrange
267  * multiplier, such as the virial and free-energy derivatives.
268  *
269  * \note When \p dd = nullptr, an empty reference is returned.
270  */
271 gmx::ArrayRef<const int> dd_constraints_nlocalatoms(const gmx_domdec_t *dd);
272
273 /* In domdec_top.c */
274
275 /*! \brief Print error output when interactions are missing */
276 [[ noreturn ]]
277 void dd_print_missing_interactions(const gmx::MDLogger  &mdlog,
278                                    t_commrec            *cr,
279                                    int                   local_count,
280                                    const gmx_mtop_t     *top_global,
281                                    const gmx_localtop_t *top_local,
282                                    const rvec           *x,
283                                    const matrix          box);
284
285 /*! \brief Generate and store the reverse topology */
286 void dd_make_reverse_top(FILE *fplog,
287                          gmx_domdec_t *dd, const gmx_mtop_t *mtop,
288                          const gmx_vsite_t *vsite,
289                          const t_inputrec *ir, gmx_bool bBCheck);
290
291 /*! \brief Generate the local topology and virtual site data */
292 void dd_make_local_top(struct gmx_domdec_t       *dd,
293                        struct gmx_domdec_zones_t *zones,
294                        int                        npbcdim,
295                        matrix                     box,
296                        rvec                       cellsize_min,
297                        const ivec                 npulse,
298                        t_forcerec                *fr,
299                        rvec                      *cgcm_or_x,
300                        const gmx_mtop_t          &top,
301                        gmx_localtop_t            *ltop);
302
303 /*! \brief Sort ltop->ilist when we are doing free energy. */
304 void dd_sort_local_top(gmx_domdec_t *dd, const t_mdatoms *mdatoms,
305                        gmx_localtop_t *ltop);
306
307 /*! \brief Initialize local topology
308  *
309  * \param[in] top_global Reference to global topology.
310  * \param[in,out] top Pointer to new local topology
311  */
312 void dd_init_local_top(const gmx_mtop_t &top_global,
313                        gmx_localtop_t   *top);
314
315 /*! \brief Construct local state */
316 void dd_init_local_state(struct gmx_domdec_t *dd,
317                          const t_state *state_global, t_state *local_state);
318
319 /*! \brief Generate a list of links between atoms that are linked by bonded interactions
320  *
321  * Also stores whether atoms are linked in \p cginfo_mb.
322  */
323 t_blocka *makeBondedLinks(const gmx_mtop_t *mtop,
324                           cginfo_mb_t      *cginfo_mb);
325
326 /*! \brief Calculate the maximum distance involved in 2-body and multi-body bonded interactions */
327 void dd_bonded_cg_distance(const gmx::MDLogger &mdlog,
328                            const gmx_mtop_t *mtop,
329                            const t_inputrec *ir,
330                            const rvec *x, const matrix box,
331                            gmx_bool bBCheck,
332                            real *r_2b, real *r_mb);
333
334 #endif