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
41 extern int ddglatnr(gmx_domdec_t
*dd
,int i
);
42 /* Returns the global topology atom number belonging to local atom index i.
43 * This function is intended for writing ascii output
44 * and returns atom numbers starting at 1.
45 * When dd=NULL returns i+1.
48 extern t_block
*dd_charge_groups_global(gmx_domdec_t
*dd
);
49 /* Return a block struct for the charge groups of the whole system */
51 extern bool dd_filled_nsgrid_home(gmx_domdec_t
*dd
);
52 /* Is the ns grid already filled with the home particles? */
54 extern void dd_store_state(gmx_domdec_t
*dd
,t_state
*state
);
55 /* Store the global cg indices of the home cgs in state,
56 * so it can be reset, even after a new DD partitioning.
59 extern gmx_domdec_zones_t
*domdec_zones(gmx_domdec_t
*dd
);
61 extern void dd_get_ns_ranges(gmx_domdec_t
*dd
,int icg
,
62 int *jcg0
,int *jcg1
,ivec shift0
,ivec shift1
);
64 extern int dd_natoms_vsite(gmx_domdec_t
*dd
);
66 extern void dd_get_constraint_range(gmx_domdec_t
*dd
,
67 int *at_start
,int *at_end
);
69 extern real
dd_cutoff_mbody(gmx_domdec_t
*dd
);
71 extern real
dd_cutoff_twobody(gmx_domdec_t
*dd
);
73 extern bool gmx_pmeonlynode(t_commrec
*cr
,int nodeid
);
74 /* Return if nodeid in cr->mpi_comm_mysim is a PME-only node */
76 extern void get_pme_ddnodes(t_commrec
*cr
,int pmenodeid
,
77 int *nmy_ddnodes
,int **my_ddnodes
,int *node_peer
);
78 /* Returns the set of DD nodes that communicate with pme node cr->nodeid */
80 extern int dd_pme_maxshift0(gmx_domdec_t
*dd
);
81 /* Returns the maximum shift for coordinate communication in PME, dim 0 */
83 extern int dd_pme_maxshift1(gmx_domdec_t
*dd
);
84 /* Returns the maximum shift for coordinate communication in PME, dim 1 */
86 extern void make_dd_communicators(FILE *fplog
,t_commrec
*cr
,int dd_node_order
);
89 init_domain_decomposition(FILE *fplog
,
93 real comm_distance_min
,real rconstr
,
94 const char *dlb_opt
,real dlb_scale
,
95 const char *sizex
,const char *sizey
,const char *sizez
,
96 gmx_mtop_t
*mtop
,t_inputrec
*ir
,
101 extern void dd_init_bondeds(FILE *fplog
,
102 gmx_domdec_t
*dd
,gmx_mtop_t
*mtop
,
103 gmx_vsite_t
*vsite
,gmx_constr_t constr
,
104 t_inputrec
*ir
,bool bBCheck
,cginfo_mb_t
*cginfo_mb
);
105 /* Initialize data structures for bonded interactions */
107 extern void set_dd_parameters(FILE *fplog
,gmx_domdec_t
*dd
,real dlb_scale
,
108 t_inputrec
*ir
,t_forcerec
*fr
,
110 /* Set DD grid dimensions and limits,
111 * should be called after calling dd_init_bondeds.
114 extern void setup_dd_grid(FILE *fplog
,gmx_domdec_t
*dd
);
116 extern void dd_collect_vec(gmx_domdec_t
*dd
,
117 t_state
*state_local
,rvec
*lv
,rvec
*v
);
119 extern void dd_collect_state(gmx_domdec_t
*dd
,
120 t_state
*state_local
,t_state
*state
);
122 enum { ddCyclStep
, ddCyclPPduringPME
, ddCyclF
, ddCyclPME
, ddCyclNr
};
124 extern void dd_cycles_add(gmx_domdec_t
*dd
,float cycles
,int ddCycl
);
125 /* Add the wallcycle count to the DD counter */
127 extern void dd_force_flop_start(gmx_domdec_t
*dd
,t_nrnb
*nrnb
);
128 /* Start the force flop count */
130 extern void dd_force_flop_stop(gmx_domdec_t
*dd
,t_nrnb
*nrnb
);
131 /* Stop the force flop count */
133 extern void dd_move_x(gmx_domdec_t
*dd
,matrix box
,rvec x
[]);
134 /* Communicate the coordinates to the neighboring cells and do pbc. */
136 extern void dd_move_f(gmx_domdec_t
*dd
,rvec f
[],rvec
*fshift
);
137 /* Sum the forces over the neighboring cells.
138 * When fshift!=NULL the shift forces are updated to obtain
139 * the correct virial from the single sum including f.
142 extern void dd_atom_spread_real(gmx_domdec_t
*dd
,real v
[]);
143 /* Communicate a real for each atom to the neighboring cells. */
145 extern void dd_atom_sum_real(gmx_domdec_t
*dd
,real v
[]);
146 /* Sum the contributions to a real for each atom over the neighboring cells. */
148 extern void dd_partition_system(FILE *fplog
,
149 gmx_large_int_t step
,
153 t_state
*state_global
,
154 gmx_mtop_t
*top_global
,
156 t_state
*state_local
,
159 gmx_localtop_t
*top_local
,
162 gmx_shellfc_t shellfc
,
165 gmx_wallcycle_t wcycle
,
167 /* Partition the system over the nodes.
168 * step is only used for printing error messages.
169 * If bMasterState==TRUE then state_global from the master node is used,
170 * else state_local is redistributed between the nodes.
171 * When f!=NULL, *f will be reallocated to the size of state_local.
174 extern void reset_dd_statistics_counters(gmx_domdec_t
*dd
);
175 /* Reset all the statistics and counters for total run counting */
177 extern void print_dd_statistics(t_commrec
*cr
,t_inputrec
*ir
,FILE *fplog
);
179 /* In domdec_con.c */
181 extern void dd_move_f_vsites(gmx_domdec_t
*dd
,rvec
*f
,rvec
*fshift
);
183 extern void dd_clear_f_vsites(gmx_domdec_t
*dd
,rvec
*f
);
185 extern void dd_move_x_constraints(gmx_domdec_t
*dd
,matrix box
,
187 /* Move x0 and also x1 if x1!=NULL */
189 extern void dd_move_x_vsites(gmx_domdec_t
*dd
,matrix box
,rvec
*x
);
191 extern int *dd_constraints_nlocalatoms(gmx_domdec_t
*dd
);
193 extern void dd_clear_local_constraint_indices(gmx_domdec_t
*dd
);
195 extern void dd_clear_local_vsite_indices(gmx_domdec_t
*dd
);
197 extern int dd_make_local_vsites(gmx_domdec_t
*dd
,int at_start
,t_ilist
*lil
);
199 extern int dd_make_local_constraints(gmx_domdec_t
*dd
,int at_start
,
201 gmx_constr_t constr
,int nrec
,
204 extern void init_domdec_constraints(gmx_domdec_t
*dd
,
205 int natoms
,gmx_mtop_t
*mtop
,
206 gmx_constr_t constr
);
208 extern void init_domdec_vsites(gmx_domdec_t
*dd
,int natoms
);
211 /* In domdec_top.c */
213 extern void dd_print_missing_interactions(FILE *fplog
,t_commrec
*cr
,
214 int local_count
, gmx_mtop_t
*top_global
, t_state
*state_local
);
216 extern void dd_make_reverse_top(FILE *fplog
,
217 gmx_domdec_t
*dd
,gmx_mtop_t
*mtop
,
218 gmx_vsite_t
*vsite
,gmx_constr_t constr
,
219 t_inputrec
*ir
,bool bBCheck
);
221 extern void dd_make_local_cgs(gmx_domdec_t
*dd
,t_block
*lcgs
);
223 extern void dd_make_local_top(FILE *fplog
,
224 gmx_domdec_t
*dd
,gmx_domdec_zones_t
*zones
,
225 int npbcdim
,matrix box
,
226 rvec cellsize_min
,ivec npulse
,
227 t_forcerec
*fr
,gmx_vsite_t
*vsite
,
228 gmx_mtop_t
*top
,gmx_localtop_t
*ltop
);
230 extern gmx_localtop_t
*dd_init_local_top(gmx_mtop_t
*top_global
);
232 extern void dd_init_local_state(gmx_domdec_t
*dd
,
233 t_state
*state_global
,t_state
*local_state
);
235 extern t_blocka
*make_charge_group_links(gmx_mtop_t
*mtop
,gmx_domdec_t
*dd
,
236 cginfo_mb_t
*cginfo_mb
);
238 extern void dd_bonded_cg_distance(FILE *fplog
,
239 gmx_domdec_t
*dd
,gmx_mtop_t
*mtop
,
240 t_inputrec
*ir
,rvec
*x
,matrix box
,
242 real
*r_2b
,real
*r_mb
);
244 extern void write_dd_pdb(const char *fn
,gmx_large_int_t step
,const char *title
,
247 int natoms
,rvec x
[],matrix box
);
248 /* Dump a pdb file with the current DD home + communicated atoms.
249 * When natoms=-1, dump all known atoms.
253 /* In domdec_setup.c */
255 extern real
comm_box_frac(ivec dd_nc
,real cutoff
,gmx_ddbox_t
*ddbox
);
256 /* Returns the volume fraction of the system communicated by each node */
258 extern real
dd_choose_grid(FILE *fplog
,
259 t_commrec
*cr
,gmx_domdec_t
*dd
,t_inputrec
*ir
,
260 gmx_mtop_t
*mtop
,matrix box
,gmx_ddbox_t
*ddbox
,
261 bool bDynLoadBal
,real dlb_scale
,
262 real cellsize_limit
,real cutoff_dd
,
263 bool bInterCGBondeds
,bool bInterCGMultiBody
);
264 /* Determines the optimal DD cell setup dd->nc and possibly npmenodes
266 * On the master node returns the actual cellsize limit used.
270 /* In domdec_box.c */
272 extern void set_ddbox(gmx_domdec_t
*dd
,bool bMasterState
,t_commrec
*cr_sum
,
273 t_inputrec
*ir
,matrix box
,
274 bool bCalcUnboundedSize
,t_block
*cgs
,rvec
*x
,
277 extern void set_ddbox_cr(t_commrec
*cr
,ivec
*dd_nc
,
278 t_inputrec
*ir
,matrix box
,t_block
*cgs
,rvec
*x
,
285 #endif /* _domdec_h */