Propose a post-submit matrix
[gromacs.git] / src / gromacs / mdlib / force.h
blob3fe505cf0b4f3397ac061371fdc115d1a1b5d8fb
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015, 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.
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.
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.
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.
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.
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.
37 #ifndef GMX_MDLIB_FORCE_H
38 #define GMX_MDLIB_FORCE_H
40 #include "gromacs/mdlib/force_flags.h"
41 #include "gromacs/mdlib/vsite.h"
42 #include "gromacs/mdtypes/fcdata.h"
43 #include "gromacs/mdtypes/forcerec.h"
44 #include "gromacs/timing/wallcycle.h"
46 struct gmx_edsam;
47 struct gmx_gpu_info_t;
48 struct gmx_groups_t;
49 struct gmx_vsite_t;
50 struct history_t;
51 struct nonbonded_verlet_t;
52 struct t_blocka;
53 struct t_commrec;
54 struct t_fcdata;
55 struct t_graph;
56 struct t_grpopts;
57 struct t_inputrec;
58 struct t_lambda;
59 struct t_mdatoms;
60 struct t_nrnb;
61 struct t_pbc;
63 void calc_vir(int nxf, rvec x[], rvec f[], tensor vir,
64 gmx_bool bScrewPBC, matrix box);
65 /* Calculate virial for nxf atoms, and add it to vir */
67 void f_calc_vir(int i0, int i1, rvec x[], rvec f[], tensor vir,
68 t_graph *g, rvec shift_vec[]);
69 /* Calculate virial taking periodicity into account */
71 real RF_excl_correction(const t_forcerec *fr, t_graph *g,
72 const t_mdatoms *mdatoms, const t_blocka *excl,
73 rvec x[], rvec f[], rvec *fshift, const t_pbc *pbc,
74 real lambda, real *dvdlambda);
75 /* Calculate the reaction-field energy correction for this node:
76 * epsfac q_i q_j (k_rf r_ij^2 - c_rf)
77 * and force correction for all excluded pairs, including self pairs.
80 void calc_rffac(FILE *fplog, int eel, real eps_r, real eps_rf,
81 real Rc, real Temp,
82 real zsq, matrix box,
83 real *kappa, real *krf, real *crf);
84 /* Determine the reaction-field constants */
86 void init_generalized_rf(FILE *fplog,
87 const gmx_mtop_t *mtop, const t_inputrec *ir,
88 t_forcerec *fr);
89 /* Initialize the generalized reaction field parameters */
92 /* In wall.c */
93 void make_wall_tables(FILE *fplog,
94 const t_inputrec *ir, const char *tabfn,
95 const gmx_groups_t *groups,
96 t_forcerec *fr);
98 real do_walls(t_inputrec *ir, t_forcerec *fr, matrix box, t_mdatoms *md,
99 rvec x[], rvec f[], real lambda, real Vlj[], t_nrnb *nrnb);
101 gmx_bool can_use_allvsall(const t_inputrec *ir,
102 gmx_bool bPrintNote, t_commrec *cr, FILE *fp);
103 /* Returns if we can use all-vs-all loops.
104 * If bPrintNote==TRUE, prints a note, if necessary, to stderr
105 * and fp (if !=NULL) on the master node.
109 gmx_bool nbnxn_gpu_acceleration_supported(FILE *fplog,
110 const t_commrec *cr,
111 const t_inputrec *ir,
112 gmx_bool bRerunMD);
113 /* Return if GPU acceleration is supported with the given settings.
115 * If the return value is FALSE and fplog/cr != NULL, prints a fallback
116 * message to fplog/stderr.
119 gmx_bool nbnxn_simd_supported(FILE *fplog,
120 const t_commrec *cr,
121 const t_inputrec *ir);
122 /* Return if CPU SIMD support exists for the given inputrec
123 * If the return value is FALSE and fplog/cr != NULL, prints a fallback
124 * message to fplog/stderr.
127 gmx_bool uses_simple_tables(int cutoff_scheme,
128 nonbonded_verlet_t *nbv,
129 int group);
130 /* Returns whether simple tables (i.e. not for use with GPUs) are used
131 * with the type of kernel indicated.
134 void init_enerdata(int ngener, int n_lambda, gmx_enerdata_t *enerd);
135 /* Intializes the energy storage struct */
137 void destroy_enerdata(gmx_enerdata_t *enerd);
138 /* Free all memory associated with enerd */
140 void reset_foreign_enerdata(gmx_enerdata_t *enerd);
141 /* Resets only the foreign energy data */
143 void reset_enerdata(gmx_enerdata_t *enerd);
144 /* Resets the energy data */
146 void sum_epot(gmx_grppairener_t *grpp, real *epot);
147 /* Locally sum the non-bonded potential energy terms */
149 void sum_dhdl(gmx_enerdata_t *enerd, real *lambda, t_lambda *fepvals);
150 /* Sum the free energy contributions */
152 /* Compute the average C6 and C12 params for LJ corrections */
153 void set_avcsixtwelve(FILE *fplog, t_forcerec *fr,
154 const gmx_mtop_t *mtop);
156 void do_force(FILE *log, t_commrec *cr,
157 t_inputrec *inputrec,
158 gmx_int64_t step, struct t_nrnb *nrnb, gmx_wallcycle_t wcycle,
159 gmx_localtop_t *top,
160 gmx_groups_t *groups,
161 matrix box, rvec x[], history_t *hist,
162 rvec f[],
163 tensor vir_force,
164 t_mdatoms *mdatoms,
165 gmx_enerdata_t *enerd, t_fcdata *fcd,
166 real *lambda, struct t_graph *graph,
167 t_forcerec *fr,
168 gmx_vsite_t *vsite, rvec mu_tot,
169 double t, FILE *field, struct gmx_edsam *ed,
170 gmx_bool bBornRadii,
171 int flags);
173 /* Communicate coordinates (if parallel).
174 * Do neighbor searching (if necessary).
175 * Calculate forces.
176 * Communicate forces (if parallel).
177 * Spread forces for vsites (if present).
179 * f is always required.
182 void ns(FILE *fplog,
183 t_forcerec *fr,
184 matrix box,
185 gmx_groups_t *groups,
186 gmx_localtop_t *top,
187 t_mdatoms *md,
188 t_commrec *cr,
189 t_nrnb *nrnb,
190 gmx_bool bFillGrid);
191 /* Call the neighborsearcher */
193 void do_force_lowlevel(t_forcerec *fr,
194 t_inputrec *ir,
195 t_idef *idef,
196 t_commrec *cr,
197 t_nrnb *nrnb,
198 gmx_wallcycle_t wcycle,
199 t_mdatoms *md,
200 rvec x[],
201 history_t *hist,
202 rvec f_shortrange[],
203 gmx_enerdata_t *enerd,
204 t_fcdata *fcd,
205 gmx_localtop_t *top,
206 gmx_genborn_t *born,
207 gmx_bool bBornRadii,
208 matrix box,
209 t_lambda *fepvals,
210 real *lambda,
211 t_graph *graph,
212 t_blocka *excl,
213 rvec mu_tot[2],
214 int flags,
215 float *cycles_pme);
216 /* Call all the force routines */
218 void free_gpu_resources(const t_forcerec *fr,
219 const t_commrec *cr,
220 const gmx_gpu_info_t *gpu_info,
221 const gmx_gpu_opt_t *gpu_opt);
223 #endif