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