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 void get_pme_nnodes(const gmx_domdec_t *dd,
64 int *npmenodes_x, int *npmenodes_y);
65 /* Get the number of PME nodes along x and y, can be called with dd=NULL */
67 gmx_bool gmx_pmeonlynode(t_commrec *cr, int nodeid);
68 /* Return if nodeid in cr->mpi_comm_mysim is a PME-only node */
70 void get_pme_ddnodes(t_commrec *cr, int pmenodeid,
71 int *nmy_ddnodes, int **my_ddnodes, int *node_peer);
72 /* Returns the set of DD nodes that communicate with pme node cr->nodeid */
74 int dd_pme_maxshift_x(gmx_domdec_t *dd);
75 /* Returns the maximum shift for coordinate communication in PME, dim x */
77 int dd_pme_maxshift_y(gmx_domdec_t *dd);
78 /* Returns the maximum shift for coordinate communication in PME, dim y */
80 void make_dd_communicators(FILE *fplog, t_commrec *cr, int dd_node_order);
83 init_domain_decomposition(FILE *fplog,
87 real comm_distance_min, real rconstr,
88 const char *dlb_opt, real dlb_scale,
89 const char *sizex, const char *sizey, const char *sizez,
90 gmx_mtop_t *mtop, t_inputrec *ir,
93 int *npme_x, int *npme_y);
95 void dd_init_bondeds(FILE *fplog,
96 gmx_domdec_t *dd, gmx_mtop_t *mtop,
98 t_inputrec *ir, gmx_bool bBCheck, cginfo_mb_t *cginfo_mb);
99 /* Initialize data structures for bonded interactions */
101 gmx_bool dd_bonded_molpbc(gmx_domdec_t *dd, int ePBC);
102 /* Returns if we need to do pbc for calculating bonded interactions */
104 void set_dd_parameters(FILE *fplog, gmx_domdec_t *dd, real dlb_scale,
107 /* Set DD grid dimensions and limits,
108 * should be called after calling dd_init_bondeds.
111 gmx_bool change_dd_cutoff(t_commrec *cr, t_state *state, t_inputrec *ir,
113 /* Change the DD non-bonded communication cut-off.
114 * This could fail when trying to increase the cut-off,
115 * then FALSE will be returned and the cut-off is not modified.
118 void change_dd_dlb_cutoff_limit(t_commrec *cr);
119 /* Domain boundary changes due to the DD dynamic load balancing can limit
120 * the cut-off distance that can be set in change_dd_cutoff. This function
121 * limits the DLB such that using the currently set cut-off should still be
122 * possible after subsequently setting a shorter cut-off with change_dd_cutoff.
125 void dd_setup_dlb_resource_sharing(t_commrec *cr,
126 const gmx_hw_info_t *hwinfo,
127 const gmx_hw_opt_t *hw_opt);
128 /* When domains (PP MPI ranks) share a GPU, the individual GPU wait times
129 * are meaningless, as it depends on the order in which tasks on the same
130 * GPU finish. Therefore there wait times need to be averaged over the ranks
131 * sharing the same GPU. This function sets up the communication for that.
134 void setup_dd_grid(FILE *fplog, gmx_domdec_t *dd);
136 void dd_collect_vec(gmx_domdec_t *dd,
137 t_state *state_local, rvec *lv, rvec *v);
139 void dd_collect_state(gmx_domdec_t *dd,
140 t_state *state_local, t_state *state);
143 ddCyclStep, ddCyclPPduringPME, ddCyclF, ddCyclWaitGPU, ddCyclPME, ddCyclNr
146 void dd_cycles_add(gmx_domdec_t *dd, float cycles, int ddCycl);
147 /* Add the wallcycle count to the DD counter */
149 void dd_force_flop_start(gmx_domdec_t *dd, t_nrnb *nrnb);
150 /* Start the force flop count */
152 void dd_force_flop_stop(gmx_domdec_t *dd, t_nrnb *nrnb);
153 /* Stop the force flop count */
155 float dd_pme_f_ratio(gmx_domdec_t *dd);
156 /* Return the PME/PP force load ratio, or -1 if nothing was measured.
157 * Should only be called on the DD master node.
160 void dd_move_x(gmx_domdec_t *dd, matrix box, rvec x[]);
161 /* Communicate the coordinates to the neighboring cells and do pbc. */
163 void dd_move_f(gmx_domdec_t *dd, rvec f[], rvec *fshift);
164 /* Sum the forces over the neighboring cells.
165 * When fshift!=NULL the shift forces are updated to obtain
166 * the correct virial from the single sum including f.
169 void dd_atom_spread_real(gmx_domdec_t *dd, real v[]);
170 /* Communicate a real for each atom to the neighboring cells. */
172 void dd_atom_sum_real(gmx_domdec_t *dd, real v[]);
173 /* Sum the contributions to a real for each atom over the neighboring cells. */
175 void dd_partition_system(FILE *fplog,
176 gmx_large_int_t step,
178 gmx_bool bMasterState,
180 t_state *state_global,
181 gmx_mtop_t *top_global,
183 t_state *state_local,
186 gmx_localtop_t *top_local,
189 gmx_shellfc_t shellfc,
192 gmx_wallcycle_t wcycle,
194 /* Partition the system over the nodes.
195 * step is only used for printing error messages.
196 * If bMasterState==TRUE then state_global from the master node is used,
197 * else state_local is redistributed between the nodes.
198 * When f!=NULL, *f will be reallocated to the size of state_local.
201 void reset_dd_statistics_counters(gmx_domdec_t *dd);
202 /* Reset all the statistics and counters for total run counting */
204 void print_dd_statistics(t_commrec *cr, t_inputrec *ir, FILE *fplog);
206 /* In domdec_con.c */
208 void dd_move_f_vsites(gmx_domdec_t *dd, rvec *f, rvec *fshift);
210 void dd_clear_f_vsites(gmx_domdec_t *dd, rvec *f);
212 void dd_move_x_constraints(gmx_domdec_t *dd, matrix box,
214 /* Move x0 and also x1 if x1!=NULL */
216 void dd_move_x_vsites(gmx_domdec_t *dd, matrix box, rvec *x);
218 int *dd_constraints_nlocalatoms(gmx_domdec_t *dd);
220 void dd_clear_local_constraint_indices(gmx_domdec_t *dd);
222 void dd_clear_local_vsite_indices(gmx_domdec_t *dd);
224 int dd_make_local_vsites(gmx_domdec_t *dd, int at_start, t_ilist *lil);
226 int dd_make_local_constraints(gmx_domdec_t *dd, int at_start,
227 const gmx_mtop_t *mtop,
229 gmx_constr_t constr, int nrec,
232 void init_domdec_constraints(gmx_domdec_t *dd,
235 void init_domdec_vsites(gmx_domdec_t *dd, int n_intercg_vsite);
238 /* In domdec_top.c */
240 void dd_print_missing_interactions(FILE *fplog, t_commrec *cr,
241 int local_count, gmx_mtop_t *top_global, t_state *state_local);
243 void dd_make_reverse_top(FILE *fplog,
244 gmx_domdec_t *dd, gmx_mtop_t *mtop,
246 t_inputrec *ir, gmx_bool bBCheck);
248 void dd_make_local_cgs(gmx_domdec_t *dd, t_block *lcgs);
250 void dd_make_local_top(gmx_domdec_t *dd, gmx_domdec_zones_t *zones,
251 int npbcdim, matrix box,
252 rvec cellsize_min, ivec npulse,
256 gmx_mtop_t *top, gmx_localtop_t *ltop);
258 void dd_sort_local_top(gmx_domdec_t *dd, t_mdatoms *mdatoms,
259 gmx_localtop_t *ltop);
260 /* Sort ltop->ilist when we are doing free energy. */
262 gmx_localtop_t *dd_init_local_top(gmx_mtop_t *top_global);
264 void dd_init_local_state(gmx_domdec_t *dd,
265 t_state *state_global, t_state *local_state);
267 t_blocka *make_charge_group_links(gmx_mtop_t *mtop, gmx_domdec_t *dd,
268 cginfo_mb_t *cginfo_mb);
270 void dd_bonded_cg_distance(FILE *fplog, gmx_mtop_t *mtop,
271 t_inputrec *ir, rvec *x, matrix box,
273 real *r_2b, real *r_mb);
275 void write_dd_pdb(const char *fn, gmx_large_int_t step, const char *title,
278 int natoms, rvec x[], matrix box);
279 /* Dump a pdb file with the current DD home + communicated atoms.
280 * When natoms=-1, dump all known atoms.
284 /* In domdec_setup.c */
286 real comm_box_frac(ivec dd_nc, real cutoff, gmx_ddbox_t *ddbox);
287 /* Returns the volume fraction of the system that is communicated */
289 real dd_choose_grid(FILE *fplog,
290 t_commrec *cr, gmx_domdec_t *dd, t_inputrec *ir,
291 gmx_mtop_t *mtop, matrix box, gmx_ddbox_t *ddbox,
292 gmx_bool bDynLoadBal, real dlb_scale,
293 real cellsize_limit, real cutoff_dd,
294 gmx_bool bInterCGBondeds);
295 /* Determines the optimal DD cell setup dd->nc and possibly npmenodes
297 * On the master node returns the actual cellsize limit used.
301 /* In domdec_box.c */
303 void set_ddbox(gmx_domdec_t *dd, gmx_bool bMasterState, t_commrec *cr_sum,
304 t_inputrec *ir, matrix box,
305 gmx_bool bCalcUnboundedSize, t_block *cgs, rvec *x,
308 void set_ddbox_cr(t_commrec *cr, ivec *dd_nc,
309 t_inputrec *ir, matrix box, t_block *cgs, rvec *x,
316 #endif /* _domdec_h */