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 RangePartitioning;
93 /*! \brief Returns the global topology atom number belonging to local atom index i.
95 * This function is intended for writing ASCII output
96 * and returns atom numbers starting at 1.
97 * When dd=NULL returns i+1.
99 int ddglatnr(const gmx_domdec_t* dd, int i);
101 /*! \brief Returns a list of update group partitioning for each molecule type or empty when update groups are not used */
102 gmx::ArrayRef<const gmx::RangePartitioning> getUpdateGroupingPerMoleculetype(const gmx_domdec_t& dd);
104 /*! \brief Store the global cg indices of the home cgs in state,
106 * This means it can be reset, even after a new DD partitioning.
108 void dd_store_state(struct gmx_domdec_t* dd, t_state* state);
110 /*! \brief Returns a pointer to the gmx_domdec_zones_t struct */
111 struct gmx_domdec_zones_t* domdec_zones(struct gmx_domdec_t* dd);
113 /*! \brief Returns the range for atoms in zones*/
114 int dd_numAtomsZones(const gmx_domdec_t& dd);
116 /*! \brief Returns the number of home atoms */
117 int dd_numHomeAtoms(const gmx_domdec_t& dd);
119 /*! \brief Returns the atom range in the local state for atoms that need to be present in mdatoms */
120 int dd_natoms_mdatoms(const gmx_domdec_t* dd);
122 /*! \brief Returns the atom range in the local state for atoms involved in virtual sites */
123 int dd_natoms_vsite(const gmx_domdec_t* dd);
125 /*! \brief Sets the atom range for atom in the local state for atoms received in constraints communication */
126 void dd_get_constraint_range(const gmx_domdec_t* dd, int* at_start, int* at_end);
128 /*! \libinternal \brief Struct for passing around the number of PME domains */
131 int x; //!< The number of PME domains along dimension x
132 int y; //!< The number of PME domains along dimension y
135 /*! \brief Returns the number of PME domains, can be called with dd=NULL */
136 NumPmeDomains getNumPmeDomains(const gmx_domdec_t* dd);
138 /*! \brief Returns the set of DD ranks that communicate with pme node cr->nodeid */
139 std::vector<int> get_pme_ddranks(const t_commrec* cr, int pmenodeid);
141 /*! \brief Returns the maximum shift for coordinate communication in PME, dim x */
142 int dd_pme_maxshift_x(const gmx_domdec_t* dd);
144 /*! \brief Returns the maximum shift for coordinate communication in PME, dim y */
145 int dd_pme_maxshift_y(const gmx_domdec_t* dd);
147 /*! \brief Return whether constraints, not including settles, cross domain boundaries */
148 bool ddHaveSplitConstraints(const gmx_domdec_t& dd);
150 /*! \brief Return whether the DD has a single dimension with a single pulse
152 * The GPU halo exchange code requires a 1D single-pulse DD, and its
153 * setup code can use the returned value to understand what it should
155 bool is1DAnd1PulseDD(const gmx_domdec_t& dd);
157 /*! \brief Initialize data structures for bonded interactions */
158 void dd_init_bondeds(FILE* fplog,
160 const gmx_mtop_t* mtop,
161 const gmx_vsite_t* vsite,
162 const t_inputrec* ir,
164 cginfo_mb_t* cginfo_mb);
166 /*! \brief Returns whether molecules are always whole, i.e. not broken by PBC */
167 bool dd_moleculesAreAlwaysWhole(const gmx_domdec_t& dd);
169 /*! \brief Returns if we need to do pbc for calculating bonded interactions */
170 gmx_bool dd_bonded_molpbc(const gmx_domdec_t* dd, int ePBC);
172 /*! \brief Change the DD non-bonded communication cut-off.
174 * This could fail when trying to increase the cut-off,
175 * then FALSE will be returned and the cut-off is not modified.
177 * \param[in] cr Communication recrod
178 * \param[in] box Box matrix, used for computing the dimensions of the system
179 * \param[in] x Position vector, used for computing the dimensions of the system
180 * \param[in] cutoffRequested The requested atom to atom cut-off distance, usually the pair-list cutoff distance
182 gmx_bool change_dd_cutoff(t_commrec* cr, const matrix box, gmx::ArrayRef<const gmx::RVec> x, real cutoffRequested);
184 /*! \brief Set up communication for averaging GPU wait times over domains
186 * When domains (PP MPI ranks) share a GPU, the individual GPU wait times
187 * are meaningless, as it depends on the order in which tasks on the same
188 * GPU finish. Therefore there wait times need to be averaged over the ranks
189 * sharing the same GPU. This function sets up the communication for that.
191 void dd_setup_dlb_resource_sharing(const t_commrec* cr, int gpu_id);
193 /*! \brief Cycle counter indices used internally in the domain decomposition */
204 /*! \brief Add the wallcycle count to the DD counter */
205 void dd_cycles_add(const gmx_domdec_t* dd, float cycles, int ddCycl);
207 /*! \brief Communicate the coordinates to the neighboring cells and do pbc. */
208 void dd_move_x(struct gmx_domdec_t* dd, const matrix box, gmx::ArrayRef<gmx::RVec> x, gmx_wallcycle* wcycle);
210 /*! \brief Sum the forces over the neighboring cells.
212 * When fshift!=NULL the shift forces are updated to obtain
213 * the correct virial from the single sum including f.
215 void dd_move_f(struct gmx_domdec_t* dd, gmx::ForceWithShiftForces* forceWithShiftForces, gmx_wallcycle* wcycle);
217 /*! \brief Communicate a real for each atom to the neighboring cells. */
218 void dd_atom_spread_real(struct gmx_domdec_t* dd, real v[]);
220 /*! \brief Sum the contributions to a real for each atom over the neighboring cells. */
221 void dd_atom_sum_real(struct gmx_domdec_t* dd, real v[]);
223 /*! \brief Reset all the statistics and counters for total run counting */
224 void reset_dd_statistics_counters(struct gmx_domdec_t* dd);
226 /* In domdec_con.c */
228 /*! \brief Communicates the virtual site forces, reduces the shift forces when \p fshift != NULL */
229 void dd_move_f_vsites(struct gmx_domdec_t* dd, rvec* f, rvec* fshift);
231 /*! \brief Clears the forces for virtual sites */
232 void dd_clear_f_vsites(struct gmx_domdec_t* dd, rvec* f);
234 /*! \brief Move x0 and also x1 if x1!=NULL. bX1IsCoord tells if to do PBC on x1 */
235 void dd_move_x_constraints(struct gmx_domdec_t* dd, const matrix box, rvec* x0, rvec* x1, gmx_bool bX1IsCoord);
237 /*! \brief Communicates the coordinates involved in virtual sites */
238 void dd_move_x_vsites(struct gmx_domdec_t* dd, const matrix box, rvec* x);
240 /*! \brief Returns the local atom count array for all constraints
242 * The local atom count for a constraint, possible values 2/1/0, is needed
243 * to avoid not/double-counting contributions linked to the Lagrange
244 * multiplier, such as the virial and free-energy derivatives.
246 * \note When \p dd = nullptr, an empty reference is returned.
248 gmx::ArrayRef<const int> dd_constraints_nlocalatoms(const gmx_domdec_t* dd);
250 /* In domdec_top.c */
252 /*! \brief Print error output when interactions are missing */
253 [[noreturn]] void dd_print_missing_interactions(const gmx::MDLogger& mdlog,
256 const gmx_mtop_t* top_global,
257 const gmx_localtop_t* top_local,
261 /*! \brief Generate and store the reverse topology */
262 void dd_make_reverse_top(FILE* fplog,
264 const gmx_mtop_t* mtop,
265 const gmx_vsite_t* vsite,
266 const t_inputrec* ir,
269 /*! \brief Generate the local topology and virtual site data */
270 void dd_make_local_top(struct gmx_domdec_t* dd,
271 struct gmx_domdec_zones_t* zones,
278 const gmx_mtop_t& top,
279 gmx_localtop_t* ltop);
281 /*! \brief Sort ltop->ilist when we are doing free energy. */
282 void dd_sort_local_top(gmx_domdec_t* dd, const t_mdatoms* mdatoms, gmx_localtop_t* ltop);
284 /*! \brief Initialize local topology
286 * \param[in] top_global Reference to global topology.
287 * \param[in,out] top Pointer to new local topology
289 void dd_init_local_top(const gmx_mtop_t& top_global, gmx_localtop_t* top);
291 /*! \brief Construct local state */
292 void dd_init_local_state(struct gmx_domdec_t* dd, const t_state* state_global, t_state* local_state);
294 /*! \brief Generate a list of links between atoms that are linked by bonded interactions
296 * Also stores whether atoms are linked in \p cginfo_mb.
298 t_blocka* makeBondedLinks(const gmx_mtop_t* mtop, cginfo_mb_t* cginfo_mb);
300 /*! \brief Calculate the maximum distance involved in 2-body and multi-body bonded interactions */
301 void dd_bonded_cg_distance(const gmx::MDLogger& mdlog,
302 const gmx_mtop_t* mtop,
303 const t_inputrec* ir,