Removed refs to charmm files in share/top/Makefile.am
[gromacs/qmmm-gamess-us.git] / include / domdec.h
blobb1813fa107cd0543ccc9dc52470486b37fcd565f
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 *
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
15 * And Hey:
16 * Gnomes, ROck Monsters And Chili Sauce
19 #ifdef HAVE_CONFIG_H
20 #include <config.h>
21 #endif
23 #ifndef _domdec_h
24 #define _domdec_h
26 #include "typedefs.h"
27 #include "vsite.h"
28 #include "genborn.h"
30 #ifdef GMX_LIB_MPI
31 #include <mpi.h>
32 #endif
33 #ifdef GMX_THREADS
34 #include "tmpi.h"
35 #endif
37 extern 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 extern t_block *dd_charge_groups_global(gmx_domdec_t *dd);
45 /* Return a block struct for the charge groups of the whole system */
47 extern bool dd_filled_nsgrid_home(gmx_domdec_t *dd);
48 /* Is the ns grid already filled with the home particles? */
50 extern 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 extern gmx_domdec_zones_t *domdec_zones(gmx_domdec_t *dd);
57 extern void dd_get_ns_ranges(gmx_domdec_t *dd,int icg,
58 int *jcg0,int *jcg1,ivec shift0,ivec shift1);
60 extern int dd_natoms_vsite(gmx_domdec_t *dd);
62 extern void dd_get_constraint_range(gmx_domdec_t *dd,
63 int *at_start,int *at_end);
65 extern real dd_cutoff_mbody(gmx_domdec_t *dd);
67 extern real dd_cutoff_twobody(gmx_domdec_t *dd);
69 extern bool gmx_pmeonlynode(t_commrec *cr,int nodeid);
70 /* Return if nodeid in cr->mpi_comm_mysim is a PME-only node */
72 extern 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 extern int dd_pme_maxshift0(gmx_domdec_t *dd);
77 /* Returns the maximum shift for coordinate communication in PME, dim 0 */
79 extern int dd_pme_maxshift1(gmx_domdec_t *dd);
80 /* Returns the maximum shift for coordinate communication in PME, dim 1 */
82 extern void make_dd_communicators(FILE *fplog,t_commrec *cr,int dd_node_order);
84 extern gmx_domdec_t *
85 init_domain_decomposition(FILE *fplog,
86 t_commrec *cr,
87 unsigned long Flags,
88 ivec nc,
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,
93 matrix box,rvec *x,
94 gmx_ddbox_t *ddbox,
95 int *npme_major);
97 extern 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,bool bBCheck,cginfo_mb_t *cginfo_mb);
101 /* Initialize data structures for bonded interactions */
103 extern void set_dd_parameters(FILE *fplog,gmx_domdec_t *dd,real dlb_scale,
104 t_inputrec *ir,t_forcerec *fr,
105 gmx_ddbox_t *ddbox);
106 /* Set DD grid dimensions and limits,
107 * should be called after calling dd_init_bondeds.
110 extern void setup_dd_grid(FILE *fplog,gmx_domdec_t *dd);
112 extern void dd_collect_vec(gmx_domdec_t *dd,
113 t_state *state_local,rvec *lv,rvec *v);
115 extern void dd_collect_state(gmx_domdec_t *dd,
116 t_state *state_local,t_state *state);
118 enum { ddCyclStep, ddCyclPPduringPME, ddCyclF, ddCyclPME, ddCyclNr };
120 extern void dd_cycles_add(gmx_domdec_t *dd,float cycles,int ddCycl);
121 /* Add the wallcycle count to the DD counter */
123 extern void dd_force_flop_start(gmx_domdec_t *dd,t_nrnb *nrnb);
124 /* Start the force flop count */
126 extern void dd_force_flop_stop(gmx_domdec_t *dd,t_nrnb *nrnb);
127 /* Stop the force flop count */
129 extern void dd_move_x(gmx_domdec_t *dd,matrix box,rvec x[]);
130 /* Communicate the coordinates to the neighboring cells and do pbc. */
132 extern 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 extern void dd_atom_spread_real(gmx_domdec_t *dd,real v[]);
139 /* Communicate a real for each atom to the neighboring cells. */
141 extern 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 extern void dd_partition_system(FILE *fplog,
145 gmx_large_int_t step,
146 t_commrec *cr,
147 bool bMasterState,
148 int nstglobalcomm,
149 t_state *state_global,
150 gmx_mtop_t *top_global,
151 t_inputrec *ir,
152 t_state *state_local,
153 rvec **f,
154 t_mdatoms *mdatoms,
155 gmx_localtop_t *top_local,
156 t_forcerec *fr,
157 gmx_vsite_t *vsite,
158 gmx_shellfc_t shellfc,
159 gmx_constr_t constr,
160 t_nrnb *nrnb,
161 gmx_wallcycle_t wcycle,
162 bool bVerbose);
163 /* Partition the system over the nodes.
164 * step is only used for printing error messages.
165 * If bMasterState==TRUE then state_global from the master node is used,
166 * else state_local is redistributed between the nodes.
167 * When f!=NULL, *f will be reallocated to the size of state_local.
170 extern void reset_dd_statistics_counters(gmx_domdec_t *dd);
171 /* Reset all the statistics and counters for total run counting */
173 extern void print_dd_statistics(t_commrec *cr,t_inputrec *ir,FILE *fplog);
175 /* In domdec_con.c */
177 extern void dd_move_f_vsites(gmx_domdec_t *dd,rvec *f,rvec *fshift);
179 extern void dd_clear_f_vsites(gmx_domdec_t *dd,rvec *f);
181 extern void dd_move_x_constraints(gmx_domdec_t *dd,matrix box,
182 rvec *x0,rvec *x1);
183 /* Move x0 and also x1 if x1!=NULL */
185 extern void dd_move_x_vsites(gmx_domdec_t *dd,matrix box,rvec *x);
187 extern int *dd_constraints_nlocalatoms(gmx_domdec_t *dd);
189 extern void dd_clear_local_constraint_indices(gmx_domdec_t *dd);
191 extern void dd_clear_local_vsite_indices(gmx_domdec_t *dd);
193 extern int dd_make_local_vsites(gmx_domdec_t *dd,int at_start,t_ilist *lil);
195 extern int dd_make_local_constraints(gmx_domdec_t *dd,int at_start,
196 gmx_mtop_t *mtop,
197 gmx_constr_t constr,int nrec,
198 t_ilist *il_local);
200 extern void init_domdec_constraints(gmx_domdec_t *dd,
201 int natoms,gmx_mtop_t *mtop,
202 gmx_constr_t constr);
204 extern void init_domdec_vsites(gmx_domdec_t *dd,int natoms);
207 /* In domdec_top.c */
209 extern void dd_print_missing_interactions(FILE *fplog,t_commrec *cr,
210 int local_count, gmx_mtop_t *top_global, t_state *state_local);
212 extern void dd_make_reverse_top(FILE *fplog,
213 gmx_domdec_t *dd,gmx_mtop_t *mtop,
214 gmx_vsite_t *vsite,gmx_constr_t constr,
215 t_inputrec *ir,bool bBCheck);
217 extern void dd_make_local_cgs(gmx_domdec_t *dd,t_block *lcgs);
219 extern void dd_make_local_top(FILE *fplog,
220 gmx_domdec_t *dd,gmx_domdec_zones_t *zones,
221 int npbcdim,matrix box,
222 rvec cellsize_min,ivec npulse,
223 t_forcerec *fr,gmx_vsite_t *vsite,
224 gmx_mtop_t *top,gmx_localtop_t *ltop);
226 extern gmx_localtop_t *dd_init_local_top(gmx_mtop_t *top_global);
228 extern void dd_init_local_state(gmx_domdec_t *dd,
229 t_state *state_global,t_state *local_state);
231 extern t_blocka *make_charge_group_links(gmx_mtop_t *mtop,gmx_domdec_t *dd,
232 cginfo_mb_t *cginfo_mb);
234 extern void dd_bonded_cg_distance(FILE *fplog,
235 gmx_domdec_t *dd,gmx_mtop_t *mtop,
236 t_inputrec *ir,rvec *x,matrix box,
237 bool bBCheck,
238 real *r_2b,real *r_mb);
240 extern void write_dd_pdb(const char *fn,gmx_large_int_t step,const char *title,
241 gmx_mtop_t *mtop,
242 t_commrec *cr,
243 int natoms,rvec x[],matrix box);
244 /* Dump a pdb file with the current DD home + communicated atoms.
245 * When natoms=-1, dump all known atoms.
249 /* In domdec_setup.c */
251 extern real comm_box_frac(ivec dd_nc,real cutoff,gmx_ddbox_t *ddbox);
252 /* Returns the volume fraction of the system communicated by each node */
254 extern real dd_choose_grid(FILE *fplog,
255 t_commrec *cr,gmx_domdec_t *dd,t_inputrec *ir,
256 gmx_mtop_t *mtop,matrix box,gmx_ddbox_t *ddbox,
257 bool bDynLoadBal,real dlb_scale,
258 real cellsize_limit,real cutoff_dd,
259 bool bInterCGBondeds,bool bInterCGMultiBody);
260 /* Determines the optimal DD cell setup dd->nc and possibly npmenodes
261 * for the system.
262 * On the master node returns the actual cellsize limit used.
266 /* In domdec_box.c */
268 extern void set_ddbox(gmx_domdec_t *dd,bool bMasterState,t_commrec *cr_sum,
269 t_inputrec *ir,matrix box,
270 bool bCalcUnboundedSize,t_block *cgs,rvec *x,
271 gmx_ddbox_t *ddbox);
273 extern void set_ddbox_cr(t_commrec *cr,ivec *dd_nc,
274 t_inputrec *ir,matrix box,t_block *cgs,rvec *x,
275 gmx_ddbox_t *ddbox);
278 #endif /* _domdec_h */