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;
90 class LocalAtomSetManager;
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 Return a block struct for the charge groups of the whole system */
102 t_block *dd_charge_groups_global(struct 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 Sets the j-charge-group range for i-charge-group \p icg */
114 void dd_get_ns_ranges(const gmx_domdec_t *dd, int icg,
115 int *jcg0, int *jcg1, ivec shift0, ivec shift1);
117 /*! \brief Returns the number of home atoms */
118 int dd_numHomeAtoms(const gmx_domdec_t &dd);
120 /*! \brief Returns the atom range in the local state for atoms that need to be present in mdatoms */
121 int dd_natoms_mdatoms(const gmx_domdec_t *dd);
123 /*! \brief Returns the atom range in the local state for atoms involved in virtual sites */
124 int dd_natoms_vsite(const gmx_domdec_t *dd);
126 /*! \brief Sets the atom range for atom in the local state for atoms received in constraints communication */
127 void dd_get_constraint_range(const gmx_domdec_t *dd,
128 int *at_start, int *at_end);
130 /*! \libinternal \brief Struct for passing around the number of PME domains */
133 int x; //!< The number of PME domains along dimension x
134 int y; //!< The number of PME domains along dimension y
137 /*! \brief Returns the number of PME domains, can be called with dd=NULL */
138 NumPmeDomains getNumPmeDomains(const gmx_domdec_t *dd);
140 /*! \brief Returns the set of DD ranks that communicate with pme node cr->nodeid */
141 std::vector<int> get_pme_ddranks(const t_commrec *cr, int pmenodeid);
143 /*! \brief Returns the maximum shift for coordinate communication in PME, dim x */
144 int dd_pme_maxshift_x(const gmx_domdec_t *dd);
146 /*! \brief Returns the maximum shift for coordinate communication in PME, dim y */
147 int dd_pme_maxshift_y(const gmx_domdec_t *dd);
149 /*! \brief The options for the domain decomposition MPI task ordering. */
150 enum class DdRankOrder
152 select, //!< First value (needed to cope with command-line parsing)
153 interleave, //!< Interleave the PP and PME ranks
154 pp_pme, //!< First all PP ranks, all PME rank at the end
155 cartesian, //!< Use Cartesian communicators for PP, PME and PP-PME
156 nr //!< The number of options
159 /*! \brief The options for the dynamic load balancing. */
162 select, //!< First value (needed to cope with command-line parsing)
163 turnOnWhenUseful, //!< Turn on DLB when we think it would improve performance
164 no, //!< Never turn on DLB
165 yes, //!< Turn on DLB from the start and keep it on
166 nr //!< The number of options
169 /*! \libinternal \brief Structure containing all (command line) options for the domain decomposition */
172 //! If true, check that all bonded interactions have been assigned to exactly one domain/rank.
173 gmx_bool checkBondedInteractions = TRUE;
174 //! If true, don't communicate all atoms between the non-bonded cut-off and the larger bonded cut-off, but only those that have non-local bonded interactions. This significantly reduces the communication volume.
175 gmx_bool useBondedCommunication = TRUE;
176 //! The domain decomposition grid cell count, 0 means let domdec choose based on the number of ranks.
178 //! The number of separate PME ranks requested, -1 = auto.
179 int numPmeRanks = -1;
180 //! Ordering of the PP and PME ranks, values from enum above.
181 DdRankOrder rankOrder = DdRankOrder::interleave;
182 //! The minimum communication range, used for extended the communication range for bonded interactions (nm).
183 real minimumCommunicationRange = 0;
184 //! Communication range for atom involved in constraints (P-LINCS) (nm).
185 real constraintCommunicationRange = 0;
186 //! Dynamic load balancing option, values from enum above.
187 DlbOption dlbOption = DlbOption::turnOnWhenUseful;
188 /*! \brief Fraction in (0,1) by whose reciprocal the initial
189 * DD cell size will be increased in order to provide a margin
190 * in which dynamic load balancing can act, while preserving
191 * the minimum cell size. */
192 real dlbScaling = 0.8;
193 //! String containing a vector of the relative sizes in the x direction of the corresponding DD cells.
194 const char *cellSizeX = nullptr;
195 //! String containing a vector of the relative sizes in the y direction of the corresponding DD cells.
196 const char *cellSizeY = nullptr;
197 //! String containing a vector of the relative sizes in the z direction of the corresponding DD cells.
198 const char *cellSizeZ = nullptr;
201 /*! \brief Initialized the domain decomposition, chooses the DD grid and PME ranks, return the DD struct */
203 init_domain_decomposition(const gmx::MDLogger &mdlog,
205 const DomdecOptions &options,
206 const gmx::MdrunOptions &mdrunOptions,
207 const gmx_mtop_t *mtop,
208 const t_inputrec *ir,
210 gmx::ArrayRef<const gmx::RVec> xGlobal,
211 gmx::LocalAtomSetManager *atomSets);
213 /*! \brief Initialize data structures for bonded interactions */
214 void dd_init_bondeds(FILE *fplog,
216 const gmx_mtop_t *mtop,
217 const gmx_vsite_t *vsite,
218 const t_inputrec *ir,
220 cginfo_mb_t *cginfo_mb);
222 /*! \brief Returns if we need to do pbc for calculating bonded interactions */
223 gmx_bool dd_bonded_molpbc(const gmx_domdec_t *dd, int ePBC);
225 /*! \brief Change the DD non-bonded communication cut-off.
227 * This could fail when trying to increase the cut-off,
228 * then FALSE will be returned and the cut-off is not modified.
230 * \param[in] cr Communication recrod
231 * \param[in] state State, used for computing the dimensions of the system
232 * \param[in] cutoffRequested The requested atom to atom cut-off distance, usually the pair-list cutoff distance
234 gmx_bool change_dd_cutoff(t_commrec *cr,
235 const t_state &state,
236 real cutoffRequested);
238 /*! \brief Set up communication for averaging GPU wait times over domains
240 * When domains (PP MPI ranks) share a GPU, the individual GPU wait times
241 * are meaningless, as it depends on the order in which tasks on the same
242 * GPU finish. Therefore there wait times need to be averaged over the ranks
243 * sharing the same GPU. This function sets up the communication for that.
245 void dd_setup_dlb_resource_sharing(t_commrec *cr,
248 /*! \brief Cycle counter indices used internally in the domain decomposition */
250 ddCyclStep, ddCyclPPduringPME, ddCyclF, ddCyclWaitGPU, ddCyclPME, ddCyclNr
253 /*! \brief Add the wallcycle count to the DD counter */
254 void dd_cycles_add(const gmx_domdec_t *dd, float cycles, int ddCycl);
256 /*! \brief Communicate the coordinates to the neighboring cells and do pbc. */
257 void dd_move_x(struct gmx_domdec_t *dd,
259 gmx::ArrayRef<gmx::RVec> x,
260 gmx_wallcycle *wcycle);
262 /*! \brief Sum the forces over the neighboring cells.
264 * When fshift!=NULL the shift forces are updated to obtain
265 * the correct virial from the single sum including f.
267 void dd_move_f(struct gmx_domdec_t *dd,
268 gmx::ArrayRef<gmx::RVec> f,
270 gmx_wallcycle *wcycle);
272 /*! \brief Communicate a real for each atom to the neighboring cells. */
273 void dd_atom_spread_real(struct gmx_domdec_t *dd, real v[]);
275 /*! \brief Sum the contributions to a real for each atom over the neighboring cells. */
276 void dd_atom_sum_real(struct gmx_domdec_t *dd, real v[]);
278 /*! \brief Reset all the statistics and counters for total run counting */
279 void reset_dd_statistics_counters(struct gmx_domdec_t *dd);
281 /* In domdec_con.c */
283 /*! \brief Communicates the virtual site forces, reduces the shift forces when \p fshift != NULL */
284 void dd_move_f_vsites(struct gmx_domdec_t *dd, rvec *f, rvec *fshift);
286 /*! \brief Clears the forces for virtual sites */
287 void dd_clear_f_vsites(struct gmx_domdec_t *dd, rvec *f);
289 /*! \brief Move x0 and also x1 if x1!=NULL. bX1IsCoord tells if to do PBC on x1 */
290 void dd_move_x_constraints(struct gmx_domdec_t *dd, matrix box,
291 rvec *x0, rvec *x1, gmx_bool bX1IsCoord);
293 /*! \brief Communicates the coordinates involved in virtual sites */
294 void dd_move_x_vsites(struct gmx_domdec_t *dd, const matrix box, rvec *x);
296 /*! \brief Returns the local atom count array for all constraints
298 * The local atom count for a constraint, possible values 2/1/0, is needed
299 * to avoid not/double-counting contributions linked to the Lagrange
300 * multiplier, such as the virial and free-energy derivatives.
302 * \note When \p dd = nullptr, an empty reference is returned.
304 gmx::ArrayRef<const int> dd_constraints_nlocalatoms(const gmx_domdec_t *dd);
306 /* In domdec_top.c */
308 /*! \brief Print error output when interactions are missing */
310 void dd_print_missing_interactions(const gmx::MDLogger &mdlog,
313 const gmx_mtop_t *top_global,
314 const gmx_localtop_t *top_local,
315 const t_state *state_local);
317 /*! \brief Generate and store the reverse topology */
318 void dd_make_reverse_top(FILE *fplog,
319 gmx_domdec_t *dd, const gmx_mtop_t *mtop,
320 const gmx_vsite_t *vsite,
321 const t_inputrec *ir, gmx_bool bBCheck);
323 /*! \brief Store the local charge group index in \p lcgs */
324 void dd_make_local_cgs(struct gmx_domdec_t *dd, t_block *lcgs);
326 /*! \brief Generate the local topology and virtual site data */
327 void dd_make_local_top(struct gmx_domdec_t *dd,
328 struct gmx_domdec_zones_t *zones,
336 const gmx_mtop_t &top,
337 gmx_localtop_t *ltop);
339 /*! \brief Sort ltop->ilist when we are doing free energy. */
340 void dd_sort_local_top(gmx_domdec_t *dd, const t_mdatoms *mdatoms,
341 gmx_localtop_t *ltop);
343 /*! \brief Initialize local topology
345 * \param[in] top_global Reference to global topology.
346 * \param[in,out] top Pointer to new local topology
348 void dd_init_local_top(const gmx_mtop_t &top_global,
349 gmx_localtop_t *top);
351 /*! \brief Construct local state */
352 void dd_init_local_state(struct gmx_domdec_t *dd,
353 const t_state *state_global, t_state *local_state);
355 /*! \brief Generate a list of links between charge groups that are linked by bonded interactions */
356 t_blocka *make_charge_group_links(const gmx_mtop_t *mtop, gmx_domdec_t *dd,
357 cginfo_mb_t *cginfo_mb);
359 /*! \brief Calculate the maximum distance involved in 2-body and multi-body bonded interactions */
360 void dd_bonded_cg_distance(const gmx::MDLogger &mdlog,
361 const gmx_mtop_t *mtop,
362 const t_inputrec *ir,
363 const rvec *x, const matrix box,
365 real *r_2b, real *r_mb);
368 /* In domdec_setup.c */
370 /*! \brief Returns the volume fraction of the system that is communicated */
371 real comm_box_frac(const ivec dd_nc, real cutoff, const gmx_ddbox_t *ddbox);
373 /*! \brief Determines the optimal DD cell setup dd->nc and possibly npmenodes
376 * On the master node returns the actual cellsize limit used.
378 real dd_choose_grid(const gmx::MDLogger &mdlog,
379 t_commrec *cr, gmx_domdec_t *dd,
380 const t_inputrec *ir,
381 const gmx_mtop_t *mtop,
382 const matrix box, const gmx_ddbox_t *ddbox,
384 gmx_bool bDynLoadBal, real dlb_scale,
385 real cellsize_limit, real cutoff_dd,
386 gmx_bool bInterCGBondeds);