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
30 int ddglatnr(gmx_domdec_t
*dd
,int i
);
31 /* Returns the global topology atom number belonging to local atom index i.
32 * This function is intended for writing ascii output
33 * and returns atom numbers starting at 1.
34 * When dd=NULL returns i+1.
37 t_block
*dd_charge_groups_global(gmx_domdec_t
*dd
);
38 /* Return a block struct for the charge groups of the whole system */
40 gmx_bool
dd_filled_nsgrid_home(gmx_domdec_t
*dd
);
41 /* Is the ns grid already filled with the home particles? */
43 void dd_store_state(gmx_domdec_t
*dd
,t_state
*state
);
44 /* Store the global cg indices of the home cgs in state,
45 * so it can be reset, even after a new DD partitioning.
48 gmx_domdec_zones_t
*domdec_zones(gmx_domdec_t
*dd
);
50 void dd_get_ns_ranges(gmx_domdec_t
*dd
,int icg
,
51 int *jcg0
,int *jcg1
,ivec shift0
,ivec shift1
);
53 int dd_natoms_vsite(gmx_domdec_t
*dd
);
55 void dd_get_constraint_range(gmx_domdec_t
*dd
,
56 int *at_start
,int *at_end
);
58 real
dd_cutoff_mbody(gmx_domdec_t
*dd
);
60 real
dd_cutoff_twobody(gmx_domdec_t
*dd
);
62 gmx_bool
gmx_pmeonlynode(t_commrec
*cr
,int nodeid
);
63 /* Return if nodeid in cr->mpi_comm_mysim is a PME-only node */
65 void get_pme_ddnodes(t_commrec
*cr
,int pmenodeid
,
66 int *nmy_ddnodes
,int **my_ddnodes
,int *node_peer
);
67 /* Returns the set of DD nodes that communicate with pme node cr->nodeid */
69 int dd_pme_maxshift_x(gmx_domdec_t
*dd
);
70 /* Returns the maximum shift for coordinate communication in PME, dim x */
72 int dd_pme_maxshift_y(gmx_domdec_t
*dd
);
73 /* Returns the maximum shift for coordinate communication in PME, dim y */
75 void make_dd_communicators(FILE *fplog
,t_commrec
*cr
,int dd_node_order
);
78 init_domain_decomposition(FILE *fplog
,
82 real comm_distance_min
,real rconstr
,
83 const char *dlb_opt
,real dlb_scale
,
84 const char *sizex
,const char *sizey
,const char *sizez
,
85 gmx_mtop_t
*mtop
,t_inputrec
*ir
,
88 int *npme_x
, int *npme_y
);
90 void dd_init_bondeds(FILE *fplog
,
91 gmx_domdec_t
*dd
,gmx_mtop_t
*mtop
,
92 gmx_vsite_t
*vsite
,gmx_constr_t constr
,
93 t_inputrec
*ir
,gmx_bool bBCheck
,cginfo_mb_t
*cginfo_mb
);
94 /* Initialize data structures for bonded interactions */
96 void set_dd_parameters(FILE *fplog
,gmx_domdec_t
*dd
,real dlb_scale
,
97 t_inputrec
*ir
,t_forcerec
*fr
,
99 /* Set DD grid dimensions and limits,
100 * should be called after calling dd_init_bondeds.
103 void setup_dd_grid(FILE *fplog
,gmx_domdec_t
*dd
);
105 void dd_collect_vec(gmx_domdec_t
*dd
,
106 t_state
*state_local
,rvec
*lv
,rvec
*v
);
108 void dd_collect_state(gmx_domdec_t
*dd
,
109 t_state
*state_local
,t_state
*state
);
111 enum { ddCyclStep
, ddCyclPPduringPME
, ddCyclF
, ddCyclPME
, ddCyclNr
};
113 void dd_cycles_add(gmx_domdec_t
*dd
,float cycles
,int ddCycl
);
114 /* Add the wallcycle count to the DD counter */
116 void dd_force_flop_start(gmx_domdec_t
*dd
,t_nrnb
*nrnb
);
117 /* Start the force flop count */
119 void dd_force_flop_stop(gmx_domdec_t
*dd
,t_nrnb
*nrnb
);
120 /* Stop the force flop count */
122 void dd_move_x(gmx_domdec_t
*dd
,matrix box
,rvec x
[]);
123 /* Communicate the coordinates to the neighboring cells and do pbc. */
125 void dd_move_f(gmx_domdec_t
*dd
,rvec f
[],rvec
*fshift
);
126 /* Sum the forces over the neighboring cells.
127 * When fshift!=NULL the shift forces are updated to obtain
128 * the correct virial from the single sum including f.
131 void dd_atom_spread_real(gmx_domdec_t
*dd
,real v
[]);
132 /* Communicate a real for each atom to the neighboring cells. */
134 void dd_atom_sum_real(gmx_domdec_t
*dd
,real v
[]);
135 /* Sum the contributions to a real for each atom over the neighboring cells. */
137 void dd_partition_system(FILE *fplog
,
138 gmx_large_int_t step
,
140 gmx_bool bMasterState
,
142 t_state
*state_global
,
143 gmx_mtop_t
*top_global
,
145 t_state
*state_local
,
148 gmx_localtop_t
*top_local
,
151 gmx_shellfc_t shellfc
,
154 gmx_wallcycle_t wcycle
,
156 /* Partition the system over the nodes.
157 * step is only used for printing error messages.
158 * If bMasterState==TRUE then state_global from the master node is used,
159 * else state_local is redistributed between the nodes.
160 * When f!=NULL, *f will be reallocated to the size of state_local.
163 void reset_dd_statistics_counters(gmx_domdec_t
*dd
);
164 /* Reset all the statistics and counters for total run counting */
166 void print_dd_statistics(t_commrec
*cr
,t_inputrec
*ir
,FILE *fplog
);
168 /* In domdec_con.c */
170 void dd_move_f_vsites(gmx_domdec_t
*dd
,rvec
*f
,rvec
*fshift
);
172 void dd_clear_f_vsites(gmx_domdec_t
*dd
,rvec
*f
);
174 void dd_move_x_constraints(gmx_domdec_t
*dd
,matrix box
,
176 /* Move x0 and also x1 if x1!=NULL */
178 void dd_move_x_vsites(gmx_domdec_t
*dd
,matrix box
,rvec
*x
);
180 int *dd_constraints_nlocalatoms(gmx_domdec_t
*dd
);
182 void dd_clear_local_constraint_indices(gmx_domdec_t
*dd
);
184 void dd_clear_local_vsite_indices(gmx_domdec_t
*dd
);
186 int dd_make_local_vsites(gmx_domdec_t
*dd
,int at_start
,t_ilist
*lil
);
188 int dd_make_local_constraints(gmx_domdec_t
*dd
,int at_start
,
190 gmx_constr_t constr
,int nrec
,
193 void init_domdec_constraints(gmx_domdec_t
*dd
,
194 int natoms
,gmx_mtop_t
*mtop
,
195 gmx_constr_t constr
);
197 void init_domdec_vsites(gmx_domdec_t
*dd
,int natoms
);
200 /* In domdec_top.c */
202 void dd_print_missing_interactions(FILE *fplog
,t_commrec
*cr
,
203 int local_count
, gmx_mtop_t
*top_global
, t_state
*state_local
);
205 void dd_make_reverse_top(FILE *fplog
,
206 gmx_domdec_t
*dd
,gmx_mtop_t
*mtop
,
207 gmx_vsite_t
*vsite
,gmx_constr_t constr
,
208 t_inputrec
*ir
,gmx_bool bBCheck
);
210 void dd_make_local_cgs(gmx_domdec_t
*dd
,t_block
*lcgs
);
212 void dd_make_local_top(FILE *fplog
,
213 gmx_domdec_t
*dd
,gmx_domdec_zones_t
*zones
,
214 int npbcdim
,matrix box
,
215 rvec cellsize_min
,ivec npulse
,
216 t_forcerec
*fr
,gmx_vsite_t
*vsite
,
217 gmx_mtop_t
*top
,gmx_localtop_t
*ltop
);
219 void dd_sort_local_top(gmx_domdec_t
*dd
,t_mdatoms
*mdatoms
,
220 gmx_localtop_t
*ltop
);
221 /* Sort ltop->ilist when we are doing free energy. */
223 gmx_localtop_t
*dd_init_local_top(gmx_mtop_t
*top_global
);
225 void dd_init_local_state(gmx_domdec_t
*dd
,
226 t_state
*state_global
,t_state
*local_state
);
228 t_blocka
*make_charge_group_links(gmx_mtop_t
*mtop
,gmx_domdec_t
*dd
,
229 cginfo_mb_t
*cginfo_mb
);
231 void dd_bonded_cg_distance(FILE *fplog
,
232 gmx_domdec_t
*dd
,gmx_mtop_t
*mtop
,
233 t_inputrec
*ir
,rvec
*x
,matrix box
,
235 real
*r_2b
,real
*r_mb
);
237 void write_dd_pdb(const char *fn
,gmx_large_int_t step
,const char *title
,
240 int natoms
,rvec x
[],matrix box
);
241 /* Dump a pdb file with the current DD home + communicated atoms.
242 * When natoms=-1, dump all known atoms.
246 /* In domdec_setup.c */
248 real
comm_box_frac(ivec dd_nc
,real cutoff
,gmx_ddbox_t
*ddbox
);
249 /* Returns the volume fraction of the system that is communicated */
251 real
dd_choose_grid(FILE *fplog
,
252 t_commrec
*cr
,gmx_domdec_t
*dd
,t_inputrec
*ir
,
253 gmx_mtop_t
*mtop
,matrix box
,gmx_ddbox_t
*ddbox
,
254 gmx_bool bDynLoadBal
,real dlb_scale
,
255 real cellsize_limit
,real cutoff_dd
,
256 gmx_bool bInterCGBondeds
,gmx_bool bInterCGMultiBody
);
257 /* Determines the optimal DD cell setup dd->nc and possibly npmenodes
259 * On the master node returns the actual cellsize limit used.
263 /* In domdec_box.c */
265 void set_ddbox(gmx_domdec_t
*dd
,gmx_bool bMasterState
,t_commrec
*cr_sum
,
266 t_inputrec
*ir
,matrix box
,
267 gmx_bool bCalcUnboundedSize
,t_block
*cgs
,rvec
*x
,
270 void set_ddbox_cr(t_commrec
*cr
,ivec
*dd_nc
,
271 t_inputrec
*ir
,matrix box
,t_block
*cgs
,rvec
*x
,
278 #endif /* _domdec_h */