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
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>
56 * \ingroup module_domdec
59 #ifndef GMX_DOMDEC_DOMDEC_H
60 #define GMX_DOMDEC_DOMDEC_H
64 #include "gromacs/math/vectypes.h"
65 #include "gromacs/utility/arrayref.h"
66 #include "gromacs/utility/basedefinitions.h"
67 #include "gromacs/utility/real.h"
72 struct gmx_domdec_zones_t;
73 struct gmx_localtop_t;
88 class ForceWithShiftForces;
90 class LocalAtomSetManager;
91 class RangePartitioning;
96 /*! \brief Returns the global topology atom number belonging to local atom index i.
98 * This function is intended for writing ASCII output
99 * and returns atom numbers starting at 1.
100 * When dd=NULL returns i+1.
102 int ddglatnr(const gmx_domdec_t *dd, int i);
104 /*! \brief Return a block struct for the charge groups of the whole system */
105 t_block *dd_charge_groups_global(struct gmx_domdec_t *dd);
107 /*! \brief Returns a list of update group partitioning for each molecule type or empty when update groups are not used */
108 gmx::ArrayRef<const gmx::RangePartitioning> getUpdateGroupingPerMoleculetype(const gmx_domdec_t &dd);
110 /*! \brief Store the global cg indices of the home cgs in state,
112 * This means it can be reset, even after a new DD partitioning.
114 void dd_store_state(struct gmx_domdec_t *dd, t_state *state);
116 /*! \brief Returns a pointer to the gmx_domdec_zones_t struct */
117 struct gmx_domdec_zones_t *domdec_zones(struct gmx_domdec_t *dd);
119 /*! \brief Sets the j-charge-group range for i-charge-group \p icg */
120 void dd_get_ns_ranges(const gmx_domdec_t *dd, int icg,
121 int *jcg0, int *jcg1, ivec shift0, ivec shift1);
123 /*! \brief Returns the number of home atoms */
124 int dd_numHomeAtoms(const gmx_domdec_t &dd);
126 /*! \brief Returns the atom range in the local state for atoms that need to be present in mdatoms */
127 int dd_natoms_mdatoms(const gmx_domdec_t *dd);
129 /*! \brief Returns the atom range in the local state for atoms involved in virtual sites */
130 int dd_natoms_vsite(const gmx_domdec_t *dd);
132 /*! \brief Sets the atom range for atom in the local state for atoms received in constraints communication */
133 void dd_get_constraint_range(const gmx_domdec_t *dd,
134 int *at_start, int *at_end);
136 /*! \libinternal \brief Struct for passing around the number of PME domains */
139 int x; //!< The number of PME domains along dimension x
140 int y; //!< The number of PME domains along dimension y
143 /*! \brief Returns the number of PME domains, can be called with dd=NULL */
144 NumPmeDomains getNumPmeDomains(const gmx_domdec_t *dd);
146 /*! \brief Returns the set of DD ranks that communicate with pme node cr->nodeid */
147 std::vector<int> get_pme_ddranks(const t_commrec *cr, int pmenodeid);
149 /*! \brief Returns the maximum shift for coordinate communication in PME, dim x */
150 int dd_pme_maxshift_x(const gmx_domdec_t *dd);
152 /*! \brief Returns the maximum shift for coordinate communication in PME, dim y */
153 int dd_pme_maxshift_y(const gmx_domdec_t *dd);
155 /*! \brief Initialized the domain decomposition, chooses the DD grid and PME ranks, return the DD struct */
157 init_domain_decomposition(const gmx::MDLogger &mdlog,
159 const gmx::DomdecOptions &options,
160 const gmx::MdrunOptions &mdrunOptions,
161 const gmx_mtop_t *mtop,
162 const t_inputrec *ir,
164 gmx::ArrayRef<const gmx::RVec> xGlobal,
165 gmx::LocalAtomSetManager *atomSets);
167 /*! \brief Initialize data structures for bonded interactions */
168 void dd_init_bondeds(FILE *fplog,
170 const gmx_mtop_t *mtop,
171 const gmx_vsite_t *vsite,
172 const t_inputrec *ir,
174 cginfo_mb_t *cginfo_mb);
176 /*! \brief Returns if we need to do pbc for calculating bonded interactions */
177 gmx_bool dd_bonded_molpbc(const gmx_domdec_t *dd, int ePBC);
179 /*! \brief Change the DD non-bonded communication cut-off.
181 * This could fail when trying to increase the cut-off,
182 * then FALSE will be returned and the cut-off is not modified.
184 * \param[in] cr Communication recrod
185 * \param[in] state State, used for computing the dimensions of the system
186 * \param[in] cutoffRequested The requested atom to atom cut-off distance, usually the pair-list cutoff distance
188 gmx_bool change_dd_cutoff(t_commrec *cr,
189 const t_state &state,
190 real cutoffRequested);
192 /*! \brief Set up communication for averaging GPU wait times over domains
194 * When domains (PP MPI ranks) share a GPU, the individual GPU wait times
195 * are meaningless, as it depends on the order in which tasks on the same
196 * GPU finish. Therefore there wait times need to be averaged over the ranks
197 * sharing the same GPU. This function sets up the communication for that.
199 void dd_setup_dlb_resource_sharing(t_commrec *cr,
202 /*! \brief Cycle counter indices used internally in the domain decomposition */
204 ddCyclStep, ddCyclPPduringPME, ddCyclF, ddCyclWaitGPU, ddCyclPME, ddCyclNr
207 /*! \brief Add the wallcycle count to the DD counter */
208 void dd_cycles_add(const gmx_domdec_t *dd, float cycles, int ddCycl);
210 /*! \brief Communicate the coordinates to the neighboring cells and do pbc. */
211 void dd_move_x(struct gmx_domdec_t *dd,
213 gmx::ArrayRef<gmx::RVec> x,
214 gmx_wallcycle *wcycle);
216 /*! \brief Sum the forces over the neighboring cells.
218 * When fshift!=NULL the shift forces are updated to obtain
219 * the correct virial from the single sum including f.
221 void dd_move_f(struct gmx_domdec_t *dd,
222 gmx::ForceWithShiftForces *forceWithShiftForces,
223 gmx_wallcycle *wcycle);
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[]);
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[]);
231 /*! \brief Reset all the statistics and counters for total run counting */
232 void reset_dd_statistics_counters(struct gmx_domdec_t *dd);
234 /* In domdec_con.c */
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);
239 /*! \brief Clears the forces for virtual sites */
240 void dd_clear_f_vsites(struct gmx_domdec_t *dd, rvec *f);
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);
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);
249 /*! \brief Returns the local atom count array for all constraints
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.
255 * \note When \p dd = nullptr, an empty reference is returned.
257 gmx::ArrayRef<const int> dd_constraints_nlocalatoms(const gmx_domdec_t *dd);
259 /* In domdec_top.c */
261 /*! \brief Print error output when interactions are missing */
263 void dd_print_missing_interactions(const gmx::MDLogger &mdlog,
266 const gmx_mtop_t *top_global,
267 const gmx_localtop_t *top_local,
271 /*! \brief Generate and store the reverse topology */
272 void dd_make_reverse_top(FILE *fplog,
273 gmx_domdec_t *dd, const gmx_mtop_t *mtop,
274 const gmx_vsite_t *vsite,
275 const t_inputrec *ir, gmx_bool bBCheck);
277 /*! \brief Generate the local topology and virtual site data */
278 void dd_make_local_top(struct gmx_domdec_t *dd,
279 struct gmx_domdec_zones_t *zones,
286 const gmx_mtop_t &top,
287 gmx_localtop_t *ltop);
289 /*! \brief Sort ltop->ilist when we are doing free energy. */
290 void dd_sort_local_top(gmx_domdec_t *dd, const t_mdatoms *mdatoms,
291 gmx_localtop_t *ltop);
293 /*! \brief Initialize local topology
295 * \param[in] top_global Reference to global topology.
296 * \param[in,out] top Pointer to new local topology
298 void dd_init_local_top(const gmx_mtop_t &top_global,
299 gmx_localtop_t *top);
301 /*! \brief Construct local state */
302 void dd_init_local_state(struct gmx_domdec_t *dd,
303 const t_state *state_global, t_state *local_state);
305 /*! \brief Generate a list of links between charge groups that are linked by bonded interactions */
306 t_blocka *make_charge_group_links(const gmx_mtop_t *mtop, gmx_domdec_t *dd,
307 cginfo_mb_t *cginfo_mb);
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, const matrix box,
315 real *r_2b, real *r_mb);
318 /* In domdec_setup.c */
320 /*! \brief Returns the volume fraction of the system that is communicated */
321 real comm_box_frac(const ivec dd_nc, real cutoff, const gmx_ddbox_t *ddbox);
323 /*! \brief Determines the optimal DD cell setup dd->nc and possibly npmenodes
326 * On the master node returns the actual cellsize limit used.
328 real dd_choose_grid(const gmx::MDLogger &mdlog,
329 t_commrec *cr, gmx_domdec_t *dd,
330 const t_inputrec *ir,
331 const gmx_mtop_t *mtop,
332 const matrix box, const gmx_ddbox_t *ddbox,
334 gmx_bool bDynLoadBal, real dlb_scale,
335 real cellsize_limit, real cutoff_dd,
336 gmx_bool bInterCGBondeds);