PK added extern C for rmpbc and atomprop
[gromacs.git] / include / force.h
blob21c8c49f2d4ad901b5c141dad7a8df0c605eb678
1 /*
2 *
3 * This source code is part of
4 *
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
9 * VERSION 3.2.0
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * Gromacs Runs On Most of All Computer Systems
36 #ifndef _force_h
37 #define _force_h
39 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
44 #include "typedefs.h"
45 #include "pbc.h"
46 #include "network.h"
47 #include "tgroup.h"
48 #include "vsite.h"
49 #include "genborn.h"
51 static const char *sepdvdlformat=" %-30s V %12.5e dVdl %12.5e\n";
53 extern void calc_vir(FILE *fplog,int nxf,rvec x[],rvec f[],tensor vir,
54 bool bScrewPBC,matrix box);
55 /* Calculate virial for nxf atoms, and add it to vir */
57 extern void f_calc_vir(FILE *fplog,int i0,int i1,rvec x[],rvec f[],tensor vir,
58 t_graph *g,rvec shift_vec[]);
59 /* Calculate virial taking periodicity into account */
61 extern real RF_excl_correction(FILE *fplog,
62 const t_forcerec *fr,t_graph *g,
63 const t_mdatoms *mdatoms,const t_blocka *excl,
64 rvec x[],rvec f[],rvec *fshift,const t_pbc *pbc,
65 real lambda,real *dvdlambda);
66 /* Calculate the reaction-field energy correction for this node:
67 * epsfac q_i q_j (k_rf r_ij^2 - c_rf)
68 * and force correction for all excluded pairs, including self pairs.
71 extern void calc_rffac(FILE *fplog,int eel,real eps_r,real eps_rf,
72 real Rc,real Temp,
73 real zsq,matrix box,
74 real *kappa,real *krf,real *crf);
75 /* Determine the reaction-field constants */
77 extern void init_generalized_rf(FILE *fplog,
78 const gmx_mtop_t *mtop,const t_inputrec *ir,
79 t_forcerec *fr);
80 /* Initialize the generalized reaction field parameters */
83 /* In wall.c */
84 extern void make_wall_tables(FILE *fplog,
85 const t_inputrec *ir,const char *tabfn,
86 const gmx_groups_t *groups,
87 t_forcerec *fr);
89 extern real do_walls(t_inputrec *ir,t_forcerec *fr,matrix box,t_mdatoms *md,
90 rvec x[],rvec f[],real lambda,real Vlj[],t_nrnb *nrnb);
94 extern t_forcerec *mk_forcerec(void);
96 #define GMX_MAKETABLES_FORCEUSER (1<<0)
97 #define GMX_MAKETABLES_14ONLY (1<<1)
99 extern t_forcetable make_tables(FILE *fp,const t_forcerec *fr,
100 bool bVerbose,const char *fn,
101 real rtab,int flags);
102 /* Return tables for inner loops. When bVerbose the tables are printed
103 * to .xvg files
106 extern bondedtable_t make_bonded_table(FILE *fplog,char *fn,int angle);
107 /* Return a table for bonded interactions,
108 * angle should be: bonds 0, angles 1, dihedrals 2
111 /* Return a table for GB calculations */
112 extern t_forcetable make_gb_table(FILE *out,const t_forcerec *fr,
113 const char *fn,
114 real rtab);
116 extern void pr_forcerec(FILE *fplog,t_forcerec *fr,t_commrec *cr);
118 extern void
119 forcerec_set_ranges(t_forcerec *fr,
120 int ncg_home,int ncg_force,
121 int natoms_force,int natoms_f_novirsum);
122 /* Set the number of cg's and atoms for the force calculation */
124 extern void init_forcerec(FILE *fplog,
125 t_forcerec *fr,
126 t_fcdata *fcd,
127 const t_inputrec *ir,
128 const gmx_mtop_t *mtop,
129 const t_commrec *cr,
130 matrix box,
131 bool bMolEpot,
132 const char *tabfn,
133 const char *tabpfn,
134 const char *tabbfn,
135 bool bNoSolvOpt,
136 real print_force);
137 /* The Force rec struct must be created with mk_forcerec
138 * The booleans have the following meaning:
139 * bSetQ: Copy the charges [ only necessary when they change ]
140 * bMolEpot: Use the free energy stuff per molecule
141 * print_force >= 0: print forces for atoms with force >= print_force
144 extern void init_enerdata(int ngener,int n_flambda,gmx_enerdata_t *enerd);
145 /* Intializes the energy storage struct */
147 extern void destroy_enerdata(gmx_enerdata_t *enerd);
148 /* Free all memory associated with enerd */
150 extern void reset_enerdata(t_grpopts *opts,
151 t_forcerec *fr,bool bNS,
152 gmx_enerdata_t *enerd,
153 bool bMaster);
154 /* Resets the energy data, if bNS=TRUE also zeros the long-range part */
156 extern void sum_epot(t_grpopts *opts,gmx_enerdata_t *enerd);
157 /* Locally sum the non-bonded potential energy terms */
159 extern void sum_dhdl(gmx_enerdata_t *enerd,double lambda,t_inputrec *ir);
160 /* Sum the free energy contributions */
162 extern void update_forcerec(FILE *fplog,t_forcerec *fr,matrix box);
163 /* Updates parameters in the forcerec that are time dependent */
165 /* Compute the average C6 and C12 params for LJ corrections */
166 extern void set_avcsixtwelve(FILE *fplog,t_forcerec *fr,
167 const gmx_mtop_t *mtop);
169 /* The state has changed */
170 #define GMX_FORCE_STATECHANGED (1<<0)
171 /* Do neighbor searching */
172 #define GMX_FORCE_NS (1<<1)
173 /* Calculate bonded energies/forces */
174 #define GMX_FORCE_BONDED (1<<2)
175 /* Calculate non-bonded energies/forces */
176 #define GMX_FORCE_NONBONDED (1<<3)
177 /* Calculate forces (not only energies) */
178 #define GMX_FORCE_FORCES (1<<4)
179 /* Calculate the virial */
180 #define GMX_FORCE_VIRIAL (1<<5)
181 /* Calculate dHdl */
182 #define GMX_FORCE_DHDL (1<<6)
183 /* Normally one want all energy terms and forces */
184 #define GMX_FORCE_ALLFORCES (GMX_FORCE_BONDED | GMX_FORCE_NONBONDED | GMX_FORCE_FORCES)
186 extern void do_force(FILE *log,t_commrec *cr,
187 t_inputrec *inputrec,
188 gmx_step_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
189 gmx_localtop_t *top,
190 gmx_mtop_t *mtop,
191 gmx_groups_t *groups,
192 matrix box,rvec x[],history_t *hist,
193 rvec f[],
194 tensor vir_force,
195 t_mdatoms *mdatoms,
196 gmx_enerdata_t *enerd,t_fcdata *fcd,
197 real lambda,t_graph *graph,
198 t_forcerec *fr,gmx_vsite_t *vsite,rvec mu_tot,
199 double t,FILE *field,gmx_edsam_t ed,
200 bool bBornRadii,
201 int flags);
202 /* Communicate coordinates (if parallel).
203 * Do neighbor searching (if necessary).
204 * Calculate forces.
205 * Communicate forces (if parallel).
206 * Spread forces for vsites (if present).
208 * f is always required.
211 extern void ns(FILE *fplog,
212 t_forcerec *fr,
213 rvec x[],
214 rvec f[],
215 matrix box,
216 gmx_groups_t *groups,
217 t_grpopts *opts,
218 gmx_localtop_t *top,
219 t_mdatoms *md,
220 t_commrec *cr,
221 t_nrnb *nrnb,
222 real lambda,
223 real *dvdlambda,
224 gmx_grppairener_t *grppener,
225 bool bFillGrid,
226 bool bDoForces);
227 /* Call the neighborsearcher */
229 extern void do_force_lowlevel(FILE *fplog,
230 gmx_step_t step,
231 t_forcerec *fr,
232 t_inputrec *ir,
233 t_idef *idef,
234 t_commrec *cr,
235 t_nrnb *nrnb,
236 gmx_wallcycle_t wcycle,
237 t_mdatoms *md,
238 t_grpopts *opts,
239 rvec x[],
240 history_t *hist,
241 rvec f[],
242 gmx_enerdata_t *enerd,
243 t_fcdata *fcd,
244 gmx_mtop_t *mtop,
245 gmx_localtop_t *top,
246 gmx_genborn_t *born,
247 t_atomtypes *atype,
248 bool bBornRadii,
249 matrix box,
250 real lambda,
251 t_graph *graph,
252 t_blocka *excl,
253 rvec mu_tot[2],
254 int flags,
255 float *cycles_pme);
256 /* Call all the force routines */
259 #endif /* _force_h */