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