Update copyright statements and change license to LGPL
[alexxy/gromacs.git] / include / domdec.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * This file is part of Gromacs        Copyright (c) 1991-2008
5  * David van der Spoel, Erik Lindahl, Berk Hess, University of Groningen.
6  * Copyright (c) 2012, by the GROMACS development team, led by
7  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
8  * others, as listed in the AUTHORS file in the top-level source
9  * directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37
38 #ifndef _domdec_h
39 #define _domdec_h
40 #include "visibility.h"
41 #include "typedefs.h"
42 #include "types/commrec.h"
43 #include "vsite.h"
44 #include "genborn.h"
45
46 #ifdef __cplusplus
47 extern "C" {
48 #endif
49
50 int ddglatnr(gmx_domdec_t *dd,int i);
51 /* Returns the global topology atom number belonging to local atom index i.
52  * This function is intended for writing ascii output
53  * and returns atom numbers starting at 1.
54  * When dd=NULL returns i+1.
55  */
56
57 t_block *dd_charge_groups_global(gmx_domdec_t *dd);
58 /* Return a block struct for the charge groups of the whole system */
59
60 gmx_bool dd_filled_nsgrid_home(gmx_domdec_t *dd);
61 /* Is the ns grid already filled with the home particles? */
62
63 void dd_store_state(gmx_domdec_t *dd,t_state *state);
64 /* Store the global cg indices of the home cgs in state,
65  * so it can be reset, even after a new DD partitioning.
66  */
67
68 gmx_domdec_zones_t *domdec_zones(gmx_domdec_t *dd);
69
70 void dd_get_ns_ranges(gmx_domdec_t *dd,int icg,
71                              int *jcg0,int *jcg1,ivec shift0,ivec shift1);
72
73 int dd_natoms_vsite(gmx_domdec_t *dd);
74
75 void dd_get_constraint_range(gmx_domdec_t *dd,
76                                     int *at_start,int *at_end);
77
78 real dd_cutoff_mbody(gmx_domdec_t *dd);
79
80 real dd_cutoff_twobody(gmx_domdec_t *dd);
81
82 gmx_bool gmx_pmeonlynode(t_commrec *cr,int nodeid);
83 /* Return if nodeid in cr->mpi_comm_mysim is a PME-only node */
84
85 void get_pme_ddnodes(t_commrec *cr,int pmenodeid,
86                             int *nmy_ddnodes,int **my_ddnodes,int *node_peer);
87 /* Returns the set of DD nodes that communicate with pme node cr->nodeid */
88
89 int dd_pme_maxshift_x(gmx_domdec_t *dd);
90 /* Returns the maximum shift for coordinate communication in PME, dim x */
91
92 int dd_pme_maxshift_y(gmx_domdec_t *dd);
93 /* Returns the maximum shift for coordinate communication in PME, dim y */
94
95 GMX_LIBMD_EXPORT
96 void make_dd_communicators(FILE *fplog,t_commrec *cr,int dd_node_order);
97
98 GMX_LIBMD_EXPORT
99 gmx_domdec_t *
100 init_domain_decomposition(FILE *fplog,
101                           t_commrec *cr,
102                           unsigned long Flags,
103                           ivec nc,
104                           real comm_distance_min,real rconstr,
105                           const char *dlb_opt,real dlb_scale,
106                           const char *sizex,const char *sizey,const char *sizez,
107                           gmx_mtop_t *mtop,t_inputrec *ir,
108                           matrix box,rvec *x,
109                           gmx_ddbox_t *ddbox,
110                           int *npme_x, int *npme_y);
111
112 GMX_LIBMD_EXPORT
113 void dd_init_bondeds(FILE *fplog,
114                             gmx_domdec_t *dd,gmx_mtop_t *mtop,
115                             gmx_vsite_t *vsite,gmx_constr_t constr,
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 GMX_LIBMD_EXPORT
123 void set_dd_parameters(FILE *fplog,gmx_domdec_t *dd,real dlb_scale,
124                               t_inputrec *ir,t_forcerec *fr,
125                               gmx_ddbox_t *ddbox);
126 /* Set DD grid dimensions and limits,
127  * should be called after calling dd_init_bondeds.
128  */
129
130 GMX_LIBMD_EXPORT
131 gmx_bool change_dd_cutoff(t_commrec *cr,t_state *state,t_inputrec *ir,
132                           real cutoff_req );
133 /* Change the DD non-bonded communication cut-off.
134  * This could fail when trying to increase the cut-off,
135  * then FALSE will be returned and the cut-off is not modified.
136  */
137
138 GMX_LIBMD_EXPORT
139 void setup_dd_grid(FILE *fplog,gmx_domdec_t *dd);
140
141 void dd_collect_vec(gmx_domdec_t *dd,
142                            t_state *state_local,rvec *lv,rvec *v);
143
144 GMX_LIBMD_EXPORT
145 void dd_collect_state(gmx_domdec_t *dd,
146                              t_state *state_local,t_state *state);
147
148 enum { ddCyclStep, ddCyclPPduringPME, ddCyclF, ddCyclPME, ddCyclNr };
149
150 GMX_LIBMD_EXPORT
151 void dd_cycles_add(gmx_domdec_t *dd,float cycles,int ddCycl);
152 /* Add the wallcycle count to the DD counter */
153
154 void dd_force_flop_start(gmx_domdec_t *dd,t_nrnb *nrnb);
155 /* Start the force flop count */
156
157 void dd_force_flop_stop(gmx_domdec_t *dd,t_nrnb *nrnb);
158 /* Stop the force flop count */
159
160 GMX_LIBMD_EXPORT
161 float dd_pme_f_ratio(gmx_domdec_t *dd);
162 /* Return the PME/PP force load ratio, or -1 if nothing was measured.
163  * Should only be called on the DD master node.
164  */
165
166 void dd_move_x(gmx_domdec_t *dd,matrix box,rvec x[]);
167 /* Communicate the coordinates to the neighboring cells and do pbc. */
168
169 void dd_move_f(gmx_domdec_t *dd,rvec f[],rvec *fshift);
170 /* Sum the forces over the neighboring cells.
171  * When fshift!=NULL the shift forces are updated to obtain
172  * the correct virial from the single sum including f.
173  */
174
175 void dd_atom_spread_real(gmx_domdec_t *dd,real v[]);
176 /* Communicate a real for each atom to the neighboring cells. */
177
178 void dd_atom_sum_real(gmx_domdec_t *dd,real v[]);
179 /* Sum the contributions to a real for each atom over the neighboring cells. */
180
181 GMX_LIBMD_EXPORT
182 void dd_partition_system(FILE            *fplog,
183                                 gmx_large_int_t      step,
184                                 t_commrec       *cr,
185                                 gmx_bool            bMasterState,
186                                 int             nstglobalcomm,
187                                 t_state         *state_global,
188                                 gmx_mtop_t      *top_global,
189                                 t_inputrec      *ir,
190                                 t_state         *state_local,
191                                 rvec            **f,
192                                 t_mdatoms       *mdatoms,
193                                 gmx_localtop_t  *top_local,
194                                 t_forcerec      *fr,
195                                 gmx_vsite_t     *vsite,
196                                 gmx_shellfc_t   shellfc,
197                                 gmx_constr_t    constr,
198                                 t_nrnb          *nrnb,
199                                 gmx_wallcycle_t wcycle,
200                                 gmx_bool            bVerbose);
201 /* Partition the system over the nodes.
202  * step is only used for printing error messages.
203  * If bMasterState==TRUE then state_global from the master node is used,
204  * else state_local is redistributed between the nodes.
205  * When f!=NULL, *f will be reallocated to the size of state_local.
206  */
207
208 GMX_LIBMD_EXPORT
209 void reset_dd_statistics_counters(gmx_domdec_t *dd);
210 /* Reset all the statistics and counters for total run counting */
211
212 void print_dd_statistics(t_commrec *cr,t_inputrec *ir,FILE *fplog);
213
214 /* In domdec_con.c */
215
216 void dd_move_f_vsites(gmx_domdec_t *dd,rvec *f,rvec *fshift);
217
218 void dd_clear_f_vsites(gmx_domdec_t *dd,rvec *f);
219
220 void dd_move_x_constraints(gmx_domdec_t *dd,matrix box,
221                                   rvec *x0,rvec *x1);
222 /* Move x0 and also x1 if x1!=NULL */
223
224 void dd_move_x_vsites(gmx_domdec_t *dd,matrix box,rvec *x);
225
226 int *dd_constraints_nlocalatoms(gmx_domdec_t *dd);
227
228 void dd_clear_local_constraint_indices(gmx_domdec_t *dd);
229
230 void dd_clear_local_vsite_indices(gmx_domdec_t *dd);
231
232 int dd_make_local_vsites(gmx_domdec_t *dd,int at_start,t_ilist *lil);
233
234 int dd_make_local_constraints(gmx_domdec_t *dd,int at_start,
235                               const gmx_mtop_t *mtop,
236                               const int *cginfo,
237                               gmx_constr_t constr,int nrec,
238                               t_ilist *il_local);
239
240 void init_domdec_constraints(gmx_domdec_t *dd,
241                              gmx_mtop_t *mtop,
242                              gmx_constr_t constr);
243
244 void init_domdec_vsites(gmx_domdec_t *dd,int n_intercg_vsite);
245
246
247 /* In domdec_top.c */
248
249 void dd_print_missing_interactions(FILE *fplog,t_commrec *cr,
250                                           int local_count,  gmx_mtop_t *top_global, t_state *state_local);
251
252 void dd_make_reverse_top(FILE *fplog,
253                                 gmx_domdec_t *dd,gmx_mtop_t *mtop,
254                                 gmx_vsite_t *vsite,gmx_constr_t constr,
255                                 t_inputrec *ir,gmx_bool bBCheck);
256
257 void dd_make_local_cgs(gmx_domdec_t *dd,t_block *lcgs);
258
259 void dd_make_local_top(FILE *fplog,
260                        gmx_domdec_t *dd,gmx_domdec_zones_t *zones,
261                        int npbcdim,matrix box,
262                        rvec cellsize_min,ivec npulse,
263                        t_forcerec *fr,
264                        rvec *cgcm_or_x,
265                        gmx_vsite_t *vsite,
266                        gmx_mtop_t *top,gmx_localtop_t *ltop);
267
268 void dd_sort_local_top(gmx_domdec_t *dd,t_mdatoms *mdatoms,
269                               gmx_localtop_t *ltop);
270 /* Sort ltop->ilist when we are doing free energy. */
271
272 GMX_LIBMD_EXPORT
273 gmx_localtop_t *dd_init_local_top(gmx_mtop_t *top_global);
274
275 GMX_LIBMD_EXPORT
276 void dd_init_local_state(gmx_domdec_t *dd,
277                                 t_state *state_global,t_state *local_state);
278
279 t_blocka *make_charge_group_links(gmx_mtop_t *mtop,gmx_domdec_t *dd,
280                                          cginfo_mb_t *cginfo_mb);
281
282 void dd_bonded_cg_distance(FILE *fplog,
283                                   gmx_domdec_t *dd,gmx_mtop_t *mtop,
284                                   t_inputrec *ir,rvec *x,matrix box,
285                                   gmx_bool bBCheck,
286                                   real *r_2b,real *r_mb);
287
288 void write_dd_pdb(const char *fn,gmx_large_int_t step,const char *title,
289                          gmx_mtop_t *mtop,
290                          t_commrec *cr,
291                          int natoms,rvec x[],matrix box);
292 /* Dump a pdb file with the current DD home + communicated atoms.
293  * When natoms=-1, dump all known atoms.
294  */
295
296
297 /* In domdec_setup.c */
298
299 real comm_box_frac(ivec dd_nc,real cutoff,gmx_ddbox_t *ddbox);
300 /* Returns the volume fraction of the system that is communicated */
301
302 real dd_choose_grid(FILE *fplog,
303                            t_commrec *cr,gmx_domdec_t *dd,t_inputrec *ir,
304                            gmx_mtop_t *mtop,matrix box,gmx_ddbox_t *ddbox,
305                            gmx_bool bDynLoadBal,real dlb_scale,
306                            real cellsize_limit,real cutoff_dd,
307                            gmx_bool bInterCGBondeds,gmx_bool bInterCGMultiBody);
308 /* Determines the optimal DD cell setup dd->nc and possibly npmenodes
309  * for the system.
310  * On the master node returns the actual cellsize limit used.
311  */
312
313
314 /* In domdec_box.c */
315
316 void set_ddbox(gmx_domdec_t *dd,gmx_bool bMasterState,t_commrec *cr_sum,
317                       t_inputrec *ir,matrix box,
318                       gmx_bool bCalcUnboundedSize,t_block *cgs,rvec *x,
319                       gmx_ddbox_t *ddbox);
320
321 void set_ddbox_cr(t_commrec *cr,ivec *dd_nc,
322                          t_inputrec *ir,matrix box,t_block *cgs,rvec *x,
323                          gmx_ddbox_t *ddbox);
324
325 #ifdef __cplusplus
326 }
327 #endif
328
329 #endif  /* _domdec_h */