Move constraint and bonded filtering info into DDSystemInfo
[gromacs.git] / src / gromacs / domdec / domdec.h
blob948d6ec55660560d801570207d59ec7c99542851
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
36 /*! \defgroup module_domdec Spatial domain decomposition (for parallelization over MPI)
37 * \ingroup group_mdrun
39 * \brief Manages the decomposition of the simulation volume over MPI
40 * ranks to try to distribute work evenly with minimal communication
41 * overheads.
43 * \todo Get domdec stuff out of mdtypes/commrec.h
45 * \author Berk Hess <hess@kth.se>
49 /*! \libinternal \file
51 * \brief This file declares functions for mdrun to call to manage the
52 * details of its domain decomposition.
54 * \author Berk Hess <hess@kth.se>
55 * \inlibraryapi
56 * \ingroup module_domdec
59 #ifndef GMX_DOMDEC_DOMDEC_H
60 #define GMX_DOMDEC_DOMDEC_H
62 #include <vector>
64 #include "gromacs/math/vectypes.h"
65 #include "gromacs/utility/arrayref.h"
66 #include "gromacs/utility/basedefinitions.h"
67 #include "gromacs/utility/real.h"
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;
87 namespace gmx
89 class ForceWithShiftForces;
90 class MDLogger;
91 class LocalAtomSetManager;
92 class RangePartitioning;
93 struct DomdecOptions;
94 struct MdrunOptions;
95 } // namespace
97 /*! \brief Returns the global topology atom number belonging to local atom index i.
99 * This function is intended for writing ASCII output
100 * and returns atom numbers starting at 1.
101 * When dd=NULL returns i+1.
103 int ddglatnr(const gmx_domdec_t *dd, int i);
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);
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);
111 /*! \brief Store the global cg indices of the home cgs in state,
113 * This means it can be reset, even after a new DD partitioning.
115 void dd_store_state(struct gmx_domdec_t *dd, t_state *state);
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);
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);
124 /*! \brief Returns the number of home atoms */
125 int dd_numHomeAtoms(const gmx_domdec_t &dd);
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);
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);
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);
137 /*! \libinternal \brief Struct for passing around the number of PME domains */
138 struct NumPmeDomains
140 int x; //!< The number of PME domains along dimension x
141 int y; //!< The number of PME domains along dimension y
144 /*! \brief Returns the number of PME domains, can be called with dd=NULL */
145 NumPmeDomains getNumPmeDomains(const gmx_domdec_t *dd);
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);
150 /*! \brief Returns the maximum shift for coordinate communication in PME, dim x */
151 int dd_pme_maxshift_x(const gmx_domdec_t *dd);
153 /*! \brief Returns the maximum shift for coordinate communication in PME, dim y */
154 int dd_pme_maxshift_y(const gmx_domdec_t *dd);
156 /*! \brief Return whether constraints, not including settles, cross domain boundaries */
157 bool ddHaveSplitConstraints(const gmx_domdec_t &dd);
159 /*! \brief Initialized the domain decomposition, chooses the DD grid and PME ranks, return the DD struct */
160 gmx_domdec_t *
161 init_domain_decomposition(const gmx::MDLogger &mdlog,
162 t_commrec *cr,
163 const gmx::DomdecOptions &options,
164 const gmx::MdrunOptions &mdrunOptions,
165 const gmx_mtop_t *mtop,
166 const t_inputrec *ir,
167 const matrix box,
168 gmx::ArrayRef<const gmx::RVec> xGlobal,
169 gmx::LocalAtomSetManager *atomSets);
171 /*! \brief Initialize data structures for bonded interactions */
172 void dd_init_bondeds(FILE *fplog,
173 gmx_domdec_t *dd,
174 const gmx_mtop_t *mtop,
175 const gmx_vsite_t *vsite,
176 const t_inputrec *ir,
177 gmx_bool bBCheck,
178 cginfo_mb_t *cginfo_mb);
180 /*! \brief Returns if we need to do pbc for calculating bonded interactions */
181 gmx_bool dd_bonded_molpbc(const gmx_domdec_t *dd, int ePBC);
183 /*! \brief Change the DD non-bonded communication cut-off.
185 * This could fail when trying to increase the cut-off,
186 * then FALSE will be returned and the cut-off is not modified.
188 * \param[in] cr Communication recrod
189 * \param[in] state State, used for computing the dimensions of the system
190 * \param[in] cutoffRequested The requested atom to atom cut-off distance, usually the pair-list cutoff distance
192 gmx_bool change_dd_cutoff(t_commrec *cr,
193 const t_state &state,
194 real cutoffRequested);
196 /*! \brief Set up communication for averaging GPU wait times over domains
198 * When domains (PP MPI ranks) share a GPU, the individual GPU wait times
199 * are meaningless, as it depends on the order in which tasks on the same
200 * GPU finish. Therefore there wait times need to be averaged over the ranks
201 * sharing the same GPU. This function sets up the communication for that.
203 void dd_setup_dlb_resource_sharing(t_commrec *cr,
204 int gpu_id);
206 /*! \brief Cycle counter indices used internally in the domain decomposition */
207 enum {
208 ddCyclStep, ddCyclPPduringPME, ddCyclF, ddCyclWaitGPU, ddCyclPME, ddCyclNr
211 /*! \brief Add the wallcycle count to the DD counter */
212 void dd_cycles_add(const gmx_domdec_t *dd, float cycles, int ddCycl);
214 /*! \brief Communicate the coordinates to the neighboring cells and do pbc. */
215 void dd_move_x(struct gmx_domdec_t *dd,
216 const matrix box,
217 gmx::ArrayRef<gmx::RVec> x,
218 gmx_wallcycle *wcycle);
220 /*! \brief Sum the forces over the neighboring cells.
222 * When fshift!=NULL the shift forces are updated to obtain
223 * the correct virial from the single sum including f.
225 void dd_move_f(struct gmx_domdec_t *dd,
226 gmx::ForceWithShiftForces *forceWithShiftForces,
227 gmx_wallcycle *wcycle);
229 /*! \brief Communicate a real for each atom to the neighboring cells. */
230 void dd_atom_spread_real(struct gmx_domdec_t *dd, real v[]);
232 /*! \brief Sum the contributions to a real for each atom over the neighboring cells. */
233 void dd_atom_sum_real(struct gmx_domdec_t *dd, real v[]);
235 /*! \brief Reset all the statistics and counters for total run counting */
236 void reset_dd_statistics_counters(struct gmx_domdec_t *dd);
238 /* In domdec_con.c */
240 /*! \brief Communicates the virtual site forces, reduces the shift forces when \p fshift != NULL */
241 void dd_move_f_vsites(struct gmx_domdec_t *dd, rvec *f, rvec *fshift);
243 /*! \brief Clears the forces for virtual sites */
244 void dd_clear_f_vsites(struct gmx_domdec_t *dd, rvec *f);
246 /*! \brief Move x0 and also x1 if x1!=NULL. bX1IsCoord tells if to do PBC on x1 */
247 void dd_move_x_constraints(struct gmx_domdec_t *dd, matrix box,
248 rvec *x0, rvec *x1, gmx_bool bX1IsCoord);
250 /*! \brief Communicates the coordinates involved in virtual sites */
251 void dd_move_x_vsites(struct gmx_domdec_t *dd, const matrix box, rvec *x);
253 /*! \brief Returns the local atom count array for all constraints
255 * The local atom count for a constraint, possible values 2/1/0, is needed
256 * to avoid not/double-counting contributions linked to the Lagrange
257 * multiplier, such as the virial and free-energy derivatives.
259 * \note When \p dd = nullptr, an empty reference is returned.
261 gmx::ArrayRef<const int> dd_constraints_nlocalatoms(const gmx_domdec_t *dd);
263 /* In domdec_top.c */
265 /*! \brief Print error output when interactions are missing */
266 [[ noreturn ]]
267 void dd_print_missing_interactions(const gmx::MDLogger &mdlog,
268 t_commrec *cr,
269 int local_count,
270 const gmx_mtop_t *top_global,
271 const gmx_localtop_t *top_local,
272 const rvec *x,
273 const matrix box);
275 /*! \brief Generate and store the reverse topology */
276 void dd_make_reverse_top(FILE *fplog,
277 gmx_domdec_t *dd, const gmx_mtop_t *mtop,
278 const gmx_vsite_t *vsite,
279 const t_inputrec *ir, gmx_bool bBCheck);
281 /*! \brief Generate the local topology and virtual site data */
282 void dd_make_local_top(struct gmx_domdec_t *dd,
283 struct gmx_domdec_zones_t *zones,
284 int npbcdim,
285 matrix box,
286 rvec cellsize_min,
287 const ivec npulse,
288 t_forcerec *fr,
289 rvec *cgcm_or_x,
290 const gmx_mtop_t &top,
291 gmx_localtop_t *ltop);
293 /*! \brief Sort ltop->ilist when we are doing free energy. */
294 void dd_sort_local_top(gmx_domdec_t *dd, const t_mdatoms *mdatoms,
295 gmx_localtop_t *ltop);
297 /*! \brief Initialize local topology
299 * \param[in] top_global Reference to global topology.
300 * \param[in,out] top Pointer to new local topology
302 void dd_init_local_top(const gmx_mtop_t &top_global,
303 gmx_localtop_t *top);
305 /*! \brief Construct local state */
306 void dd_init_local_state(struct gmx_domdec_t *dd,
307 const t_state *state_global, t_state *local_state);
309 /*! \brief Generate a list of links between charge groups that are linked by bonded interactions */
310 t_blocka *make_charge_group_links(const gmx_mtop_t *mtop, gmx_domdec_t *dd,
311 cginfo_mb_t *cginfo_mb);
313 /*! \brief Calculate the maximum distance involved in 2-body and multi-body bonded interactions */
314 void dd_bonded_cg_distance(const gmx::MDLogger &mdlog,
315 const gmx_mtop_t *mtop,
316 const t_inputrec *ir,
317 const rvec *x, const matrix box,
318 gmx_bool bBCheck,
319 real *r_2b, real *r_mb);
322 /* In domdec_setup.c */
324 /*! \brief Returns the volume fraction of the system that is communicated */
325 real comm_box_frac(const ivec dd_nc, real cutoff, const gmx_ddbox_t *ddbox);
327 struct DDSetup
329 int numPmeRanks = 0;
330 ivec numDomains = { 0, 0, 0 };
331 real cellsizeLimit = 0;
334 /*! \brief Determines the optimal DD cell setup dd->nc and possibly npmenodes
335 * for the system.
337 DDSetup
338 dd_choose_grid(const gmx::MDLogger &mdlog,
339 const t_commrec *cr,
340 const t_inputrec *ir,
341 const gmx_mtop_t *mtop,
342 const matrix box, const gmx_ddbox_t *ddbox,
343 int numPmeRanksRequested,
344 gmx_bool bDynLoadBal, real dlb_scale,
345 const DDSystemInfo &systemInfo);
347 #endif