Used IWYU to partially clean up some includes
[alexxy/gromacs.git] / src / gromacs / legacyheaders / 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, 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
36 #ifndef _domdec_h
37 #define _domdec_h
38
39 #include <stdio.h>
40
41 #include "gromacs/legacyheaders/vsite.h"
42 #include "gromacs/legacyheaders/types/commrec_fwd.h"
43 #include "gromacs/legacyheaders/types/constr.h"
44 #include "gromacs/legacyheaders/types/forcerec.h"
45 #include "gromacs/legacyheaders/types/hw_info.h"
46 #include "gromacs/legacyheaders/types/inputrec.h"
47 #include "gromacs/legacyheaders/types/mdatom.h"
48 #include "gromacs/legacyheaders/types/nrnb.h"
49 #include "gromacs/legacyheaders/types/shellfc.h"
50 #include "gromacs/legacyheaders/types/state.h"
51 #include "gromacs/math/vectypes.h"
52 #include "gromacs/timing/wallcycle.h"
53 #include "gromacs/topology/block.h"
54 #include "gromacs/topology/idef.h"
55 #include "gromacs/topology/topology.h"
56 #include "gromacs/utility/basedefinitions.h"
57 #include "gromacs/utility/real.h"
58
59 #ifdef __cplusplus
60 extern "C" {
61 #endif
62
63 int ddglatnr(gmx_domdec_t *dd, int i);
64 /* Returns the global topology atom number belonging to local atom index i.
65  * This function is intended for writing ascii output
66  * and returns atom numbers starting at 1.
67  * When dd=NULL returns i+1.
68  */
69
70 t_block *dd_charge_groups_global(gmx_domdec_t *dd);
71 /* Return a block struct for the charge groups of the whole system */
72
73 gmx_bool dd_filled_nsgrid_home(gmx_domdec_t *dd);
74 /* Is the ns grid already filled with the home particles? */
75
76 void dd_store_state(gmx_domdec_t *dd, t_state *state);
77 /* Store the global cg indices of the home cgs in state,
78  * so it can be reset, even after a new DD partitioning.
79  */
80
81 gmx_domdec_zones_t *domdec_zones(gmx_domdec_t *dd);
82
83 void dd_get_ns_ranges(gmx_domdec_t *dd, int icg,
84                       int *jcg0, int *jcg1, ivec shift0, ivec shift1);
85
86 int dd_natoms_vsite(gmx_domdec_t *dd);
87
88 void dd_get_constraint_range(gmx_domdec_t *dd,
89                              int *at_start, int *at_end);
90
91 real dd_cutoff_mbody(gmx_domdec_t *dd);
92
93 real dd_cutoff_twobody(gmx_domdec_t *dd);
94
95 void get_pme_nnodes(const gmx_domdec_t *dd,
96                     int *npmenodes_x, int *npmenodes_y);
97 /* Get the number of PME nodes along x and y, can be called with dd=NULL */
98
99 gmx_bool gmx_pmeonlynode(t_commrec *cr, int nodeid);
100 /* Return if nodeid in cr->mpi_comm_mysim is a PME-only node */
101
102 void get_pme_ddnodes(t_commrec *cr, int pmenodeid,
103                      int *nmy_ddnodes, int **my_ddnodes, int *node_peer);
104 /* Returns the set of DD nodes that communicate with pme node cr->nodeid */
105
106 int dd_pme_maxshift_x(gmx_domdec_t *dd);
107 /* Returns the maximum shift for coordinate communication in PME, dim x */
108
109 int dd_pme_maxshift_y(gmx_domdec_t *dd);
110 /* Returns the maximum shift for coordinate communication in PME, dim y */
111
112 void make_dd_communicators(FILE *fplog, t_commrec *cr, int dd_node_order);
113
114 gmx_domdec_t *
115 init_domain_decomposition(FILE *fplog,
116                           t_commrec *cr,
117                           unsigned long Flags,
118                           ivec nc,
119                           real comm_distance_min, real rconstr,
120                           const char *dlb_opt, real dlb_scale,
121                           const char *sizex, const char *sizey, const char *sizez,
122                           gmx_mtop_t *mtop, t_inputrec *ir,
123                           matrix box, rvec *x,
124                           gmx_ddbox_t *ddbox,
125                           int *npme_x, int *npme_y);
126
127 void dd_init_bondeds(FILE *fplog,
128                      gmx_domdec_t *dd, gmx_mtop_t *mtop,
129                      gmx_vsite_t *vsite,
130                      t_inputrec *ir, gmx_bool bBCheck, cginfo_mb_t *cginfo_mb);
131 /* Initialize data structures for bonded interactions */
132
133 gmx_bool dd_bonded_molpbc(gmx_domdec_t *dd, int ePBC);
134 /* Returns if we need to do pbc for calculating bonded interactions */
135
136 void set_dd_parameters(FILE *fplog, gmx_domdec_t *dd, real dlb_scale,
137                        t_inputrec *ir,
138                        gmx_ddbox_t *ddbox);
139 /* Set DD grid dimensions and limits,
140  * should be called after calling dd_init_bondeds.
141  */
142
143 gmx_bool change_dd_cutoff(t_commrec *cr, t_state *state, t_inputrec *ir,
144                           real cutoff_req );
145 /* Change the DD non-bonded communication cut-off.
146  * This could fail when trying to increase the cut-off,
147  * then FALSE will be returned and the cut-off is not modified.
148  */
149
150 void change_dd_dlb_cutoff_limit(t_commrec *cr);
151 /* Domain boundary changes due to the DD dynamic load balancing can limit
152  * the cut-off distance that can be set in change_dd_cutoff. This function
153  * limits the DLB such that using the currently set cut-off should still be
154  * possible after subsequently setting a shorter cut-off with change_dd_cutoff.
155  */
156
157 gmx_bool dd_dlb_is_locked(const gmx_domdec_t *dd);
158 /* Return if the DLB lock is set */
159
160 void dd_dlb_set_lock(gmx_domdec_t *dd, gmx_bool bValue);
161 /* Set a lock such that with DLB=auto DLB can (not) get turned on */
162
163 void dd_setup_dlb_resource_sharing(t_commrec           *cr,
164                                    const gmx_hw_info_t *hwinfo,
165                                    const gmx_hw_opt_t  *hw_opt);
166 /* When domains (PP MPI ranks) share a GPU, the individual GPU wait times
167  * are meaningless, as it depends on the order in which tasks on the same
168  * GPU finish. Therefore there wait times need to be averaged over the ranks
169  * sharing the same GPU. This function sets up the communication for that.
170  */
171
172 void setup_dd_grid(FILE *fplog, gmx_domdec_t *dd);
173
174 void dd_collect_vec(gmx_domdec_t *dd,
175                     t_state *state_local, rvec *lv, rvec *v);
176
177 void dd_collect_state(gmx_domdec_t *dd,
178                       t_state *state_local, t_state *state);
179
180 enum {
181     ddCyclStep, ddCyclPPduringPME, ddCyclF, ddCyclWaitGPU, ddCyclPME, ddCyclNr
182 };
183
184 void dd_cycles_add(gmx_domdec_t *dd, float cycles, int ddCycl);
185 /* Add the wallcycle count to the DD counter */
186
187 void dd_force_flop_start(gmx_domdec_t *dd, t_nrnb *nrnb);
188 /* Start the force flop count */
189
190 void dd_force_flop_stop(gmx_domdec_t *dd, t_nrnb *nrnb);
191 /* Stop the force flop count */
192
193 float dd_pme_f_ratio(gmx_domdec_t *dd);
194 /* Return the PME/PP force load ratio, or -1 if nothing was measured.
195  * Should only be called on the DD master node.
196  */
197
198 void dd_move_x(gmx_domdec_t *dd, matrix box, rvec x[]);
199 /* Communicate the coordinates to the neighboring cells and do pbc. */
200
201 void dd_move_f(gmx_domdec_t *dd, rvec f[], rvec *fshift);
202 /* Sum the forces over the neighboring cells.
203  * When fshift!=NULL the shift forces are updated to obtain
204  * the correct virial from the single sum including f.
205  */
206
207 void dd_atom_spread_real(gmx_domdec_t *dd, real v[]);
208 /* Communicate a real for each atom to the neighboring cells. */
209
210 void dd_atom_sum_real(gmx_domdec_t *dd, real v[]);
211 /* Sum the contributions to a real for each atom over the neighboring cells. */
212
213 void dd_partition_system(FILE                *fplog,
214                          gmx_int64_t          step,
215                          t_commrec           *cr,
216                          gmx_bool             bMasterState,
217                          int                  nstglobalcomm,
218                          t_state             *state_global,
219                          gmx_mtop_t          *top_global,
220                          t_inputrec          *ir,
221                          t_state             *state_local,
222                          rvec               **f,
223                          t_mdatoms           *mdatoms,
224                          gmx_localtop_t      *top_local,
225                          t_forcerec          *fr,
226                          gmx_vsite_t         *vsite,
227                          gmx_shellfc_t        shellfc,
228                          gmx_constr_t         constr,
229                          t_nrnb              *nrnb,
230                          gmx_wallcycle_t      wcycle,
231                          gmx_bool             bVerbose);
232 /* Partition the system over the nodes.
233  * step is only used for printing error messages.
234  * If bMasterState==TRUE then state_global from the master node is used,
235  * else state_local is redistributed between the nodes.
236  * When f!=NULL, *f will be reallocated to the size of state_local.
237  */
238
239 void reset_dd_statistics_counters(gmx_domdec_t *dd);
240 /* Reset all the statistics and counters for total run counting */
241
242 void print_dd_statistics(t_commrec *cr, t_inputrec *ir, FILE *fplog);
243
244 /* In domdec_con.c */
245
246 void dd_move_f_vsites(gmx_domdec_t *dd, rvec *f, rvec *fshift);
247
248 void dd_clear_f_vsites(gmx_domdec_t *dd, rvec *f);
249
250 void dd_move_x_constraints(gmx_domdec_t *dd, matrix box,
251                            rvec *x0, rvec *x1, gmx_bool bX1IsCoord);
252 /* Move x0 and also x1 if x1!=NULL. bX1IsCoord tells if to do PBC on x1 */
253
254 void dd_move_x_vsites(gmx_domdec_t *dd, matrix box, rvec *x);
255
256 int *dd_constraints_nlocalatoms(gmx_domdec_t *dd);
257
258 void dd_clear_local_constraint_indices(gmx_domdec_t *dd);
259
260 void dd_clear_local_vsite_indices(gmx_domdec_t *dd);
261
262 int dd_make_local_vsites(gmx_domdec_t *dd, int at_start, t_ilist *lil);
263
264 int dd_make_local_constraints(gmx_domdec_t *dd, int at_start,
265                               const gmx_mtop_t *mtop,
266                               const int *cginfo,
267                               gmx_constr_t constr, int nrec,
268                               t_ilist *il_local);
269
270 void init_domdec_constraints(gmx_domdec_t *dd,
271                              gmx_mtop_t   *mtop);
272
273 void init_domdec_vsites(gmx_domdec_t *dd, int n_intercg_vsite);
274
275
276 /* In domdec_top.c */
277
278 void dd_print_missing_interactions(FILE *fplog, t_commrec *cr,
279                                    int local_count,  gmx_mtop_t *top_global, t_state *state_local);
280
281 void dd_make_reverse_top(FILE *fplog,
282                          gmx_domdec_t *dd, gmx_mtop_t *mtop,
283                          gmx_vsite_t *vsite,
284                          t_inputrec *ir, gmx_bool bBCheck);
285
286 void dd_make_local_cgs(gmx_domdec_t *dd, t_block *lcgs);
287
288 void dd_make_local_top(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
289                        int npbcdim, matrix box,
290                        rvec cellsize_min, ivec npulse,
291                        t_forcerec *fr,
292                        rvec *cgcm_or_x,
293                        gmx_vsite_t *vsite,
294                        gmx_mtop_t *top, gmx_localtop_t *ltop);
295
296 void dd_sort_local_top(gmx_domdec_t *dd, t_mdatoms *mdatoms,
297                        gmx_localtop_t *ltop);
298 /* Sort ltop->ilist when we are doing free energy. */
299
300 gmx_localtop_t *dd_init_local_top(gmx_mtop_t *top_global);
301
302 void dd_init_local_state(gmx_domdec_t *dd,
303                          t_state *state_global, t_state *local_state);
304
305 t_blocka *make_charge_group_links(gmx_mtop_t *mtop, gmx_domdec_t *dd,
306                                   cginfo_mb_t *cginfo_mb);
307
308 void dd_bonded_cg_distance(FILE *fplog, gmx_mtop_t *mtop,
309                            t_inputrec *ir, rvec *x, matrix box,
310                            gmx_bool bBCheck,
311                            real *r_2b, real *r_mb);
312
313 void write_dd_pdb(const char *fn, gmx_int64_t step, const char *title,
314                   gmx_mtop_t *mtop,
315                   t_commrec *cr,
316                   int natoms, rvec x[], matrix box);
317 /* Dump a pdb file with the current DD home + communicated atoms.
318  * When natoms=-1, dump all known atoms.
319  */
320
321
322 /* In domdec_setup.c */
323
324 real comm_box_frac(ivec dd_nc, real cutoff, gmx_ddbox_t *ddbox);
325 /* Returns the volume fraction of the system that is communicated */
326
327 real dd_choose_grid(FILE *fplog,
328                     t_commrec *cr, gmx_domdec_t *dd, t_inputrec *ir,
329                     gmx_mtop_t *mtop, matrix box, gmx_ddbox_t *ddbox,
330                     gmx_bool bDynLoadBal, real dlb_scale,
331                     real cellsize_limit, real cutoff_dd,
332                     gmx_bool bInterCGBondeds);
333 /* Determines the optimal DD cell setup dd->nc and possibly npmenodes
334  * for the system.
335  * On the master node returns the actual cellsize limit used.
336  */
337
338
339 /* In domdec_box.c */
340
341 void set_ddbox(gmx_domdec_t *dd, gmx_bool bMasterState, t_commrec *cr_sum,
342                t_inputrec *ir, matrix box,
343                gmx_bool bCalcUnboundedSize, t_block *cgs, rvec *x,
344                gmx_ddbox_t *ddbox);
345
346 void set_ddbox_cr(t_commrec *cr, ivec *dd_nc,
347                   t_inputrec *ir, matrix box, t_block *cgs, rvec *x,
348                   gmx_ddbox_t *ddbox);
349
350 #ifdef __cplusplus
351 }
352 #endif
353
354 #endif  /* _domdec_h */