d8aae0b18cd68b66282d2d5413da50216c4deb57
[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 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 gmx_vsite_t;
76 struct t_block;
77 struct t_blocka;
78 struct t_commrec;
79 struct t_forcerec;
80 struct t_inputrec;
81 struct t_mdatoms;
82 struct t_nrnb;
83 struct gmx_wallcycle;
84 class t_state;
85
86 namespace gmx
87 {
88 class MDLogger;
89 class LocalAtomSetManager;
90 class RangePartitioning;
91 struct DomdecOptions;
92 struct MdrunOptions;
93 } // namespace
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 Return a block struct for the charge groups of the whole system */
104 t_block *dd_charge_groups_global(struct gmx_domdec_t *dd);
105
106 /*! \brief Returns a list of update group partitioning for each molecule type or empty when update groups are not used */
107 gmx::ArrayRef<const gmx::RangePartitioning> getUpdateGroupingPerMoleculetype(const gmx_domdec_t &dd);
108
109 /*! \brief Store the global cg indices of the home cgs in state,
110  *
111  * This means it can be reset, even after a new DD partitioning.
112  */
113 void dd_store_state(struct gmx_domdec_t *dd, t_state *state);
114
115 /*! \brief Returns a pointer to the gmx_domdec_zones_t struct */
116 struct gmx_domdec_zones_t *domdec_zones(struct gmx_domdec_t *dd);
117
118 /*! \brief Sets the j-charge-group range for i-charge-group \p icg */
119 void dd_get_ns_ranges(const gmx_domdec_t *dd, int icg,
120                       int *jcg0, int *jcg1, ivec shift0, ivec shift1);
121
122 /*! \brief Returns the number of home atoms */
123 int dd_numHomeAtoms(const gmx_domdec_t &dd);
124
125 /*! \brief Returns the atom range in the local state for atoms that need to be present in mdatoms */
126 int dd_natoms_mdatoms(const gmx_domdec_t *dd);
127
128 /*! \brief Returns the atom range in the local state for atoms involved in virtual sites */
129 int dd_natoms_vsite(const gmx_domdec_t *dd);
130
131 /*! \brief Sets the atom range for atom in the local state for atoms received in constraints communication */
132 void dd_get_constraint_range(const gmx_domdec_t *dd,
133                              int *at_start, int *at_end);
134
135 /*! \libinternal \brief Struct for passing around the number of PME domains */
136 struct NumPmeDomains
137 {
138     int x; //!< The number of PME domains along dimension x
139     int y; //!< The number of PME domains along dimension y
140 };
141
142 /*! \brief Returns the number of PME domains, can be called with dd=NULL */
143 NumPmeDomains getNumPmeDomains(const gmx_domdec_t *dd);
144
145 /*! \brief Returns the set of DD ranks that communicate with pme node cr->nodeid */
146 std::vector<int> get_pme_ddranks(const t_commrec *cr, int pmenodeid);
147
148 /*! \brief Returns the maximum shift for coordinate communication in PME, dim x */
149 int dd_pme_maxshift_x(const gmx_domdec_t *dd);
150
151 /*! \brief Returns the maximum shift for coordinate communication in PME, dim y */
152 int dd_pme_maxshift_y(const gmx_domdec_t *dd);
153
154 /*! \brief Initialized the domain decomposition, chooses the DD grid and PME ranks, return the DD struct */
155 gmx_domdec_t *
156 init_domain_decomposition(const gmx::MDLogger            &mdlog,
157                           t_commrec                      *cr,
158                           const gmx::DomdecOptions       &options,
159                           const gmx::MdrunOptions        &mdrunOptions,
160                           const gmx_mtop_t               *mtop,
161                           const t_inputrec               *ir,
162                           const matrix                    box,
163                           gmx::ArrayRef<const gmx::RVec>  xGlobal,
164                           gmx::LocalAtomSetManager       *atomSets);
165
166 /*! \brief Initialize data structures for bonded interactions */
167 void dd_init_bondeds(FILE              *fplog,
168                      gmx_domdec_t      *dd,
169                      const gmx_mtop_t  *mtop,
170                      const gmx_vsite_t *vsite,
171                      const t_inputrec  *ir,
172                      gmx_bool           bBCheck,
173                      cginfo_mb_t       *cginfo_mb);
174
175 /*! \brief Returns if we need to do pbc for calculating bonded interactions */
176 gmx_bool dd_bonded_molpbc(const gmx_domdec_t *dd, int ePBC);
177
178 /*! \brief Change the DD non-bonded communication cut-off.
179  *
180  * This could fail when trying to increase the cut-off,
181  * then FALSE will be returned and the cut-off is not modified.
182  *
183  * \param[in] cr               Communication recrod
184  * \param[in] state            State, 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,
188                           const t_state &state,
189                           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(t_commrec           *cr,
199                                    int                  gpu_id);
200
201 /*! \brief Cycle counter indices used internally in the domain decomposition */
202 enum {
203     ddCyclStep, ddCyclPPduringPME, ddCyclF, ddCyclWaitGPU, ddCyclPME, ddCyclNr
204 };
205
206 /*! \brief Add the wallcycle count to the DD counter */
207 void dd_cycles_add(const gmx_domdec_t *dd, float cycles, int ddCycl);
208
209 /*! \brief Communicate the coordinates to the neighboring cells and do pbc. */
210 void dd_move_x(struct gmx_domdec_t      *dd,
211                matrix                    box,
212                gmx::ArrayRef<gmx::RVec>  x,
213                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,
221                gmx::ArrayRef<gmx::RVec>  f,
222                rvec                     *fshift,
223                gmx_wallcycle            *wcycle);
224
225 /*! \brief Communicate a real for each atom to the neighboring cells. */
226 void dd_atom_spread_real(struct gmx_domdec_t *dd, real v[]);
227
228 /*! \brief Sum the contributions to a real for each atom over the neighboring cells. */
229 void dd_atom_sum_real(struct gmx_domdec_t *dd, real v[]);
230
231 /*! \brief Reset all the statistics and counters for total run counting */
232 void reset_dd_statistics_counters(struct gmx_domdec_t *dd);
233
234 /* In domdec_con.c */
235
236 /*! \brief Communicates the virtual site forces, reduces the shift forces when \p fshift != NULL */
237 void dd_move_f_vsites(struct gmx_domdec_t *dd, rvec *f, rvec *fshift);
238
239 /*! \brief Clears the forces for virtual sites */
240 void dd_clear_f_vsites(struct gmx_domdec_t *dd, rvec *f);
241
242 /*! \brief Move x0 and also x1 if x1!=NULL. bX1IsCoord tells if to do PBC on x1 */
243 void dd_move_x_constraints(struct gmx_domdec_t *dd, matrix box,
244                            rvec *x0, rvec *x1, 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 ]]
263 void dd_print_missing_interactions(const gmx::MDLogger  &mdlog,
264                                    t_commrec            *cr,
265                                    int                   local_count,
266                                    const gmx_mtop_t     *top_global,
267                                    const gmx_localtop_t *top_local,
268                                    const t_state        *state_local);
269
270 /*! \brief Generate and store the reverse topology */
271 void dd_make_reverse_top(FILE *fplog,
272                          gmx_domdec_t *dd, const gmx_mtop_t *mtop,
273                          const gmx_vsite_t *vsite,
274                          const t_inputrec *ir, gmx_bool bBCheck);
275
276 /*! \brief Store the local charge group index in \p lcgs */
277 void dd_make_local_cgs(struct gmx_domdec_t *dd, t_block *lcgs);
278
279 /*! \brief Generate the local topology and virtual site data */
280 void dd_make_local_top(struct gmx_domdec_t       *dd,
281                        struct gmx_domdec_zones_t *zones,
282                        int                        npbcdim,
283                        matrix                     box,
284                        rvec                       cellsize_min,
285                        const ivec                 npulse,
286                        t_forcerec                *fr,
287                        rvec                      *cgcm_or_x,
288                        const gmx_mtop_t          &top,
289                        gmx_localtop_t            *ltop);
290
291 /*! \brief Sort ltop->ilist when we are doing free energy. */
292 void dd_sort_local_top(gmx_domdec_t *dd, const t_mdatoms *mdatoms,
293                        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,
301                        gmx_localtop_t   *top);
302
303 /*! \brief Construct local state */
304 void dd_init_local_state(struct gmx_domdec_t *dd,
305                          const t_state *state_global, t_state *local_state);
306
307 /*! \brief Generate a list of links between charge groups that are linked by bonded interactions */
308 t_blocka *make_charge_group_links(const gmx_mtop_t *mtop, gmx_domdec_t *dd,
309                                   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, const matrix box,
316                            gmx_bool bBCheck,
317                            real *r_2b, real *r_mb);
318
319
320 /* In domdec_setup.c */
321
322 /*! \brief Returns the volume fraction of the system that is communicated */
323 real comm_box_frac(const ivec dd_nc, real cutoff, const gmx_ddbox_t *ddbox);
324
325 /*! \brief Determines the optimal DD cell setup dd->nc and possibly npmenodes
326  * for the system.
327  *
328  * On the master node returns the actual cellsize limit used.
329  */
330 real dd_choose_grid(const gmx::MDLogger &mdlog,
331                     t_commrec *cr, gmx_domdec_t *dd,
332                     const t_inputrec *ir,
333                     const gmx_mtop_t *mtop,
334                     const matrix box, const gmx_ddbox_t *ddbox,
335                     int nPmeRanks,
336                     gmx_bool bDynLoadBal, real dlb_scale,
337                     real cellsize_limit, real cutoff_dd,
338                     gmx_bool bInterCGBondeds);
339
340 #endif