2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2005,2006,2007,2008,2009,2010,2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
35 /*! \defgroup module_domdec Spatial domain decomposition (for parallelization over MPI)
36 * \ingroup group_mdrun
38 * \brief Manages the decomposition of the simulation volume over MPI
39 * ranks to try to distribute work evenly with minimal communication
42 * \todo Get domdec stuff out of legacyheaders/types/commrec.h
44 * \author Berk Hess <hess@kth.se>
48 /*! \libinternal \file
50 * \brief This file declares functions for mdrun to call to manage the
51 * details of its domain decomposition.
53 * \author Berk Hess <hess@kth.se>
55 * \ingroup module_domdec
58 #ifndef GMX_DOMDEC_DOMDEC_H
59 #define GMX_DOMDEC_DOMDEC_H
63 #include "gromacs/legacyheaders/vsite.h"
64 #include "gromacs/legacyheaders/types/commrec_fwd.h"
65 #include "gromacs/legacyheaders/types/constr.h"
66 #include "gromacs/legacyheaders/types/forcerec.h"
67 #include "gromacs/legacyheaders/types/hw_info.h"
68 #include "gromacs/legacyheaders/types/inputrec.h"
69 #include "gromacs/legacyheaders/types/mdatom.h"
70 #include "gromacs/legacyheaders/types/nrnb.h"
71 #include "gromacs/legacyheaders/types/shellfc.h"
72 #include "gromacs/legacyheaders/types/state.h"
73 #include "gromacs/math/vectypes.h"
74 #include "gromacs/timing/wallcycle.h"
75 #include "gromacs/topology/block.h"
76 #include "gromacs/topology/idef.h"
77 #include "gromacs/topology/topology.h"
78 #include "gromacs/utility/basedefinitions.h"
79 #include "gromacs/utility/real.h"
85 /*! \brief Returns the global topology atom number belonging to local atom index i.
87 * This function is intended for writing ASCII output
88 * and returns atom numbers starting at 1.
89 * When dd=NULL returns i+1.
91 int ddglatnr(gmx_domdec_t *dd, int i);
93 /*! \brief Return a block struct for the charge groups of the whole system */
94 t_block *dd_charge_groups_global(gmx_domdec_t *dd);
96 /*! \brief Store the global cg indices of the home cgs in state,
98 * This means it can be reset, even after a new DD partitioning.
100 void dd_store_state(gmx_domdec_t *dd, t_state *state);
102 /*! \brief Returns a pointer to the gmx_domdec_zones_t struct */
103 gmx_domdec_zones_t *domdec_zones(gmx_domdec_t *dd);
105 /*! \brief Sets the j-charge-group range for i-charge-group \p icg */
106 void dd_get_ns_ranges(gmx_domdec_t *dd, int icg,
107 int *jcg0, int *jcg1, ivec shift0, ivec shift1);
109 /*! \brief Returns the atom range in the local state for atoms involved in virtual sites */
110 int dd_natoms_vsite(gmx_domdec_t *dd);
112 /*! \brief Sets the atom range for atom in the local state for atoms received in constraints communication */
113 void dd_get_constraint_range(gmx_domdec_t *dd,
114 int *at_start, int *at_end);
116 /*! \brief Get the number of PME nodes along x and y, can be called with dd=NULL */
117 void get_pme_nnodes(const gmx_domdec_t *dd,
118 int *npmenodes_x, int *npmenodes_y);
120 /*! \brief Returns the set of DD nodes that communicate with pme node cr->nodeid */
121 void get_pme_ddnodes(t_commrec *cr, int pmenodeid,
122 int *nmy_ddnodes, int **my_ddnodes, int *node_peer);
124 /*! \brief Returns the maximum shift for coordinate communication in PME, dim x */
125 int dd_pme_maxshift_x(gmx_domdec_t *dd);
127 /*! \brief Returns the maximum shift for coordinate communication in PME, dim y */
128 int dd_pme_maxshift_y(gmx_domdec_t *dd);
130 /*! \brief Generates the MPI communicators for domain decomposition */
131 void make_dd_communicators(FILE *fplog, t_commrec *cr, int dd_node_order);
133 /*! \brief Initialized the domain decomposition, chooses the DD grid and PME ranks */
135 init_domain_decomposition(FILE *fplog,
139 real comm_distance_min, real rconstr,
140 const char *dlb_opt, real dlb_scale,
141 const char *sizex, const char *sizey, const char *sizez,
142 gmx_mtop_t *mtop, t_inputrec *ir,
145 int *npme_x, int *npme_y);
147 /*! \brief Initialize data structures for bonded interactions */
148 void dd_init_bondeds(FILE *fplog,
149 gmx_domdec_t *dd, gmx_mtop_t *mtop,
151 t_inputrec *ir, gmx_bool bBCheck, cginfo_mb_t *cginfo_mb);
153 /*! \brief Returns if we need to do pbc for calculating bonded interactions */
154 gmx_bool dd_bonded_molpbc(gmx_domdec_t *dd, int ePBC);
156 /*! \brief Set DD grid dimensions and limits
158 * Should be called after calling dd_init_bondeds.
160 void set_dd_parameters(FILE *fplog, gmx_domdec_t *dd, real dlb_scale,
164 /*! \brief Change the DD non-bonded communication cut-off.
166 * This could fail when trying to increase the cut-off,
167 * then FALSE will be returned and the cut-off is not modified.
169 gmx_bool change_dd_cutoff(t_commrec *cr, t_state *state, t_inputrec *ir,
172 /*! \brief Limit DLB to preserve the option of returning to the current cut-off.
174 * Domain boundary changes due to the DD dynamic load balancing can limit
175 * the cut-off distance that can be set in change_dd_cutoff. This function
176 * limits the DLB such that using the currently set cut-off should still be
177 * possible after subsequently setting a shorter cut-off with change_dd_cutoff.
179 void change_dd_dlb_cutoff_limit(t_commrec *cr);
181 /*! \brief Return if the DLB lock is set */
182 gmx_bool dd_dlb_is_locked(const gmx_domdec_t *dd);
184 /*! \brief Set a lock such that with DLB=auto DLB can (not) get turned on */
185 void dd_dlb_set_lock(gmx_domdec_t *dd, gmx_bool bValue);
187 /*! \brief Set up communication for averaging GPU wait times over ranks
189 * When domains (PP MPI ranks) share a GPU, the individual GPU wait times
190 * are meaningless, as it depends on the order in which tasks on the same
191 * GPU finish. Therefore there wait times need to be averaged over the ranks
192 * sharing the same GPU. This function sets up the communication for that.
194 void dd_setup_dlb_resource_sharing(t_commrec *cr,
195 const gmx_hw_info_t *hwinfo,
196 const gmx_hw_opt_t *hw_opt);
198 /*! \brief Sets up the DD communication setup */
199 void setup_dd_grid(FILE *fplog, gmx_domdec_t *dd);
201 /*! \brief Collects local rvec arrays \p lv to \p v on the master rank */
202 void dd_collect_vec(gmx_domdec_t *dd,
203 t_state *state_local, rvec *lv, rvec *v);
205 /*! \brief Collects the local state \p state_local to \p state on the master rank */
206 void dd_collect_state(gmx_domdec_t *dd,
207 t_state *state_local, t_state *state);
209 /*! \brief Cycle counter indices used internally in the domain decomposition */
211 ddCyclStep, ddCyclPPduringPME, ddCyclF, ddCyclWaitGPU, ddCyclPME, ddCyclNr
214 /*! \brief Add the wallcycle count to the DD counter */
215 void dd_cycles_add(gmx_domdec_t *dd, float cycles, int ddCycl);
217 /*! \brief Start the force flop count */
218 void dd_force_flop_start(gmx_domdec_t *dd, t_nrnb *nrnb);
220 /*! \brief Stop the force flop count */
221 void dd_force_flop_stop(gmx_domdec_t *dd, t_nrnb *nrnb);
223 /*! \brief Return the PME/PP force load ratio, or -1 if nothing was measured.
225 * Should only be called on the DD master node.
227 float dd_pme_f_ratio(gmx_domdec_t *dd);
229 /*! \brief Communicate the coordinates to the neighboring cells and do pbc. */
230 void dd_move_x(gmx_domdec_t *dd, matrix box, rvec x[]);
232 /*! \brief Sum the forces over the neighboring cells.
234 * When fshift!=NULL the shift forces are updated to obtain
235 * the correct virial from the single sum including f.
237 void dd_move_f(gmx_domdec_t *dd, rvec f[], rvec *fshift);
239 /*! \brief Communicate a real for each atom to the neighboring cells. */
240 void dd_atom_spread_real(gmx_domdec_t *dd, real v[]);
242 /*! \brief Sum the contributions to a real for each atom over the neighboring cells. */
243 void dd_atom_sum_real(gmx_domdec_t *dd, real v[]);
245 /*! \brief Partition the system over the nodes.
247 * step is only used for printing error messages.
248 * If bMasterState==TRUE then state_global from the master node is used,
249 * else state_local is redistributed between the nodes.
250 * When f!=NULL, *f will be reallocated to the size of state_local.
252 void dd_partition_system(FILE *fplog,
255 gmx_bool bMasterState,
257 t_state *state_global,
258 gmx_mtop_t *top_global,
260 t_state *state_local,
263 gmx_localtop_t *top_local,
266 gmx_shellfc_t shellfc,
269 gmx_wallcycle_t wcycle,
272 /*! \brief Reset all the statistics and counters for total run counting */
273 void reset_dd_statistics_counters(gmx_domdec_t *dd);
275 /*! \brief Print statistics for domain decomposition communication */
276 void print_dd_statistics(t_commrec *cr, t_inputrec *ir, FILE *fplog);
278 /* In domdec_con.c */
280 /*! \brief Communicates the virtual site forces, reduces the shift forces when \p fshift != NULL */
281 void dd_move_f_vsites(gmx_domdec_t *dd, rvec *f, rvec *fshift);
283 /*! \brief Clears the forces for virtual sites */
284 void dd_clear_f_vsites(gmx_domdec_t *dd, rvec *f);
286 /*! \brief Move x0 and also x1 if x1!=NULL. bX1IsCoord tells if to do PBC on x1 */
287 void dd_move_x_constraints(gmx_domdec_t *dd, matrix box,
288 rvec *x0, rvec *x1, gmx_bool bX1IsCoord);
290 /*! \brief Communicates the coordinates involved in virtual sites */
291 void dd_move_x_vsites(gmx_domdec_t *dd, matrix box, rvec *x);
293 /*! \brief Returns the local atom count array for all constraints
295 * The local atom count for a constraint, possible values 2/1/0, is needed
296 * to avoid not/double-counting contributions linked to the Lagrange
297 * multiplier, such as the virial and free-energy derivatives.
299 int *dd_constraints_nlocalatoms(gmx_domdec_t *dd);
301 /* In domdec_top.c */
303 /*! \brief Print error output when interactions are missing */
304 void dd_print_missing_interactions(FILE *fplog, t_commrec *cr,
305 int local_count, gmx_mtop_t *top_global, t_state *state_local);
307 /*! \brief Generate and store the reverse topology */
308 void dd_make_reverse_top(FILE *fplog,
309 gmx_domdec_t *dd, gmx_mtop_t *mtop,
311 t_inputrec *ir, gmx_bool bBCheck);
313 /*! \brief Store the local charge group index in \p lcgs */
314 void dd_make_local_cgs(gmx_domdec_t *dd, t_block *lcgs);
316 /*! \brief Generate the local topology and virtual site data */
317 void dd_make_local_top(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
318 int npbcdim, matrix box,
319 rvec cellsize_min, ivec npulse,
323 gmx_mtop_t *top, gmx_localtop_t *ltop);
325 /*! \brief Sort ltop->ilist when we are doing free energy. */
326 void dd_sort_local_top(gmx_domdec_t *dd, t_mdatoms *mdatoms,
327 gmx_localtop_t *ltop);
329 /*! \brief Construct local topology */
330 gmx_localtop_t *dd_init_local_top(gmx_mtop_t *top_global);
332 /*! \brief Construct local state */
333 void dd_init_local_state(gmx_domdec_t *dd,
334 t_state *state_global, t_state *local_state);
336 /*! \brief Generate a list of links between charge groups that are linked by bonded interactions */
337 t_blocka *make_charge_group_links(gmx_mtop_t *mtop, gmx_domdec_t *dd,
338 cginfo_mb_t *cginfo_mb);
340 /*! \brief Calculate the maximum distance involved in 2-body and multi-body bonded interactions */
341 void dd_bonded_cg_distance(FILE *fplog, gmx_mtop_t *mtop,
342 t_inputrec *ir, rvec *x, matrix box,
344 real *r_2b, real *r_mb);
346 /*! \brief Dump a pdb file with the current DD home + communicated atoms.
348 * When natoms=-1, dump all known atoms.
350 void write_dd_pdb(const char *fn, gmx_int64_t step, const char *title,
353 int natoms, rvec x[], matrix box);
356 /* In domdec_setup.c */
358 /*! \brief Returns the volume fraction of the system that is communicated */
359 real comm_box_frac(ivec dd_nc, real cutoff, gmx_ddbox_t *ddbox);
361 /*! \brief Determines the optimal DD cell setup dd->nc and possibly npmenodes
364 * On the master node returns the actual cellsize limit used.
366 real dd_choose_grid(FILE *fplog,
367 t_commrec *cr, gmx_domdec_t *dd, t_inputrec *ir,
368 gmx_mtop_t *mtop, matrix box, gmx_ddbox_t *ddbox,
369 gmx_bool bDynLoadBal, real dlb_scale,
370 real cellsize_limit, real cutoff_dd,
371 gmx_bool bInterCGBondeds);
374 /* In domdec_box.c */
376 /*! \brief Set the box and PBC data in \p ddbox */
377 void set_ddbox(gmx_domdec_t *dd, gmx_bool bMasterState, t_commrec *cr_sum,
378 t_inputrec *ir, matrix box,
379 gmx_bool bCalcUnboundedSize, t_block *cgs, rvec *x,
382 /*! \brief Set the box and PBC data in \p ddbox */
383 void set_ddbox_cr(t_commrec *cr, ivec *dd_nc,
384 t_inputrec *ir, matrix box, t_block *cgs, rvec *x,