6f8fa8c8ef42ff0a174175b26966a353f2401783
[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 Return a block struct for the charge groups of the whole system */
106 t_block *dd_charge_groups_global(struct gmx_domdec_t *dd);
107
108 /*! \brief Returns a list of update group partitioning for each molecule type or empty when update groups are not used */
109 gmx::ArrayRef<const gmx::RangePartitioning> getUpdateGroupingPerMoleculetype(const gmx_domdec_t &dd);
110
111 /*! \brief Store the global cg indices of the home cgs in state,
112  *
113  * This means it can be reset, even after a new DD partitioning.
114  */
115 void dd_store_state(struct gmx_domdec_t *dd, t_state *state);
116
117 /*! \brief Returns a pointer to the gmx_domdec_zones_t struct */
118 struct gmx_domdec_zones_t *domdec_zones(struct gmx_domdec_t *dd);
119
120 /*! \brief Sets the j-charge-group range for i-charge-group \p icg */
121 void dd_get_ns_ranges(const gmx_domdec_t *dd, int icg,
122                       int *jcg0, int *jcg1, ivec shift0, ivec shift1);
123
124 /*! \brief Returns the number of home atoms */
125 int dd_numHomeAtoms(const gmx_domdec_t &dd);
126
127 /*! \brief Returns the atom range in the local state for atoms that need to be present in mdatoms */
128 int dd_natoms_mdatoms(const gmx_domdec_t *dd);
129
130 /*! \brief Returns the atom range in the local state for atoms involved in virtual sites */
131 int dd_natoms_vsite(const gmx_domdec_t *dd);
132
133 /*! \brief Sets the atom range for atom in the local state for atoms received in constraints communication */
134 void dd_get_constraint_range(const gmx_domdec_t *dd,
135                              int *at_start, int *at_end);
136
137 /*! \libinternal \brief Struct for passing around the number of PME domains */
138 struct NumPmeDomains
139 {
140     int x; //!< The number of PME domains along dimension x
141     int y; //!< The number of PME domains along dimension y
142 };
143
144 /*! \brief Returns the number of PME domains, can be called with dd=NULL */
145 NumPmeDomains getNumPmeDomains(const gmx_domdec_t *dd);
146
147 /*! \brief Returns the set of DD ranks that communicate with pme node cr->nodeid */
148 std::vector<int> get_pme_ddranks(const t_commrec *cr, int pmenodeid);
149
150 /*! \brief Returns the maximum shift for coordinate communication in PME, dim x */
151 int dd_pme_maxshift_x(const gmx_domdec_t *dd);
152
153 /*! \brief Returns the maximum shift for coordinate communication in PME, dim y */
154 int dd_pme_maxshift_y(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                           const gmx_mtop_t               *mtop,
163                           const t_inputrec               *ir,
164                           const matrix                    box,
165                           gmx::ArrayRef<const gmx::RVec>  xGlobal,
166                           gmx::LocalAtomSetManager       *atomSets);
167
168 /*! \brief Initialize data structures for bonded interactions */
169 void dd_init_bondeds(FILE              *fplog,
170                      gmx_domdec_t      *dd,
171                      const gmx_mtop_t  *mtop,
172                      const gmx_vsite_t *vsite,
173                      const t_inputrec  *ir,
174                      gmx_bool           bBCheck,
175                      cginfo_mb_t       *cginfo_mb);
176
177 /*! \brief Returns if we need to do pbc for calculating bonded interactions */
178 gmx_bool dd_bonded_molpbc(const gmx_domdec_t *dd, int ePBC);
179
180 /*! \brief Change the DD non-bonded communication cut-off.
181  *
182  * This could fail when trying to increase the cut-off,
183  * then FALSE will be returned and the cut-off is not modified.
184  *
185  * \param[in] cr               Communication recrod
186  * \param[in] state            State, 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,
190                           const t_state &state,
191                           real           cutoffRequested);
192
193 /*! \brief Set up communication for averaging GPU wait times over domains
194  *
195  * When domains (PP MPI ranks) share a GPU, the individual GPU wait times
196  * are meaningless, as it depends on the order in which tasks on the same
197  * GPU finish. Therefore there wait times need to be averaged over the ranks
198  * sharing the same GPU. This function sets up the communication for that.
199  */
200 void dd_setup_dlb_resource_sharing(t_commrec           *cr,
201                                    int                  gpu_id);
202
203 /*! \brief Cycle counter indices used internally in the domain decomposition */
204 enum {
205     ddCyclStep, ddCyclPPduringPME, ddCyclF, ddCyclWaitGPU, ddCyclPME, ddCyclNr
206 };
207
208 /*! \brief Add the wallcycle count to the DD counter */
209 void dd_cycles_add(const gmx_domdec_t *dd, float cycles, int ddCycl);
210
211 /*! \brief Communicate the coordinates to the neighboring cells and do pbc. */
212 void dd_move_x(struct gmx_domdec_t      *dd,
213                const matrix              box,
214                gmx::ArrayRef<gmx::RVec>  x,
215                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,
223                gmx::ForceWithShiftForces *forceWithShiftForces,
224                gmx_wallcycle             *wcycle);
225
226 /*! \brief Communicate a real for each atom to the neighboring cells. */
227 void dd_atom_spread_real(struct gmx_domdec_t *dd, real v[]);
228
229 /*! \brief Sum the contributions to a real for each atom over the neighboring cells. */
230 void dd_atom_sum_real(struct gmx_domdec_t *dd, real v[]);
231
232 /*! \brief Reset all the statistics and counters for total run counting */
233 void reset_dd_statistics_counters(struct gmx_domdec_t *dd);
234
235 /* In domdec_con.c */
236
237 /*! \brief Communicates the virtual site forces, reduces the shift forces when \p fshift != NULL */
238 void dd_move_f_vsites(struct gmx_domdec_t *dd, rvec *f, rvec *fshift);
239
240 /*! \brief Clears the forces for virtual sites */
241 void dd_clear_f_vsites(struct gmx_domdec_t *dd, rvec *f);
242
243 /*! \brief Move x0 and also x1 if x1!=NULL. bX1IsCoord tells if to do PBC on x1 */
244 void dd_move_x_constraints(struct gmx_domdec_t *dd, matrix box,
245                            rvec *x0, rvec *x1, gmx_bool bX1IsCoord);
246
247 /*! \brief Communicates the coordinates involved in virtual sites */
248 void dd_move_x_vsites(struct gmx_domdec_t *dd, const matrix box, rvec *x);
249
250 /*! \brief Returns the local atom count array for all constraints
251  *
252  * The local atom count for a constraint, possible values 2/1/0, is needed
253  * to avoid not/double-counting contributions linked to the Lagrange
254  * multiplier, such as the virial and free-energy derivatives.
255  *
256  * \note When \p dd = nullptr, an empty reference is returned.
257  */
258 gmx::ArrayRef<const int> dd_constraints_nlocalatoms(const gmx_domdec_t *dd);
259
260 /* In domdec_top.c */
261
262 /*! \brief Print error output when interactions are missing */
263 [[ noreturn ]]
264 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                                    const 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, const gmx_mtop_t *mtop,
275                          const gmx_vsite_t *vsite,
276                          const t_inputrec *ir, 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,
292                        gmx_localtop_t *ltop);
293
294 /*! \brief Initialize local topology
295  *
296  * \param[in] top_global Reference to global topology.
297  * \param[in,out] top Pointer to new local topology
298  */
299 void dd_init_local_top(const gmx_mtop_t &top_global,
300                        gmx_localtop_t   *top);
301
302 /*! \brief Construct local state */
303 void dd_init_local_state(struct gmx_domdec_t *dd,
304                          const t_state *state_global, t_state *local_state);
305
306 /*! \brief Generate a list of links between charge groups that are linked by bonded interactions */
307 t_blocka *make_charge_group_links(const gmx_mtop_t *mtop, gmx_domdec_t *dd,
308                                   cginfo_mb_t *cginfo_mb);
309
310 /*! \brief Calculate the maximum distance involved in 2-body and multi-body bonded interactions */
311 void dd_bonded_cg_distance(const gmx::MDLogger &mdlog,
312                            const gmx_mtop_t *mtop,
313                            const t_inputrec *ir,
314                            const rvec *x, const matrix box,
315                            gmx_bool bBCheck,
316                            real *r_2b, real *r_mb);
317
318
319 /* In domdec_setup.c */
320
321 /*! \brief Returns the volume fraction of the system that is communicated */
322 real comm_box_frac(const ivec dd_nc, real cutoff, const gmx_ddbox_t *ddbox);
323
324 struct DDSetup
325 {
326     int  numPmeRanks   = 0;
327     ivec numDomains    = { 0, 0, 0 };
328     real cellsizeLimit = 0;
329 };
330
331 /*! \brief Determines the optimal DD cell setup dd->nc and possibly npmenodes
332  * for the system.
333  */
334 DDSetup
335 dd_choose_grid(const gmx::MDLogger &mdlog,
336                const t_commrec *cr,
337                const t_inputrec *ir,
338                const gmx_mtop_t *mtop,
339                const matrix box, const gmx_ddbox_t *ddbox,
340                int numPmeRanksRequested,
341                gmx_bool bDynLoadBal, real dlb_scale,
342                const DDSystemInfo &systemInfo);
343
344 #endif