Improve DLB+PME tuning with GPUs
[alexxy/gromacs.git] / src / gromacs / domdec / domdec.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2005,2006,2007,2008,2009,2010,2012,2013,2014,2015, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \defgroup module_domdec Spatial domain decomposition (for parallelization over MPI)
36  * \ingroup group_mdrun
37  *
38  * \brief Manages the decomposition of the simulation volume over MPI
39  * ranks to try to distribute work evenly with minimal communication
40  * overheads.
41  *
42  * \todo Get domdec stuff out of legacyheaders/types/commrec.h
43  *
44  * \author Berk Hess <hess@kth.se>
45  *
46  */
47
48 /*! \libinternal \file
49  *
50  * \brief This file declares functions for mdrun to call to manage the
51  * details of its domain decomposition.
52  *
53  * \author Berk Hess <hess@kth.se>
54  * \inlibraryapi
55  * \ingroup module_domdec
56  */
57
58 #ifndef GMX_DOMDEC_DOMDEC_H
59 #define GMX_DOMDEC_DOMDEC_H
60
61 #include <stdio.h>
62
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"
80
81 #ifdef __cplusplus
82 extern "C" {
83 #endif
84
85 /*! \brief Returns the global topology atom number belonging to local atom index i.
86  *
87  * This function is intended for writing ASCII output
88  * and returns atom numbers starting at 1.
89  * When dd=NULL returns i+1.
90  */
91 int ddglatnr(gmx_domdec_t *dd, int i);
92
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);
95
96 /*! \brief Store the global cg indices of the home cgs in state,
97  *
98  * This means it can be reset, even after a new DD partitioning.
99  */
100 void dd_store_state(gmx_domdec_t *dd, t_state *state);
101
102 /*! \brief Returns a pointer to the gmx_domdec_zones_t struct */
103 gmx_domdec_zones_t *domdec_zones(gmx_domdec_t *dd);
104
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);
108
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);
111
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);
115
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);
119
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);
123
124 /*! \brief Returns the maximum shift for coordinate communication in PME, dim x */
125 int dd_pme_maxshift_x(gmx_domdec_t *dd);
126
127 /*! \brief Returns the maximum shift for coordinate communication in PME, dim y */
128 int dd_pme_maxshift_y(gmx_domdec_t *dd);
129
130 /*! \brief Generates the MPI communicators for domain decomposition */
131 void make_dd_communicators(FILE *fplog, t_commrec *cr, int dd_node_order);
132
133 /*! \brief Initialized the domain decomposition, chooses the DD grid and PME ranks */
134 gmx_domdec_t *
135 init_domain_decomposition(FILE *fplog,
136                           t_commrec *cr,
137                           unsigned long Flags,
138                           ivec nc,
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,
143                           matrix box, rvec *x,
144                           gmx_ddbox_t *ddbox,
145                           int *npme_x, int *npme_y);
146
147 /*! \brief Initialize data structures for bonded interactions */
148 void dd_init_bondeds(FILE *fplog,
149                      gmx_domdec_t *dd, gmx_mtop_t *mtop,
150                      gmx_vsite_t *vsite,
151                      t_inputrec *ir, gmx_bool bBCheck, cginfo_mb_t *cginfo_mb);
152
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);
155
156 /*! \brief Set DD grid dimensions and limits
157  *
158  * Should be called after calling dd_init_bondeds.
159  */
160 void set_dd_parameters(FILE *fplog, gmx_domdec_t *dd, real dlb_scale,
161                        t_inputrec *ir,
162                        gmx_ddbox_t *ddbox);
163
164 /*! \brief Change the DD non-bonded communication cut-off.
165  *
166  * This could fail when trying to increase the cut-off,
167  * then FALSE will be returned and the cut-off is not modified.
168  */
169 gmx_bool change_dd_cutoff(t_commrec *cr, t_state *state, t_inputrec *ir,
170                           real cutoff_req );
171
172 /*! \brief Limit DLB to preserve the option of returning to the current cut-off.
173  *
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  * sets/changes the DLB limit such that using the passed (pair-list) cut-off
177  * should still be possible after subsequently setting a shorter cut-off
178  * with change_dd_cutoff.
179  */
180 void set_dd_dlb_max_cutoff(t_commrec *cr, real cutoff);
181
182 /*! \brief Return if we are currently using dynamic load balancing */
183 gmx_bool dd_dlb_is_on(const gmx_domdec_t *dd);
184
185 /*! \brief Return if the DLB lock is set */
186 gmx_bool dd_dlb_is_locked(const gmx_domdec_t *dd);
187
188 /*! \brief Set a lock such that with DLB=auto DLB cannot get turned on */
189 void dd_dlb_lock(gmx_domdec_t *dd);
190
191 /*! \brief Clear a lock such that with DLB=auto DLB may get turned on later */
192 void dd_dlb_unlock(gmx_domdec_t *dd);
193
194 /*! \brief Set up communication for averaging GPU wait times over ranks
195  *
196  * When domains (PP MPI ranks) share a GPU, the individual GPU wait times
197  * are meaningless, as it depends on the order in which tasks on the same
198  * GPU finish. Therefore there wait times need to be averaged over the ranks
199  * sharing the same GPU. This function sets up the communication for that.
200  */
201 void dd_setup_dlb_resource_sharing(t_commrec           *cr,
202                                    const gmx_hw_info_t *hwinfo,
203                                    const gmx_hw_opt_t  *hw_opt);
204
205 /*! \brief Sets up the DD communication setup */
206 void setup_dd_grid(FILE *fplog, gmx_domdec_t *dd);
207
208 /*! \brief Collects local rvec arrays \p lv to \p v on the master rank */
209 void dd_collect_vec(gmx_domdec_t *dd,
210                     t_state *state_local, rvec *lv, rvec *v);
211
212 /*! \brief Collects the local state \p state_local to \p state on the master rank */
213 void dd_collect_state(gmx_domdec_t *dd,
214                       t_state *state_local, t_state *state);
215
216 /*! \brief Cycle counter indices used internally in the domain decomposition */
217 enum {
218     ddCyclStep, ddCyclPPduringPME, ddCyclF, ddCyclWaitGPU, ddCyclPME, ddCyclNr
219 };
220
221 /*! \brief Add the wallcycle count to the DD counter */
222 void dd_cycles_add(gmx_domdec_t *dd, float cycles, int ddCycl);
223
224 /*! \brief Start the force flop count */
225 void dd_force_flop_start(gmx_domdec_t *dd, t_nrnb *nrnb);
226
227 /*! \brief Stop the force flop count */
228 void dd_force_flop_stop(gmx_domdec_t *dd, t_nrnb *nrnb);
229
230 /*! \brief Return the PME/PP force load ratio, or -1 if nothing was measured.
231  *
232  * Should only be called on the DD master node.
233  */
234 float dd_pme_f_ratio(gmx_domdec_t *dd);
235
236 /*! \brief Communicate the coordinates to the neighboring cells and do pbc. */
237 void dd_move_x(gmx_domdec_t *dd, matrix box, rvec x[]);
238
239 /*! \brief Sum the forces over the neighboring cells.
240  *
241  * When fshift!=NULL the shift forces are updated to obtain
242  * the correct virial from the single sum including f.
243  */
244 void dd_move_f(gmx_domdec_t *dd, rvec f[], rvec *fshift);
245
246 /*! \brief Communicate a real for each atom to the neighboring cells. */
247 void dd_atom_spread_real(gmx_domdec_t *dd, real v[]);
248
249 /*! \brief Sum the contributions to a real for each atom over the neighboring cells. */
250 void dd_atom_sum_real(gmx_domdec_t *dd, real v[]);
251
252 /*! \brief Partition the system over the nodes.
253  *
254  * step is only used for printing error messages.
255  * If bMasterState==TRUE then state_global from the master node is used,
256  * else state_local is redistributed between the nodes.
257  * When f!=NULL, *f will be reallocated to the size of state_local.
258  */
259 void dd_partition_system(FILE                *fplog,
260                          gmx_int64_t          step,
261                          t_commrec           *cr,
262                          gmx_bool             bMasterState,
263                          int                  nstglobalcomm,
264                          t_state             *state_global,
265                          gmx_mtop_t          *top_global,
266                          t_inputrec          *ir,
267                          t_state             *state_local,
268                          rvec               **f,
269                          t_mdatoms           *mdatoms,
270                          gmx_localtop_t      *top_local,
271                          t_forcerec          *fr,
272                          gmx_vsite_t         *vsite,
273                          gmx_shellfc_t        shellfc,
274                          gmx_constr_t         constr,
275                          t_nrnb              *nrnb,
276                          gmx_wallcycle_t      wcycle,
277                          gmx_bool             bVerbose);
278
279 /*! \brief Reset all the statistics and counters for total run counting */
280 void reset_dd_statistics_counters(gmx_domdec_t *dd);
281
282 /*! \brief Print statistics for domain decomposition communication */
283 void print_dd_statistics(t_commrec *cr, t_inputrec *ir, FILE *fplog);
284
285 /* In domdec_con.c */
286
287 /*! \brief Communicates the virtual site forces, reduces the shift forces when \p fshift != NULL */
288 void dd_move_f_vsites(gmx_domdec_t *dd, rvec *f, rvec *fshift);
289
290 /*! \brief Clears the forces for virtual sites */
291 void dd_clear_f_vsites(gmx_domdec_t *dd, rvec *f);
292
293 /*! \brief Move x0 and also x1 if x1!=NULL. bX1IsCoord tells if to do PBC on x1 */
294 void dd_move_x_constraints(gmx_domdec_t *dd, matrix box,
295                            rvec *x0, rvec *x1, gmx_bool bX1IsCoord);
296
297 /*! \brief Communicates the coordinates involved in virtual sites */
298 void dd_move_x_vsites(gmx_domdec_t *dd, matrix box, rvec *x);
299
300 /*! \brief Returns the local atom count array for all constraints
301  *
302  * The local atom count for a constraint, possible values 2/1/0, is needed
303  * to avoid not/double-counting contributions linked to the Lagrange
304  * multiplier, such as the virial and free-energy derivatives.
305  */
306 int *dd_constraints_nlocalatoms(gmx_domdec_t *dd);
307
308 /* In domdec_top.c */
309
310 /*! \brief Print error output when interactions are missing */
311 void dd_print_missing_interactions(FILE *fplog, t_commrec *cr,
312                                    int local_count,  gmx_mtop_t *top_global, t_state *state_local);
313
314 /*! \brief Generate and store the reverse topology */
315 void dd_make_reverse_top(FILE *fplog,
316                          gmx_domdec_t *dd, gmx_mtop_t *mtop,
317                          gmx_vsite_t *vsite,
318                          t_inputrec *ir, gmx_bool bBCheck);
319
320 /*! \brief Store the local charge group index in \p lcgs */
321 void dd_make_local_cgs(gmx_domdec_t *dd, t_block *lcgs);
322
323 /*! \brief Generate the local topology and virtual site data */
324 void dd_make_local_top(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
325                        int npbcdim, matrix box,
326                        rvec cellsize_min, ivec npulse,
327                        t_forcerec *fr,
328                        rvec *cgcm_or_x,
329                        gmx_vsite_t *vsite,
330                        gmx_mtop_t *top, gmx_localtop_t *ltop);
331
332 /*! \brief Sort ltop->ilist when we are doing free energy. */
333 void dd_sort_local_top(gmx_domdec_t *dd, t_mdatoms *mdatoms,
334                        gmx_localtop_t *ltop);
335
336 /*! \brief Construct local topology */
337 gmx_localtop_t *dd_init_local_top(gmx_mtop_t *top_global);
338
339 /*! \brief Construct local state */
340 void dd_init_local_state(gmx_domdec_t *dd,
341                          t_state *state_global, t_state *local_state);
342
343 /*! \brief Generate a list of links between charge groups that are linked by bonded interactions */
344 t_blocka *make_charge_group_links(gmx_mtop_t *mtop, gmx_domdec_t *dd,
345                                   cginfo_mb_t *cginfo_mb);
346
347 /*! \brief Calculate the maximum distance involved in 2-body and multi-body bonded interactions */
348 void dd_bonded_cg_distance(FILE *fplog, gmx_mtop_t *mtop,
349                            t_inputrec *ir, rvec *x, matrix box,
350                            gmx_bool bBCheck,
351                            real *r_2b, real *r_mb);
352
353 /*! \brief Dump a pdb file with the current DD home + communicated atoms.
354  *
355  * When natoms=-1, dump all known atoms.
356  */
357 void write_dd_pdb(const char *fn, gmx_int64_t step, const char *title,
358                   gmx_mtop_t *mtop,
359                   t_commrec *cr,
360                   int natoms, rvec x[], matrix box);
361
362
363 /* In domdec_setup.c */
364
365 /*! \brief Returns the volume fraction of the system that is communicated */
366 real comm_box_frac(ivec dd_nc, real cutoff, gmx_ddbox_t *ddbox);
367
368 /*! \brief Determines the optimal DD cell setup dd->nc and possibly npmenodes
369  * for the system.
370  *
371  * On the master node returns the actual cellsize limit used.
372  */
373 real dd_choose_grid(FILE *fplog,
374                     t_commrec *cr, gmx_domdec_t *dd, t_inputrec *ir,
375                     gmx_mtop_t *mtop, matrix box, gmx_ddbox_t *ddbox,
376                     gmx_bool bDynLoadBal, real dlb_scale,
377                     real cellsize_limit, real cutoff_dd,
378                     gmx_bool bInterCGBondeds);
379
380
381 /* In domdec_box.c */
382
383 /*! \brief Set the box and PBC data in \p ddbox */
384 void set_ddbox(gmx_domdec_t *dd, gmx_bool bMasterState, t_commrec *cr_sum,
385                t_inputrec *ir, matrix box,
386                gmx_bool bCalcUnboundedSize, t_block *cgs, rvec *x,
387                gmx_ddbox_t *ddbox);
388
389 /*! \brief Set the box and PBC data in \p ddbox */
390 void set_ddbox_cr(t_commrec *cr, ivec *dd_nc,
391                   t_inputrec *ir, matrix box, t_block *cgs, rvec *x,
392                   gmx_ddbox_t *ddbox);
393
394 #ifdef __cplusplus
395 }
396 #endif
397
398 #endif