Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / legacyheaders / force.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37
38 #ifndef _force_h
39 #define _force_h
40
41
42 #include "typedefs.h"
43 #include "types/force_flags.h"
44 #include "network.h"
45 #include "tgroup.h"
46 #include "vsite.h"
47 #include "genborn.h"
48
49 #include "../timing/wallcycle.h"
50
51 #ifdef __cplusplus
52 extern "C" {
53 #endif
54
55 struct t_graph;
56 struct t_pbc;
57
58 void gmx_print_sepdvdl(FILE *fplog, const char *s, real v, real dvdlambda);
59
60 void calc_vir(int nxf, rvec x[], rvec f[], tensor vir,
61               gmx_bool bScrewPBC, matrix box);
62 /* Calculate virial for nxf atoms, and add it to vir */
63
64 void f_calc_vir(int i0, int i1, rvec x[], rvec f[], tensor vir,
65                 struct t_graph *g, rvec shift_vec[]);
66 /* Calculate virial taking periodicity into account */
67
68 real RF_excl_correction(const t_forcerec *fr, struct t_graph *g,
69                         const t_mdatoms *mdatoms, const t_blocka *excl,
70                         rvec x[], rvec f[], rvec *fshift, const struct t_pbc *pbc,
71                         real lambda, real *dvdlambda);
72 /* Calculate the reaction-field energy correction for this node:
73  * epsfac q_i q_j (k_rf r_ij^2 - c_rf)
74  * and force correction for all excluded pairs, including self pairs.
75  */
76
77 void calc_rffac(FILE *fplog, int eel, real eps_r, real eps_rf,
78                 real Rc, real Temp,
79                 real zsq, matrix box,
80                 real *kappa, real *krf, real *crf);
81 /* Determine the reaction-field constants */
82
83 void init_generalized_rf(FILE *fplog,
84                          const gmx_mtop_t *mtop, const t_inputrec *ir,
85                          t_forcerec *fr);
86 /* Initialize the generalized reaction field parameters */
87
88
89 /* In wall.c */
90 void make_wall_tables(FILE *fplog, const output_env_t oenv,
91                       const t_inputrec *ir, const char *tabfn,
92                       const gmx_groups_t *groups,
93                       t_forcerec *fr);
94
95 real do_walls(t_inputrec *ir, t_forcerec *fr, matrix box, t_mdatoms *md,
96               rvec x[], rvec f[], real lambda, real Vlj[], t_nrnb *nrnb);
97
98 t_forcerec *mk_forcerec(void);
99
100 #define GMX_MAKETABLES_FORCEUSER  (1<<0)
101 #define GMX_MAKETABLES_14ONLY     (1<<1)
102
103 t_forcetable make_tables(FILE *fp, const output_env_t oenv,
104                          const t_forcerec *fr, gmx_bool bVerbose,
105                          const char *fn, real rtab, int flags);
106 /* Return tables for inner loops. When bVerbose the tables are printed
107  * to .xvg files
108  */
109
110 bondedtable_t make_bonded_table(FILE *fplog, char *fn, int angle);
111 /* Return a table for bonded interactions,
112  * angle should be: bonds 0, angles 1, dihedrals 2
113  */
114
115 /* Return a table for GB calculations */
116 t_forcetable make_gb_table(const output_env_t oenv,
117                            const t_forcerec  *fr);
118
119 /* Read a table for AdResS Thermo Force calculations */
120 extern t_forcetable make_atf_table(FILE *out, const output_env_t oenv,
121                                    const t_forcerec *fr,
122                                    const char *fn,
123                                    matrix box);
124
125 void pr_forcerec(FILE *fplog, t_forcerec *fr);
126
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 */
133
134 gmx_bool can_use_allvsall(const t_inputrec *ir,
135                           gmx_bool bPrintNote, t_commrec *cr, FILE *fp);
136 /* Returns if we can use all-vs-all loops.
137  * If bPrintNote==TRUE, prints a note, if necessary, to stderr
138  * and fp (if !=NULL) on the master node.
139  */
140
141
142 gmx_bool nbnxn_acceleration_supported(FILE             *fplog,
143                                       const t_commrec  *cr,
144                                       const t_inputrec *ir,
145                                       gmx_bool          bGPU);
146 /* Return if GPU/CPU-SIMD acceleration is supported with the given inputrec
147  * with bGPU TRUE/FALSE.
148  * If the return value is FALSE and fplog/cr != NULL, prints a fallback
149  * message to fplog/stderr.
150  */
151
152 gmx_bool uses_simple_tables(int                 cutoff_scheme,
153                             nonbonded_verlet_t *nbv,
154                             int                 group);
155 /* Returns whether simple tables (i.e. not for use with GPUs) are used
156  * with the type of kernel indicated.
157  */
158
159 void init_interaction_const_tables(FILE                *fp,
160                                    interaction_const_t *ic,
161                                    gmx_bool             bSimpleTable,
162                                    real                 rtab);
163 /* Initializes the tables in the interaction constant data structure.
164  * Setting verlet_kernel_type to -1 always initializes tables for
165  * use with group kernels.
166  */
167
168 void init_forcerec(FILE              *fplog,
169                    const output_env_t oenv,
170                    t_forcerec        *fr,
171                    t_fcdata          *fcd,
172                    const t_inputrec  *ir,
173                    const gmx_mtop_t  *mtop,
174                    const t_commrec   *cr,
175                    matrix             box,
176                    const char        *tabfn,
177                    const char        *tabafn,
178                    const char        *tabpfn,
179                    const char        *tabbfn,
180                    const char        *nbpu_opt,
181                    gmx_bool           bNoSolvOpt,
182                    real               print_force);
183 /* The Force rec struct must be created with mk_forcerec
184  * The gmx_booleans have the following meaning:
185  * bSetQ:    Copy the charges [ only necessary when they change ]
186  * bMolEpot: Use the free energy stuff per molecule
187  * print_force >= 0: print forces for atoms with force >= print_force
188  */
189
190 void forcerec_set_excl_load(t_forcerec           *fr,
191                             const gmx_localtop_t *top);
192 /* Set the exclusion load for the local exclusions and possibly threads */
193
194 void init_enerdata(int ngener, int n_lambda, gmx_enerdata_t *enerd);
195 /* Intializes the energy storage struct */
196
197 void destroy_enerdata(gmx_enerdata_t *enerd);
198 /* Free all memory associated with enerd */
199
200 void reset_foreign_enerdata(gmx_enerdata_t *enerd);
201 /* Resets only the foreign energy data */
202
203 void reset_enerdata(t_forcerec *fr, gmx_bool bNS,
204                     gmx_enerdata_t *enerd,
205                     gmx_bool bMaster);
206 /* Resets the energy data, if bNS=TRUE also zeros the long-range part */
207
208 void sum_epot(gmx_grppairener_t *grpp, real *epot);
209 /* Locally sum the non-bonded potential energy terms */
210
211 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals);
212 /* Sum the free energy contributions */
213
214 void update_forcerec(t_forcerec *fr, matrix box);
215 /* Updates parameters in the forcerec that are time dependent */
216
217 /* Compute the average C6 and C12 params for LJ corrections */
218 void set_avcsixtwelve(FILE *fplog, t_forcerec *fr,
219                       const gmx_mtop_t *mtop);
220
221 extern void do_force(FILE *log, t_commrec *cr,
222                      t_inputrec *inputrec,
223                      gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
224                      gmx_localtop_t *top,
225                      gmx_groups_t *groups,
226                      matrix box, rvec x[], history_t *hist,
227                      rvec f[],
228                      tensor vir_force,
229                      t_mdatoms *mdatoms,
230                      gmx_enerdata_t *enerd, t_fcdata *fcd,
231                      real *lambda, struct t_graph *graph,
232                      t_forcerec *fr,
233                      gmx_vsite_t *vsite, rvec mu_tot,
234                      double t, FILE *field, gmx_edsam_t ed,
235                      gmx_bool bBornRadii,
236                      int flags);
237
238 /* Communicate coordinates (if parallel).
239  * Do neighbor searching (if necessary).
240  * Calculate forces.
241  * Communicate forces (if parallel).
242  * Spread forces for vsites (if present).
243  *
244  * f is always required.
245  */
246
247 void ns(FILE              *fplog,
248         t_forcerec        *fr,
249         matrix             box,
250         gmx_groups_t      *groups,
251         gmx_localtop_t    *top,
252         t_mdatoms         *md,
253         t_commrec         *cr,
254         t_nrnb            *nrnb,
255         gmx_bool           bFillGrid,
256         gmx_bool           bDoLongRangeNS);
257 /* Call the neighborsearcher */
258
259 extern void do_force_lowlevel(FILE         *fplog,
260                               gmx_int64_t   step,
261                               t_forcerec   *fr,
262                               t_inputrec   *ir,
263                               t_idef       *idef,
264                               t_commrec    *cr,
265                               t_nrnb       *nrnb,
266                               gmx_wallcycle_t wcycle,
267                               t_mdatoms    *md,
268                               rvec         x[],
269                               history_t    *hist,
270                               rvec         f_shortrange[],
271                               rvec         f_longrange[],
272                               gmx_enerdata_t *enerd,
273                               t_fcdata     *fcd,
274                               gmx_localtop_t *top,
275                               gmx_genborn_t *born,
276                               t_atomtypes  *atype,
277                               gmx_bool         bBornRadii,
278                               matrix       box,
279                               t_lambda     *fepvals,
280                               real         *lambda,
281                               struct t_graph      *graph,
282                               t_blocka     *excl,
283                               rvec         mu_tot[2],
284                               int          flags,
285                               float        *cycles_pme);
286 /* Call all the force routines */
287
288 #ifdef __cplusplus
289 }
290 #endif
291
292 #endif  /* _force_h */