automation for setting GMX_GPU & cmake GPU detection
[gromacs.git] / include / force.h
bloba91be65d17037aecb49f8ef95868abebd26e38c2
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 #include "visibility.h"
40 #include "typedefs.h"
41 #include "types/force_flags.h"
42 #include "pbc.h"
43 #include "network.h"
44 #include "tgroup.h"
45 #include "vsite.h"
46 #include "genborn.h"
49 #ifdef __cplusplus
50 extern "C" {
51 #endif
53 static const char *sepdvdlformat=" %-30s V %12.5e dVdl %12.5e\n";
55 void calc_vir(FILE *fplog,int nxf,rvec x[],rvec f[],tensor vir,
56 gmx_bool bScrewPBC,matrix box);
57 /* Calculate virial for nxf atoms, and add it to vir */
59 void f_calc_vir(FILE *fplog,int i0,int i1,rvec x[],rvec f[],tensor vir,
60 t_graph *g,rvec shift_vec[]);
61 /* Calculate virial taking periodicity into account */
63 real RF_excl_correction(FILE *fplog,
64 const t_forcerec *fr,t_graph *g,
65 const t_mdatoms *mdatoms,const t_blocka *excl,
66 rvec x[],rvec f[],rvec *fshift,const t_pbc *pbc,
67 real lambda,real *dvdlambda);
68 /* Calculate the reaction-field energy correction for this node:
69 * epsfac q_i q_j (k_rf r_ij^2 - c_rf)
70 * and force correction for all excluded pairs, including self pairs.
73 void calc_rffac(FILE *fplog,int eel,real eps_r,real eps_rf,
74 real Rc,real Temp,
75 real zsq,matrix box,
76 real *kappa,real *krf,real *crf);
77 /* Determine the reaction-field constants */
79 void init_generalized_rf(FILE *fplog,
80 const gmx_mtop_t *mtop,const t_inputrec *ir,
81 t_forcerec *fr);
82 /* Initialize the generalized reaction field parameters */
85 /* In wall.c */
86 void make_wall_tables(FILE *fplog,const output_env_t oenv,
87 const t_inputrec *ir,const char *tabfn,
88 const gmx_groups_t *groups,
89 t_forcerec *fr);
91 real do_walls(t_inputrec *ir,t_forcerec *fr,matrix box,t_mdatoms *md,
92 rvec x[],rvec f[],real lambda,real Vlj[],t_nrnb *nrnb);
94 GMX_LIBMD_EXPORT
95 t_forcerec *mk_forcerec(void);
97 #define GMX_MAKETABLES_FORCEUSER (1<<0)
98 #define GMX_MAKETABLES_14ONLY (1<<1)
100 t_forcetable make_tables(FILE *fp,const output_env_t oenv,
101 const t_forcerec *fr, gmx_bool bVerbose,
102 const char *fn, real rtab,int flags);
103 /* Return tables for inner loops. When bVerbose the tables are printed
104 * to .xvg files
107 bondedtable_t make_bonded_table(FILE *fplog,char *fn,int angle);
108 /* Return a table for bonded interactions,
109 * angle should be: bonds 0, angles 1, dihedrals 2
112 /* Return a table for GB calculations */
113 t_forcetable make_gb_table(FILE *out,const output_env_t oenv,
114 const t_forcerec *fr,
115 const char *fn,
116 real rtab);
118 /* Read a table for AdResS Thermo Force calculations */
119 extern t_forcetable make_atf_table(FILE *out,const output_env_t oenv,
120 const t_forcerec *fr,
121 const char *fn,
122 matrix box);
124 GMX_LIBMD_EXPORT
125 void pr_forcerec(FILE *fplog,t_forcerec *fr,t_commrec *cr);
127 void
128 forcerec_set_ranges(t_forcerec *fr,
129 int ncg_home,int ncg_force,
130 int natoms_force,
131 int natoms_force_constr,int natoms_f_novirsum);
132 /* Set the number of cg's and atoms for the force calculation */
134 GMX_LIBMD_EXPORT
135 gmx_bool can_use_allvsall(const t_inputrec *ir, const gmx_mtop_t *mtop,
136 gmx_bool bPrintNote,t_commrec *cr,FILE *fp);
137 /* Returns if we can use all-vs-all loops.
138 * If bPrintNote==TRUE, prints a note, if necessary, to stderr
139 * and fp (if !=NULL) on the master node.
142 GMX_LIBMD_EXPORT
143 gmx_bool uses_simple_tables(int cutoff_scheme,
144 nonbonded_verlet_t *nbv,
145 int group);
146 /* Returns whether simple tables (i.e. not for use with GPUs) are used
147 * with the type of kernel indicated.
150 GMX_LIBMD_EXPORT
151 void init_interaction_const_tables(FILE *fp,
152 interaction_const_t *ic,
153 gmx_bool bSimpleTable,
154 real rtab);
155 /* Initializes the tables in the interaction constant data structure.
156 * Setting verlet_kernel_type to -1 always initializes tables for
157 * use with group kernels.
160 void init_interaction_const(FILE *fp,
161 interaction_const_t **interaction_const,
162 const t_forcerec *fr,
163 real rtab);
164 /* Initializes the interaction constant data structure. Currently it
165 * uses forcerec as input.
168 GMX_LIBMD_EXPORT
169 void init_forcerec(FILE *fplog,
170 const output_env_t oenv,
171 t_forcerec *fr,
172 t_fcdata *fcd,
173 const t_inputrec *ir,
174 const gmx_mtop_t *mtop,
175 const t_commrec *cr,
176 matrix box,
177 gmx_bool bMolEpot,
178 const char *tabfn,
179 const char *tabafn,
180 const char *tabpfn,
181 const char *tabbfn,
182 const char *nbpu_opt,
183 gmx_bool bNoSolvOpt,
184 real print_force);
185 /* The Force rec struct must be created with mk_forcerec
186 * The gmx_booleans have the following meaning:
187 * bSetQ: Copy the charges [ only necessary when they change ]
188 * bMolEpot: Use the free energy stuff per molecule
189 * print_force >= 0: print forces for atoms with force >= print_force
192 GMX_LIBMD_EXPORT
193 void forcerec_set_excl_load(t_forcerec *fr,
194 const gmx_localtop_t *top,const t_commrec *cr);
195 /* Set the exclusion load for the local exclusions and possibly threads */
197 GMX_LIBMD_EXPORT
198 void init_enerdata(int ngener,int n_lambda,gmx_enerdata_t *enerd);
199 /* Intializes the energy storage struct */
201 void destroy_enerdata(gmx_enerdata_t *enerd);
202 /* Free all memory associated with enerd */
204 void reset_foreign_enerdata(gmx_enerdata_t *enerd);
205 /* Resets only the foreign energy data */
207 void reset_enerdata(t_grpopts *opts,
208 t_forcerec *fr,gmx_bool bNS,
209 gmx_enerdata_t *enerd,
210 gmx_bool bMaster);
211 /* Resets the energy data, if bNS=TRUE also zeros the long-range part */
213 void sum_epot(t_grpopts *opts, gmx_grppairener_t *grpp, real *epot);
214 /* Locally sum the non-bonded potential energy terms */
216 GMX_LIBMD_EXPORT
217 void sum_dhdl(gmx_enerdata_t *enerd,real *lambda,t_lambda *fepvals);
218 /* Sum the free energy contributions */
220 void update_forcerec(FILE *fplog,t_forcerec *fr,matrix box);
221 /* Updates parameters in the forcerec that are time dependent */
223 /* Compute the average C6 and C12 params for LJ corrections */
224 void set_avcsixtwelve(FILE *fplog,t_forcerec *fr,
225 const gmx_mtop_t *mtop);
227 GMX_LIBMD_EXPORT
228 extern void do_force(FILE *log,t_commrec *cr,
229 t_inputrec *inputrec,
230 gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
231 gmx_localtop_t *top,
232 gmx_mtop_t *mtop,
233 gmx_groups_t *groups,
234 matrix box,rvec x[],history_t *hist,
235 rvec f[],
236 tensor vir_force,
237 t_mdatoms *mdatoms,
238 gmx_enerdata_t *enerd,t_fcdata *fcd,
239 real *lambda,t_graph *graph,
240 t_forcerec *fr,
241 gmx_vsite_t *vsite,rvec mu_tot,
242 double t,FILE *field,gmx_edsam_t ed,
243 gmx_bool bBornRadii,
244 int flags);
246 /* Communicate coordinates (if parallel).
247 * Do neighbor searching (if necessary).
248 * Calculate forces.
249 * Communicate forces (if parallel).
250 * Spread forces for vsites (if present).
252 * f is always required.
255 GMX_LIBMD_EXPORT
256 void ns(FILE *fplog,
257 t_forcerec *fr,
258 rvec x[],
259 matrix box,
260 gmx_groups_t *groups,
261 t_grpopts *opts,
262 gmx_localtop_t *top,
263 t_mdatoms *md,
264 t_commrec *cr,
265 t_nrnb *nrnb,
266 real *lambda,
267 real *dvdlambda,
268 gmx_grppairener_t *grppener,
269 gmx_bool bFillGrid,
270 gmx_bool bDoLongRangeNS);
271 /* Call the neighborsearcher */
273 extern void do_force_lowlevel(FILE *fplog,
274 gmx_large_int_t step,
275 t_forcerec *fr,
276 t_inputrec *ir,
277 t_idef *idef,
278 t_commrec *cr,
279 t_nrnb *nrnb,
280 gmx_wallcycle_t wcycle,
281 t_mdatoms *md,
282 t_grpopts *opts,
283 rvec x[],
284 history_t *hist,
285 rvec f_shortrange[],
286 rvec f_longrange[],
287 gmx_enerdata_t *enerd,
288 t_fcdata *fcd,
289 gmx_mtop_t *mtop,
290 gmx_localtop_t *top,
291 gmx_genborn_t *born,
292 t_atomtypes *atype,
293 gmx_bool bBornRadii,
294 matrix box,
295 t_lambda *fepvals,
296 real *lambda,
297 t_graph *graph,
298 t_blocka *excl,
299 rvec mu_tot[2],
300 int flags,
301 float *cycles_pme);
302 /* Call all the force routines */
304 #ifdef __cplusplus
306 #endif
308 #endif /* _force_h */