More const correctness
[gromacs.git] / src / gromacs / mdlib / force.h
blob34075da9ce248c390fbc195a53a3483b61348db3
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,2016,2017,2018, 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/domdec/dlbtiming.h"
41 #include "gromacs/mdlib/force_flags.h"
42 #include "gromacs/mdlib/vsite.h"
43 #include "gromacs/mdtypes/fcdata.h"
44 #include "gromacs/mdtypes/forcerec.h"
45 #include "gromacs/timing/wallcycle.h"
46 #include "gromacs/utility/arrayref.h"
48 struct gmx_device_info_t;
49 struct gmx_edsam;
50 struct gmx_gpu_info_t;
51 struct gmx_groups_t;
52 struct gmx_multisim_t;
53 struct gmx_vsite_t;
54 class history_t;
55 struct nonbonded_verlet_t;
56 struct t_blocka;
57 struct t_commrec;
58 struct t_fcdata;
59 struct t_graph;
60 struct t_grpopts;
61 struct t_inputrec;
62 struct t_lambda;
63 struct t_mdatoms;
64 struct t_nrnb;
65 struct t_pbc;
67 namespace gmx
69 class ForceWithVirial;
70 class MDLogger;
71 class PhysicalNodeCommunicator;
74 void calc_vir(int nxf, rvec x[], rvec f[], tensor vir,
75 gmx_bool bScrewPBC, matrix box);
76 /* Calculate virial for nxf atoms, and add it to vir */
78 void f_calc_vir(int i0, int i1, rvec x[], rvec f[], tensor vir,
79 t_graph *g, rvec shift_vec[]);
80 /* Calculate virial taking periodicity into account */
82 real RF_excl_correction(const t_forcerec *fr,
83 const t_graph *g,
84 const t_mdatoms *mdatoms,
85 const t_blocka *excl,
86 rvec x[],
87 rvec f[],
88 rvec *fshift,
89 const t_pbc *pbc,
90 real lambda,
91 real *dvdlambda);
92 /* Calculate the reaction-field energy correction for this node:
93 * epsfac q_i q_j (k_rf r_ij^2 - c_rf)
94 * and force correction for all excluded pairs, including self pairs.
97 void calc_rffac(FILE *fplog, int eel, real eps_r, real eps_rf,
98 real Rc, real Temp,
99 real zsq, matrix box,
100 real *krf, real *crf);
101 /* Determine the reaction-field constants */
103 void init_generalized_rf(FILE *fplog,
104 const gmx_mtop_t *mtop, const t_inputrec *ir,
105 t_forcerec *fr);
106 /* Initialize the generalized reaction field parameters */
109 /* In wall.c */
110 void make_wall_tables(FILE *fplog,
111 const t_inputrec *ir, const char *tabfn,
112 const gmx_groups_t *groups,
113 t_forcerec *fr);
115 real do_walls(const t_inputrec *ir,
116 t_forcerec *fr,
117 matrix box,
118 const t_mdatoms *md,
119 const rvec x[],
120 rvec f[],
121 real lambda,
122 real Vlj[],
123 t_nrnb *nrnb);
125 gmx_bool can_use_allvsall(const t_inputrec *ir,
126 gmx_bool bPrintNote, const t_commrec *cr, FILE *fp);
127 /* Returns if we can use all-vs-all loops.
128 * If bPrintNote==TRUE, prints a note, if necessary, to stderr
129 * and fp (if !=NULL) on the master node.
132 gmx_bool nbnxn_simd_supported(const gmx::MDLogger &mdlog,
133 const t_inputrec *ir);
134 /* Return if CPU SIMD support exists for the given inputrec
135 * If the return value is FALSE and fplog/cr != NULL, prints a fallback
136 * message to fplog/stderr.
139 gmx_bool uses_simple_tables(int cutoff_scheme,
140 nonbonded_verlet_t *nbv,
141 int group);
142 /* Returns whether simple tables (i.e. not for use with GPUs) are used
143 * with the type of kernel indicated.
146 void init_enerdata(int ngener, int n_lambda, gmx_enerdata_t *enerd);
147 /* Intializes the energy storage struct */
149 void destroy_enerdata(gmx_enerdata_t *enerd);
150 /* Free all memory associated with enerd */
152 void reset_foreign_enerdata(gmx_enerdata_t *enerd);
153 /* Resets only the foreign energy data */
155 void reset_enerdata(gmx_enerdata_t *enerd);
156 /* Resets the energy data */
158 void sum_epot(gmx_grppairener_t *grpp, real *epot);
159 /* Locally sum the non-bonded potential energy terms */
161 void sum_dhdl(gmx_enerdata_t *enerd, gmx::ArrayRef<const real> lambda, t_lambda *fepvals);
162 /* Sum the free energy contributions */
164 /* Compute the average C6 and C12 params for LJ corrections */
165 void set_avcsixtwelve(FILE *fplog, t_forcerec *fr,
166 const gmx_mtop_t *mtop);
168 void do_force(FILE *log,
169 const t_commrec *cr,
170 const gmx_multisim_t *ms,
171 const t_inputrec *inputrec,
172 gmx_int64_t step,
173 t_nrnb *nrnb,
174 gmx_wallcycle_t wcycle,
175 // TODO top can be const when the group scheme no longer
176 // builds exclusions during neighbor searching within
177 // do_force_cutsGROUP.
178 gmx_localtop_t *top,
179 const gmx_groups_t *groups,
180 matrix box,
181 gmx::PaddedArrayRef<gmx::RVec> coordinates,
182 history_t *hist,
183 gmx::PaddedArrayRef<gmx::RVec> force,
184 tensor vir_force,
185 const t_mdatoms *mdatoms,
186 gmx_enerdata_t *enerd,
187 t_fcdata *fcd,
188 gmx::ArrayRef<real> lambda,
189 t_graph *graph,
190 t_forcerec *fr,
191 const gmx_vsite_t *vsite,
192 rvec mu_tot,
193 double t,
194 const gmx_edsam *ed,
195 int flags,
196 DdOpenBalanceRegionBeforeForceComputation ddOpenBalanceRegion,
197 DdCloseBalanceRegionAfterForceComputation ddCloseBalanceRegion);
199 /* Communicate coordinates (if parallel).
200 * Do neighbor searching (if necessary).
201 * Calculate forces.
202 * Communicate forces (if parallel).
203 * Spread forces for vsites (if present).
205 * f is always required.
208 void ns(FILE *fplog,
209 t_forcerec *fr,
210 matrix box,
211 const gmx_groups_t *groups,
212 gmx_localtop_t *top,
213 const t_mdatoms *md,
214 const t_commrec *cr,
215 t_nrnb *nrnb,
216 gmx_bool bFillGrid);
217 /* Call the neighborsearcher */
219 void do_force_lowlevel(t_forcerec *fr,
220 const t_inputrec *ir,
221 const t_idef *idef,
222 const t_commrec *cr,
223 const gmx_multisim_t *ms,
224 t_nrnb *nrnb,
225 gmx_wallcycle_t wcycle,
226 const t_mdatoms *md,
227 rvec x[],
228 history_t *hist,
229 rvec f_shortrange[],
230 gmx::ForceWithVirial *forceWithVirial,
231 gmx_enerdata_t *enerd,
232 t_fcdata *fcd,
233 matrix box,
234 t_lambda *fepvals,
235 real *lambda,
236 const t_graph *graph,
237 const t_blocka *excl,
238 rvec mu_tot[2],
239 int flags,
240 float *cycles_pme);
241 /* Call all the force routines */
243 void free_gpu_resources(const t_forcerec *fr,
244 const gmx::PhysicalNodeCommunicator &physicalNodeCommunicator);
246 #endif