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