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
37 int ddglatnr(gmx_domdec_t
*dd
,int i
);
38 /* Returns the global topology atom number belonging to local atom index i.
39 * This function is intended for writing ascii output
40 * and returns atom numbers starting at 1.
41 * When dd=NULL returns i+1.
44 t_block
*dd_charge_groups_global(gmx_domdec_t
*dd
);
45 /* Return a block struct for the charge groups of the whole system */
47 gmx_bool
dd_filled_nsgrid_home(gmx_domdec_t
*dd
);
48 /* Is the ns grid already filled with the home particles? */
50 void dd_store_state(gmx_domdec_t
*dd
,t_state
*state
);
51 /* Store the global cg indices of the home cgs in state,
52 * so it can be reset, even after a new DD partitioning.
55 gmx_domdec_zones_t
*domdec_zones(gmx_domdec_t
*dd
);
57 void dd_get_ns_ranges(gmx_domdec_t
*dd
,int icg
,
58 int *jcg0
,int *jcg1
,ivec shift0
,ivec shift1
);
60 int dd_natoms_vsite(gmx_domdec_t
*dd
);
62 void dd_get_constraint_range(gmx_domdec_t
*dd
,
63 int *at_start
,int *at_end
);
65 real
dd_cutoff_mbody(gmx_domdec_t
*dd
);
67 real
dd_cutoff_twobody(gmx_domdec_t
*dd
);
69 gmx_bool
gmx_pmeonlynode(t_commrec
*cr
,int nodeid
);
70 /* Return if nodeid in cr->mpi_comm_mysim is a PME-only node */
72 void get_pme_ddnodes(t_commrec
*cr
,int pmenodeid
,
73 int *nmy_ddnodes
,int **my_ddnodes
,int *node_peer
);
74 /* Returns the set of DD nodes that communicate with pme node cr->nodeid */
76 int dd_pme_maxshift_x(gmx_domdec_t
*dd
);
77 /* Returns the maximum shift for coordinate communication in PME, dim x */
79 int dd_pme_maxshift_y(gmx_domdec_t
*dd
);
80 /* Returns the maximum shift for coordinate communication in PME, dim y */
82 void make_dd_communicators(FILE *fplog
,t_commrec
*cr
,int dd_node_order
);
85 init_domain_decomposition(FILE *fplog
,
89 real comm_distance_min
,real rconstr
,
90 const char *dlb_opt
,real dlb_scale
,
91 const char *sizex
,const char *sizey
,const char *sizez
,
92 gmx_mtop_t
*mtop
,t_inputrec
*ir
,
95 int *npme_x
, int *npme_y
);
97 void dd_init_bondeds(FILE *fplog
,
98 gmx_domdec_t
*dd
,gmx_mtop_t
*mtop
,
99 gmx_vsite_t
*vsite
,gmx_constr_t constr
,
100 t_inputrec
*ir
,gmx_bool bBCheck
,cginfo_mb_t
*cginfo_mb
);
101 /* Initialize data structures for bonded interactions */
103 void set_dd_parameters(FILE *fplog
,gmx_domdec_t
*dd
,real dlb_scale
,
104 t_inputrec
*ir
,t_forcerec
*fr
,
106 /* Set DD grid dimensions and limits,
107 * should be called after calling dd_init_bondeds.
110 void setup_dd_grid(FILE *fplog
,gmx_domdec_t
*dd
);
112 void dd_collect_vec(gmx_domdec_t
*dd
,
113 t_state
*state_local
,rvec
*lv
,rvec
*v
);
115 void dd_collect_state(gmx_domdec_t
*dd
,
116 t_state
*state_local
,t_state
*state
);
118 enum { ddCyclStep
, ddCyclPPduringPME
, ddCyclF
, ddCyclPME
, ddCyclNr
};
120 void dd_cycles_add(gmx_domdec_t
*dd
,float cycles
,int ddCycl
);
121 /* Add the wallcycle count to the DD counter */
123 void dd_force_flop_start(gmx_domdec_t
*dd
,t_nrnb
*nrnb
);
124 /* Start the force flop count */
126 void dd_force_flop_stop(gmx_domdec_t
*dd
,t_nrnb
*nrnb
);
127 /* Stop the force flop count */
129 void dd_move_x(gmx_domdec_t
*dd
,matrix box
,rvec x
[]);
130 /* Communicate the coordinates to the neighboring cells and do pbc. */
132 void dd_move_f(gmx_domdec_t
*dd
,rvec f
[],rvec
*fshift
);
133 /* Sum the forces over the neighboring cells.
134 * When fshift!=NULL the shift forces are updated to obtain
135 * the correct virial from the single sum including f.
138 void dd_atom_spread_real(gmx_domdec_t
*dd
,real v
[]);
139 /* Communicate a real for each atom to the neighboring cells. */
141 void dd_atom_sum_real(gmx_domdec_t
*dd
,real v
[]);
142 /* Sum the contributions to a real for each atom over the neighboring cells. */
144 void init_local_rigid_groups(t_inputrec
*ir
, t_mdatoms
*mdatoms
, t_state
*state_local
);
145 /* Initialize the rigid body groups */
147 void get_local_rigid_groups(t_inputrec
*ir
, t_mdatoms
*mdatoms
);
148 /* Update the rigid body groups according to the current domain decomposition */
150 void dd_partition_system(FILE *fplog
,
151 gmx_large_int_t step
,
153 gmx_bool bMasterState
,
155 t_state
*state_global
,
156 gmx_mtop_t
*top_global
,
158 t_state
*state_local
,
161 gmx_localtop_t
*top_local
,
164 gmx_shellfc_t shellfc
,
167 gmx_wallcycle_t wcycle
,
169 /* Partition the system over the nodes.
170 * step is only used for printing error messages.
171 * If bMasterState==TRUE then state_global from the master node is used,
172 * else state_local is redistributed between the nodes.
173 * When f!=NULL, *f will be reallocated to the size of state_local.
176 void reset_dd_statistics_counters(gmx_domdec_t
*dd
);
177 /* Reset all the statistics and counters for total run counting */
179 void print_dd_statistics(t_commrec
*cr
,t_inputrec
*ir
,FILE *fplog
);
181 /* In domdec_con.c */
183 void dd_move_f_vsites(gmx_domdec_t
*dd
,rvec
*f
,rvec
*fshift
);
185 void dd_clear_f_vsites(gmx_domdec_t
*dd
,rvec
*f
);
187 void dd_move_x_constraints(gmx_domdec_t
*dd
,matrix box
,
189 /* Move x0 and also x1 if x1!=NULL */
191 void dd_move_x_vsites(gmx_domdec_t
*dd
,matrix box
,rvec
*x
);
193 int *dd_constraints_nlocalatoms(gmx_domdec_t
*dd
);
195 void dd_clear_local_constraint_indices(gmx_domdec_t
*dd
);
197 void dd_clear_local_vsite_indices(gmx_domdec_t
*dd
);
199 int dd_make_local_vsites(gmx_domdec_t
*dd
,int at_start
,t_ilist
*lil
);
201 int dd_make_local_constraints(gmx_domdec_t
*dd
,int at_start
,
203 gmx_constr_t constr
,int nrec
,
206 void init_domdec_constraints(gmx_domdec_t
*dd
,
207 int natoms
,gmx_mtop_t
*mtop
,
208 gmx_constr_t constr
);
210 void init_domdec_vsites(gmx_domdec_t
*dd
,int natoms
);
213 /* In domdec_top.c */
215 void dd_print_missing_interactions(FILE *fplog
,t_commrec
*cr
,
216 int local_count
, gmx_mtop_t
*top_global
, t_state
*state_local
);
218 void dd_make_reverse_top(FILE *fplog
,
219 gmx_domdec_t
*dd
,gmx_mtop_t
*mtop
,
220 gmx_vsite_t
*vsite
,gmx_constr_t constr
,
221 t_inputrec
*ir
,gmx_bool bBCheck
);
223 void dd_make_local_cgs(gmx_domdec_t
*dd
,t_block
*lcgs
);
225 void dd_make_local_top(FILE *fplog
,
226 gmx_domdec_t
*dd
,gmx_domdec_zones_t
*zones
,
227 int npbcdim
,matrix box
,
228 rvec cellsize_min
,ivec npulse
,
229 t_forcerec
*fr
,gmx_vsite_t
*vsite
,
230 gmx_mtop_t
*top
,gmx_localtop_t
*ltop
);
232 void dd_sort_local_top(gmx_domdec_t
*dd
,t_mdatoms
*mdatoms
,
233 gmx_localtop_t
*ltop
);
234 /* Sort ltop->ilist when we are doing free energy. */
236 gmx_localtop_t
*dd_init_local_top(gmx_mtop_t
*top_global
);
238 void dd_init_local_state(gmx_domdec_t
*dd
,
239 t_state
*state_global
,t_state
*local_state
);
241 t_blocka
*make_charge_group_links(gmx_mtop_t
*mtop
,gmx_domdec_t
*dd
,
242 cginfo_mb_t
*cginfo_mb
);
244 void dd_bonded_cg_distance(FILE *fplog
,
245 gmx_domdec_t
*dd
,gmx_mtop_t
*mtop
,
246 t_inputrec
*ir
,rvec
*x
,matrix box
,
248 real
*r_2b
,real
*r_mb
);
250 void write_dd_pdb(const char *fn
,gmx_large_int_t step
,const char *title
,
253 int natoms
,rvec x
[],matrix box
);
254 /* Dump a pdb file with the current DD home + communicated atoms.
255 * When natoms=-1, dump all known atoms.
259 /* In domdec_setup.c */
261 real
comm_box_frac(ivec dd_nc
,real cutoff
,gmx_ddbox_t
*ddbox
);
262 /* Returns the volume fraction of the system that is communicated */
264 real
dd_choose_grid(FILE *fplog
,
265 t_commrec
*cr
,gmx_domdec_t
*dd
,t_inputrec
*ir
,
266 gmx_mtop_t
*mtop
,matrix box
,gmx_ddbox_t
*ddbox
,
267 gmx_bool bDynLoadBal
,real dlb_scale
,
268 real cellsize_limit
,real cutoff_dd
,
269 gmx_bool bInterCGBondeds
,gmx_bool bInterCGMultiBody
);
270 /* Determines the optimal DD cell setup dd->nc and possibly npmenodes
272 * On the master node returns the actual cellsize limit used.
276 /* In domdec_box.c */
278 void set_ddbox(gmx_domdec_t
*dd
,gmx_bool bMasterState
,t_commrec
*cr_sum
,
279 t_inputrec
*ir
,matrix box
,
280 gmx_bool bCalcUnboundedSize
,t_block
*cgs
,rvec
*x
,
283 void set_ddbox_cr(t_commrec
*cr
,ivec
*dd_nc
,
284 t_inputrec
*ir
,matrix box
,t_block
*cgs
,rvec
*x
,
291 #endif /* _domdec_h */