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 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 #ifndef GMX_MDTYPES_TYPES_FORCEREC_H
39 #define GMX_MDTYPES_TYPES_FORCEREC_H
45 #include "gromacs/math/vectypes.h"
46 #include "gromacs/mdtypes/md_enums.h"
47 #include "gromacs/pbcutil/pbc.h"
48 #include "gromacs/utility/basedefinitions.h"
49 #include "gromacs/utility/real.h"
51 /* Abstract type for PME that is defined only in the routine that use them. */
53 struct nonbonded_verlet_t;
54 struct bonded_threading_t;
56 class DispersionCorrection;
62 class DeviceStreamManager;
65 class StatePropagatorDataGpu;
67 class WholeMoleculeTransform;
70 /* macros for the cginfo data in forcerec
72 * Since the tpx format support max 256 energy groups, we do the same here.
73 * Note that we thus have bits 8-14 still unused.
75 * The maximum cg size in cginfo is 63
76 * because we only have space for 6 bits in cginfo,
77 * this cg size entry is actually only read with domain decomposition.
79 #define SET_CGINFO_GID(cgi, gid) (cgi) = (((cgi) & ~255) | (gid))
80 #define GET_CGINFO_GID(cgi) ((cgi)&255)
81 #define SET_CGINFO_FEP(cgi) (cgi) = ((cgi) | (1 << 15))
82 #define GET_CGINFO_FEP(cgi) ((cgi) & (1 << 15))
83 #define SET_CGINFO_EXCL_INTER(cgi) (cgi) = ((cgi) | (1 << 17))
84 #define GET_CGINFO_EXCL_INTER(cgi) ((cgi) & (1 << 17))
85 #define SET_CGINFO_CONSTR(cgi) (cgi) = ((cgi) | (1 << 20))
86 #define GET_CGINFO_CONSTR(cgi) ((cgi) & (1 << 20))
87 #define SET_CGINFO_SETTLE(cgi) (cgi) = ((cgi) | (1 << 21))
88 #define GET_CGINFO_SETTLE(cgi) ((cgi) & (1 << 21))
89 /* This bit is only used with bBondComm in the domain decomposition */
90 #define SET_CGINFO_BOND_INTER(cgi) (cgi) = ((cgi) | (1 << 22))
91 #define GET_CGINFO_BOND_INTER(cgi) ((cgi) & (1 << 22))
92 #define SET_CGINFO_HAS_VDW(cgi) (cgi) = ((cgi) | (1 << 23))
93 #define GET_CGINFO_HAS_VDW(cgi) ((cgi) & (1 << 23))
94 #define SET_CGINFO_HAS_Q(cgi) (cgi) = ((cgi) | (1 << 24))
95 #define GET_CGINFO_HAS_Q(cgi) ((cgi) & (1 << 24))
98 /* Value to be used in mdrun for an infinite cut-off.
99 * Since we need to compare with the cut-off squared,
100 * this value should be slighlty smaller than sqrt(GMX_FLOAT_MAX).
102 #define GMX_CUTOFF_INF 1E+18
104 /* enums for the neighborlist type */
119 std::vector<int> cginfo;
123 /* Forward declaration of type for managing Ewald tables */
124 struct gmx_ewald_tab_t;
126 struct ewald_corr_thread_t;
129 { // NOLINT (clang-analyzer-optin.performance.Padding)
130 // Declare an explicit constructor and destructor, so they can be
131 // implemented in a single source file, so that not every source
132 // file that includes this one needs to understand how to find the
133 // destructors of the objects pointed to by unique_ptr members.
137 struct interaction_const_t* ic = nullptr;
140 PbcType pbcType = PbcType::Xyz;
141 //! Tells whether atoms inside a molecule can be in different periodic images,
142 // i.e. whether we need to take into account PBC when computing distances inside molecules.
143 // This determines whether PBC must be considered for e.g. bonded interactions.
144 gmx_bool bMolPBC = FALSE;
146 rvec posres_com = { 0 };
147 rvec posres_comB = { 0 };
149 gmx_bool use_simd_kernels = FALSE;
151 /* Interaction for calculated in kernels. In many cases this is similar to
152 * the electrostatics settings in the inputrecord, but the difference is that
153 * these variables always specify the actual interaction in the kernel - if
154 * we are tabulating reaction-field the inputrec will say reaction-field, but
155 * the kernel interaction will say cubic-spline-table. To be safe we also
156 * have a kernel-specific setting for the modifiers - if the interaction is
157 * tabulated we already included the inputrec modification there, so the kernel
158 * modification setting will say 'none' in that case.
160 int nbkernel_elec_interaction = 0;
161 int nbkernel_vdw_interaction = 0;
162 int nbkernel_elec_modifier = 0;
163 int nbkernel_vdw_modifier = 0;
166 * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
170 /* Charge sum for topology A/B ([0]/[1]) for Ewald corrections */
171 double qsum[2] = { 0 };
172 double q2sum[2] = { 0 };
173 double c6sum[2] = { 0 };
175 /* Dispersion correction stuff */
176 std::unique_ptr<DispersionCorrection> dispersionCorrection;
182 gmx_bool bcoultab = FALSE;
183 gmx_bool bvdwtab = FALSE;
185 t_forcetable* pairsTable = nullptr; /* for 1-4 interactions, [pairs] and [pairs_nb] */
189 real sc_alphavdw = 0;
190 real sc_alphacoul = 0;
193 real sc_sigma6_def = 0;
194 real sc_sigma6_min = 0;
196 /* Information about atom properties for the molecule blocks in the system */
197 std::vector<cginfo_mb_t> cginfo_mb;
198 /* Information about atom properties for local and non-local atoms */
199 std::vector<int> cginfo;
201 rvec* shift_vec = nullptr;
203 std::unique_ptr<gmx::WholeMoleculeTransform> wholeMoleculeTransform;
205 int cutoff_scheme = 0; /* group- or Verlet-style cutoff */
206 gmx_bool bNonbonded = FALSE; /* true if nonbonded calculations are *not* turned off */
208 /* The Nbnxm Verlet non-bonded machinery */
209 std::unique_ptr<nonbonded_verlet_t> nbv;
211 /* The wall tables (if used) */
213 t_forcetable*** wall_tab = nullptr;
215 /* The number of atoms participating in do_force_lowlevel */
216 int natoms_force = 0;
217 /* The number of atoms participating in force and constraints */
218 int natoms_force_constr = 0;
219 /* The allocation size of vectors of size natoms_force */
220 int nalloc_force = 0;
222 /* Forces that should not enter into the coord x force virial summation:
223 * PPPM/PME/Ewald/posres/ForceProviders
225 /* True when we have contributions that are directly added to the virial */
226 bool haveDirectVirialContributions = false;
227 /* Force buffer for force computation with direct virial contributions */
228 std::vector<gmx::RVec> forceBufferForDirectVirialContributions;
230 /* Data for PPPM/PME/Ewald */
231 struct gmx_pme_t* pmedata = nullptr;
232 int ljpme_combination_rule = 0;
234 /* PME/Ewald stuff */
235 struct gmx_ewald_tab_t* ewald_table = nullptr;
237 /* Shift force array for computing the virial, size SHIFTS */
238 std::vector<gmx::RVec> shiftForces;
240 /* Non bonded Parameter lists */
241 int ntype = 0; /* Number of atom types */
242 gmx_bool bBHAM = FALSE;
243 std::vector<real> nbfp;
244 real* ljpme_c6grid = nullptr; /* C6-values used on grid in LJPME */
246 /* Energy group pair flags */
247 int* egp_flags = nullptr;
249 /* Shell molecular dynamics flexible constraints */
250 real fc_stepsize = 0;
252 /* If > 0 signals Test Particle Insertion,
253 * the value is the number of atoms of the molecule to insert
254 * Only the energy difference due to the addition of the last molecule
255 * should be calculated.
259 /* Limit for printing large forces, negative is don't print */
260 real print_force = 0;
262 /* User determined parameters, copied from the inputrec */
272 /* Pointer to struct for managing threading of bonded force calculation */
273 struct bonded_threading_t* bondedThreading = nullptr;
275 /* TODO: Replace the pointer by an object once we got rid of C */
276 gmx::GpuBonded* gpuBonded = nullptr;
278 /* Ewald correction thread local virial and energy data */
280 struct ewald_corr_thread_t* ewc_t = nullptr;
282 gmx::ForceProviders* forceProviders = nullptr;
284 // The stateGpu object is created in runner, forcerec just keeps the copy of the pointer.
285 // TODO: This is not supposed to be here. StatePropagatorDataGpu should be a part of
286 // general StatePropagatorData object that is passed around
287 gmx::StatePropagatorDataGpu* stateGpu = nullptr;
288 // TODO: Should not be here. This is here only to pass the pointer around.
289 gmx::DeviceStreamManager* deviceStreamManager = nullptr;
291 //! GPU device context
292 DeviceContext* deviceContext = nullptr;
294 /* For PME-PP GPU communication */
295 std::unique_ptr<gmx::PmePpCommGpu> pmePpCommGpu;
298 /* Important: Starting with Gromacs-4.6, the values of c6 and c12 in the nbfp array have
299 * been scaled by 6.0 or 12.0 to save flops in the kernels. We have corrected this everywhere
300 * in the code, but beware if you are using these macros externally.
302 #define C6(nbfp, ntp, ai, aj) (nbfp)[2 * ((ntp) * (ai) + (aj))]
303 #define C12(nbfp, ntp, ai, aj) (nbfp)[2 * ((ntp) * (ai) + (aj)) + 1]
304 #define BHAMC(nbfp, ntp, ai, aj) (nbfp)[3 * ((ntp) * (ai) + (aj))]
305 #define BHAMA(nbfp, ntp, ai, aj) (nbfp)[3 * ((ntp) * (ai) + (aj)) + 1]
306 #define BHAMB(nbfp, ntp, ai, aj) (nbfp)[3 * ((ntp) * (ai) + (aj)) + 2]