1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This file is part of Gromacs Copyright (c) 1991-2008
5 * David van der Spoel, Erik Lindahl, Berk Hess, University of Groningen.
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.
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
16 * Gnomes, ROck Monsters And Chili Sauce
23 #include "types/commrec.h"
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.
38 t_block *dd_charge_groups_global(gmx_domdec_t *dd);
39 /* Return a block struct for the charge groups of the whole system */
41 gmx_bool dd_filled_nsgrid_home(gmx_domdec_t *dd);
42 /* Is the ns grid already filled with the home particles? */
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.
49 gmx_domdec_zones_t *domdec_zones(gmx_domdec_t *dd);
51 void dd_get_ns_ranges(gmx_domdec_t *dd,int icg,
52 int *jcg0,int *jcg1,ivec shift0,ivec shift1);
54 int dd_natoms_vsite(gmx_domdec_t *dd);
56 void dd_get_constraint_range(gmx_domdec_t *dd,
57 int *at_start,int *at_end);
59 real dd_cutoff_mbody(gmx_domdec_t *dd);
61 real dd_cutoff_twobody(gmx_domdec_t *dd);
63 gmx_bool gmx_pmeonlynode(t_commrec *cr,int nodeid);
64 /* Return if nodeid in cr->mpi_comm_mysim is a PME-only node */
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 */
70 int dd_pme_maxshift_x(gmx_domdec_t *dd);
71 /* Returns the maximum shift for coordinate communication in PME, dim x */
73 int dd_pme_maxshift_y(gmx_domdec_t *dd);
74 /* Returns the maximum shift for coordinate communication in PME, dim y */
76 void make_dd_communicators(FILE *fplog,t_commrec *cr,int dd_node_order);
79 init_domain_decomposition(FILE *fplog,
83 real comm_distance_min,real rconstr,
84 const char *dlb_opt,real dlb_scale,
85 const char *sizex,const char *sizey,const char *sizez,
86 gmx_mtop_t *mtop,t_inputrec *ir,
89 int *npme_x, int *npme_y);
91 void dd_init_bondeds(FILE *fplog,
92 gmx_domdec_t *dd,gmx_mtop_t *mtop,
93 gmx_vsite_t *vsite,gmx_constr_t constr,
94 t_inputrec *ir,gmx_bool bBCheck,cginfo_mb_t *cginfo_mb);
95 /* Initialize data structures for bonded interactions */
97 gmx_bool dd_bonded_molpbc(gmx_domdec_t *dd,int ePBC);
98 /* Returns if we need to do pbc for calculating bonded interactions */
100 void set_dd_parameters(FILE *fplog,gmx_domdec_t *dd,real dlb_scale,
101 t_inputrec *ir,t_forcerec *fr,
103 /* Set DD grid dimensions and limits,
104 * should be called after calling dd_init_bondeds.
107 gmx_bool change_dd_cutoff(t_commrec *cr,t_state *state,t_inputrec *ir,
109 /* Change the DD non-bonded communication cut-off.
110 * This could fail when trying to increase the cut-off,
111 * then FALSE will be returned and the cut-off is not modified.
114 void change_dd_dlb_cutoff_limit(t_commrec *cr);
115 /* Domain boundary changes due to the DD dynamic load balancing can limit
116 * the cut-off distance that can be set in change_dd_cutoff. This function
117 * limits the DLB such that using the currently set cut-off should still be
118 * possible after subsequently setting a shorter cut-off with change_dd_cutoff.
121 void setup_dd_grid(FILE *fplog,gmx_domdec_t *dd);
123 void dd_collect_vec(gmx_domdec_t *dd,
124 t_state *state_local,rvec *lv,rvec *v);
126 void dd_collect_state(gmx_domdec_t *dd,
127 t_state *state_local,t_state *state);
129 enum { ddCyclStep, ddCyclPPduringPME, ddCyclF, ddCyclPME, ddCyclNr };
131 void dd_cycles_add(gmx_domdec_t *dd,float cycles,int ddCycl);
132 /* Add the wallcycle count to the DD counter */
134 void dd_force_flop_start(gmx_domdec_t *dd,t_nrnb *nrnb);
135 /* Start the force flop count */
137 void dd_force_flop_stop(gmx_domdec_t *dd,t_nrnb *nrnb);
138 /* Stop the force flop count */
140 float dd_pme_f_ratio(gmx_domdec_t *dd);
141 /* Return the PME/PP force load ratio, or -1 if nothing was measured.
142 * Should only be called on the DD master node.
145 void dd_move_x(gmx_domdec_t *dd,matrix box,rvec x[]);
146 /* Communicate the coordinates to the neighboring cells and do pbc. */
148 void dd_move_f(gmx_domdec_t *dd,rvec f[],rvec *fshift);
149 /* Sum the forces over the neighboring cells.
150 * When fshift!=NULL the shift forces are updated to obtain
151 * the correct virial from the single sum including f.
154 void dd_atom_spread_real(gmx_domdec_t *dd,real v[]);
155 /* Communicate a real for each atom to the neighboring cells. */
157 void dd_atom_sum_real(gmx_domdec_t *dd,real v[]);
158 /* Sum the contributions to a real for each atom over the neighboring cells. */
160 void dd_partition_system(FILE *fplog,
161 gmx_large_int_t step,
163 gmx_bool bMasterState,
165 t_state *state_global,
166 gmx_mtop_t *top_global,
168 t_state *state_local,
171 gmx_localtop_t *top_local,
174 gmx_shellfc_t shellfc,
177 gmx_wallcycle_t wcycle,
179 /* Partition the system over the nodes.
180 * step is only used for printing error messages.
181 * If bMasterState==TRUE then state_global from the master node is used,
182 * else state_local is redistributed between the nodes.
183 * When f!=NULL, *f will be reallocated to the size of state_local.
186 void reset_dd_statistics_counters(gmx_domdec_t *dd);
187 /* Reset all the statistics and counters for total run counting */
189 void print_dd_statistics(t_commrec *cr,t_inputrec *ir,FILE *fplog);
191 /* In domdec_con.c */
193 void dd_move_f_vsites(gmx_domdec_t *dd,rvec *f,rvec *fshift);
195 void dd_clear_f_vsites(gmx_domdec_t *dd,rvec *f);
197 void dd_move_x_constraints(gmx_domdec_t *dd,matrix box,
199 /* Move x0 and also x1 if x1!=NULL */
201 void dd_move_x_vsites(gmx_domdec_t *dd,matrix box,rvec *x);
203 int *dd_constraints_nlocalatoms(gmx_domdec_t *dd);
205 void dd_clear_local_constraint_indices(gmx_domdec_t *dd);
207 void dd_clear_local_vsite_indices(gmx_domdec_t *dd);
209 int dd_make_local_vsites(gmx_domdec_t *dd,int at_start,t_ilist *lil);
211 int dd_make_local_constraints(gmx_domdec_t *dd,int at_start,
212 const gmx_mtop_t *mtop,
214 gmx_constr_t constr,int nrec,
217 void init_domdec_constraints(gmx_domdec_t *dd,
219 gmx_constr_t constr);
221 void init_domdec_vsites(gmx_domdec_t *dd,int n_intercg_vsite);
224 /* In domdec_top.c */
226 void dd_print_missing_interactions(FILE *fplog,t_commrec *cr,
227 int local_count, gmx_mtop_t *top_global, t_state *state_local);
229 void dd_make_reverse_top(FILE *fplog,
230 gmx_domdec_t *dd,gmx_mtop_t *mtop,
231 gmx_vsite_t *vsite,gmx_constr_t constr,
232 t_inputrec *ir,gmx_bool bBCheck);
234 void dd_make_local_cgs(gmx_domdec_t *dd,t_block *lcgs);
236 void dd_make_local_top(FILE *fplog,
237 gmx_domdec_t *dd,gmx_domdec_zones_t *zones,
238 int npbcdim,matrix box,
239 rvec cellsize_min,ivec npulse,
243 gmx_mtop_t *top,gmx_localtop_t *ltop);
245 void dd_sort_local_top(gmx_domdec_t *dd,t_mdatoms *mdatoms,
246 gmx_localtop_t *ltop);
247 /* Sort ltop->ilist when we are doing free energy. */
249 gmx_localtop_t *dd_init_local_top(gmx_mtop_t *top_global);
251 void dd_init_local_state(gmx_domdec_t *dd,
252 t_state *state_global,t_state *local_state);
254 t_blocka *make_charge_group_links(gmx_mtop_t *mtop,gmx_domdec_t *dd,
255 cginfo_mb_t *cginfo_mb);
257 void dd_bonded_cg_distance(FILE *fplog,
258 gmx_domdec_t *dd,gmx_mtop_t *mtop,
259 t_inputrec *ir,rvec *x,matrix box,
261 real *r_2b,real *r_mb);
263 void write_dd_pdb(const char *fn,gmx_large_int_t step,const char *title,
266 int natoms,rvec x[],matrix box);
267 /* Dump a pdb file with the current DD home + communicated atoms.
268 * When natoms=-1, dump all known atoms.
272 /* In domdec_setup.c */
274 real comm_box_frac(ivec dd_nc,real cutoff,gmx_ddbox_t *ddbox);
275 /* Returns the volume fraction of the system that is communicated */
277 real dd_choose_grid(FILE *fplog,
278 t_commrec *cr,gmx_domdec_t *dd,t_inputrec *ir,
279 gmx_mtop_t *mtop,matrix box,gmx_ddbox_t *ddbox,
280 gmx_bool bDynLoadBal,real dlb_scale,
281 real cellsize_limit,real cutoff_dd,
282 gmx_bool bInterCGBondeds,gmx_bool bInterCGMultiBody);
283 /* Determines the optimal DD cell setup dd->nc and possibly npmenodes
285 * On the master node returns the actual cellsize limit used.
289 /* In domdec_box.c */
291 void set_ddbox(gmx_domdec_t *dd,gmx_bool bMasterState,t_commrec *cr_sum,
292 t_inputrec *ir,matrix box,
293 gmx_bool bCalcUnboundedSize,t_block *cgs,rvec *x,
296 void set_ddbox_cr(t_commrec *cr,ivec *dd_nc,
297 t_inputrec *ir,matrix box,t_block *cgs,rvec *x,
304 #endif /* _domdec_h */