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,2019, 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_MDTYPES_TYPES_FORCEREC_H
38 #define GMX_MDTYPES_TYPES_FORCEREC_H
44 #include "gromacs/math/vectypes.h"
45 #include "gromacs/mdtypes/interaction_const.h"
46 #include "gromacs/mdtypes/md_enums.h"
47 #include "gromacs/utility/basedefinitions.h"
48 #include "gromacs/utility/real.h"
50 struct ForceProviders
;
52 /* Abstract type for PME that is defined only in the routine that use them. */
55 struct nonbonded_verlet_t
;
56 struct bonded_threading_t
;
65 /* macros for the cginfo data in forcerec
67 * Since the tpx format support max 256 energy groups, we do the same here.
68 * Note that we thus have bits 8-14 still unused.
70 * The maximum cg size in cginfo is 63
71 * because we only have space for 6 bits in cginfo,
72 * this cg size entry is actually only read with domain decomposition.
74 #define SET_CGINFO_GID(cgi, gid) (cgi) = (((cgi) & ~255) | (gid))
75 #define GET_CGINFO_GID(cgi) ( (cgi) & 255)
76 #define SET_CGINFO_FEP(cgi) (cgi) = ((cgi) | (1<<15))
77 #define GET_CGINFO_FEP(cgi) ( (cgi) & (1<<15))
78 #define SET_CGINFO_EXCL_INTRA(cgi) (cgi) = ((cgi) | (1<<16))
79 #define GET_CGINFO_EXCL_INTRA(cgi) ( (cgi) & (1<<16))
80 #define SET_CGINFO_EXCL_INTER(cgi) (cgi) = ((cgi) | (1<<17))
81 #define GET_CGINFO_EXCL_INTER(cgi) ( (cgi) & (1<<17))
82 #define SET_CGINFO_SOLOPT(cgi, opt) (cgi) = (((cgi) & ~(3<<18)) | ((opt)<<18))
83 #define GET_CGINFO_SOLOPT(cgi) (((cgi)>>18) & 3)
84 #define SET_CGINFO_CONSTR(cgi) (cgi) = ((cgi) | (1<<20))
85 #define GET_CGINFO_CONSTR(cgi) ( (cgi) & (1<<20))
86 #define SET_CGINFO_SETTLE(cgi) (cgi) = ((cgi) | (1<<21))
87 #define GET_CGINFO_SETTLE(cgi) ( (cgi) & (1<<21))
88 /* This bit is only used with bBondComm in the domain decomposition */
89 #define SET_CGINFO_BOND_INTER(cgi) (cgi) = ((cgi) | (1<<22))
90 #define GET_CGINFO_BOND_INTER(cgi) ( (cgi) & (1<<22))
91 #define SET_CGINFO_HAS_VDW(cgi) (cgi) = ((cgi) | (1<<23))
92 #define GET_CGINFO_HAS_VDW(cgi) ( (cgi) & (1<<23))
93 #define SET_CGINFO_HAS_Q(cgi) (cgi) = ((cgi) | (1<<24))
94 #define GET_CGINFO_HAS_Q(cgi) ( (cgi) & (1<<24))
95 #define SET_CGINFO_NATOMS(cgi, opt) (cgi) = (((cgi) & ~(63<<25)) | ((opt)<<25))
96 #define GET_CGINFO_NATOMS(cgi) (((cgi)>>25) & 63)
99 /* Value to be used in mdrun for an infinite cut-off.
100 * Since we need to compare with the cut-off squared,
101 * this value should be slighlty smaller than sqrt(GMX_FLOAT_MAX).
103 #define GMX_CUTOFF_INF 1E+18
105 /* enums for the neighborlist type */
107 enbvdwNONE
, enbvdwLJ
, enbvdwBHAM
, enbvdwTAB
, enbvdwNR
119 /* Forward declaration of type for managing Ewald tables */
120 struct gmx_ewald_tab_t
;
122 struct ewald_corr_thread_t
;
124 struct t_forcerec
{ // NOLINT (clang-analyzer-optin.performance.Padding)
125 struct interaction_const_t
*ic
= nullptr;
129 //! Tells whether atoms inside a molecule can be in different periodic images,
130 // i.e. whether we need to take into account PBC when computing distances inside molecules.
131 // This determines whether PBC must be considered for e.g. bonded interactions.
132 gmx_bool bMolPBC
= FALSE
;
134 rvec posres_com
= { 0 };
135 rvec posres_comB
= { 0 };
137 gmx_bool use_simd_kernels
= FALSE
;
139 /* Interaction for calculated in kernels. In many cases this is similar to
140 * the electrostatics settings in the inputrecord, but the difference is that
141 * these variables always specify the actual interaction in the kernel - if
142 * we are tabulating reaction-field the inputrec will say reaction-field, but
143 * the kernel interaction will say cubic-spline-table. To be safe we also
144 * have a kernel-specific setting for the modifiers - if the interaction is
145 * tabulated we already included the inputrec modification there, so the kernel
146 * modification setting will say 'none' in that case.
148 int nbkernel_elec_interaction
= 0;
149 int nbkernel_vdw_interaction
= 0;
150 int nbkernel_elec_modifier
= 0;
151 int nbkernel_vdw_modifier
= 0;
154 * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
158 /* Parameters for generalized reaction field */
162 /* Charge sum and dipole for topology A/B ([0]/[1]) for Ewald corrections */
163 double qsum
[2] = { 0 };
164 double q2sum
[2] = { 0 };
165 double c6sum
[2] = { 0 };
166 rvec mu_tot
[2] = { { 0 } };
168 /* Dispersion correction stuff */
170 int numAtomsForDispersionCorrection
= 0;
171 struct t_forcetable
*dispersionCorrectionTable
= nullptr;
173 /* The shift of the shift or user potentials */
174 real enershiftsix
= 0;
175 real enershifttwelve
= 0;
176 /* Integrated differces for energy and virial with cut-off functions */
177 real enerdiffsix
= 0;
178 real enerdifftwelve
= 0;
180 real virdifftwelve
= 0;
181 /* Constant for long range dispersion correction (average dispersion)
182 * for topology A/B ([0]/[1]) */
183 real avcsix
[2] = { 0 };
184 /* Constant for long range repulsion term. Relative difference of about
185 * 0.1 percent with 0.8 nm cutoffs. But hey, it's cheap anyway...
187 real avctwelve
[2] = { 0 };
193 gmx_bool bcoultab
= FALSE
;
194 gmx_bool bvdwtab
= FALSE
;
196 struct t_forcetable
*pairsTable
; /* for 1-4 interactions, [pairs] and [pairs_nb] */
200 real sc_alphavdw
= 0;
201 real sc_alphacoul
= 0;
204 real sc_sigma6_def
= 0;
205 real sc_sigma6_min
= 0;
210 /* solvent_opt contains the enum for the most common solvent
211 * in the system, which will be optimized.
212 * It can be set to esolNO to disable all water optimization */
215 gmx_bool bGrid
= FALSE
;
216 gmx_bool bExcl_IntraCGAll_InterCGNone
= FALSE
;
217 struct cginfo_mb_t
*cginfo_mb
= nullptr;
218 int *cginfo
= nullptr;
219 rvec
*cg_cm
= nullptr;
221 rvec
*shift_vec
= nullptr;
223 int cutoff_scheme
= 0; /* group- or Verlet-style cutoff */
224 gmx_bool bNonbonded
= FALSE
; /* true if nonbonded calculations are *not* turned off */
226 /* The Nbnxm Verlet non-bonded machinery */
227 std::unique_ptr
<nonbonded_verlet_t
> nbv
;
229 /* The wall tables (if used) */
231 struct t_forcetable
***wall_tab
= nullptr;
233 /* The number of charge groups participating in do_force_lowlevel */
235 /* The number of atoms participating in do_force_lowlevel */
236 int natoms_force
= 0;
237 /* The number of atoms participating in force and constraints */
238 int natoms_force_constr
= 0;
239 /* The allocation size of vectors of size natoms_force */
240 int nalloc_force
= 0;
242 /* Forces that should not enter into the coord x force virial summation:
243 * PPPM/PME/Ewald/posres/ForceProviders
245 /* True when we have contributions that are directly added to the virial */
246 gmx_bool haveDirectVirialContributions
= FALSE
;
247 /* TODO: Replace the pointer by an object once we got rid of C */
248 std::vector
<gmx::RVec
> *forceBufferForDirectVirialContributions
= nullptr;
250 /* Data for PPPM/PME/Ewald */
251 struct gmx_pme_t
*pmedata
= nullptr;
252 int ljpme_combination_rule
= 0;
254 /* PME/Ewald stuff */
255 struct gmx_ewald_tab_t
*ewald_table
= nullptr;
257 /* Shift force array for computing the virial */
258 rvec
*fshift
= nullptr;
260 /* Non bonded Parameter lists */
261 int ntype
= 0; /* Number of atom types */
262 gmx_bool bBHAM
= FALSE
;
263 real
*nbfp
= nullptr;
264 real
*ljpme_c6grid
= nullptr; /* C6-values used on grid in LJPME */
266 /* Energy group pair flags */
267 int *egp_flags
= nullptr;
269 /* Shell molecular dynamics flexible constraints */
270 real fc_stepsize
= 0;
272 /* If > 0 signals Test Particle Insertion,
273 * the value is the number of atoms of the molecule to insert
274 * Only the energy difference due to the addition of the last molecule
275 * should be calculated.
280 gmx_bool bQMMM
= FALSE
;
281 struct t_QMMMrec
*qr
= nullptr;
283 /* QM-MM neighborlists */
284 struct t_nblist
*QMMMlist
= nullptr;
286 /* Limit for printing large forces, negative is don't print */
287 real print_force
= 0;
289 /* User determined parameters, copied from the inputrec */
299 /* Pointer to struct for managing threading of bonded force calculation */
300 struct bonded_threading_t
*bondedThreading
= nullptr;
302 /* TODO: Replace the pointer by an object once we got rid of C */
303 gmx::GpuBonded
*gpuBonded
= nullptr;
305 /* Ewald correction thread local virial and energy data */
307 struct ewald_corr_thread_t
*ewc_t
= nullptr;
309 struct ForceProviders
*forceProviders
= nullptr;
312 /* Important: Starting with Gromacs-4.6, the values of c6 and c12 in the nbfp array have
313 * been scaled by 6.0 or 12.0 to save flops in the kernels. We have corrected this everywhere
314 * in the code, but beware if you are using these macros externally.
316 #define C6(nbfp, ntp, ai, aj) (nbfp)[2*((ntp)*(ai)+(aj))]
317 #define C12(nbfp, ntp, ai, aj) (nbfp)[2*((ntp)*(ai)+(aj))+1]
318 #define BHAMC(nbfp, ntp, ai, aj) (nbfp)[3*((ntp)*(ai)+(aj))]
319 #define BHAMA(nbfp, ntp, ai, aj) (nbfp)[3*((ntp)*(ai)+(aj))+1]
320 #define BHAMB(nbfp, ntp, ai, aj) (nbfp)[3*((ntp)*(ai)+(aj))+2]