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/arrayref.h"
49 #include "gromacs/utility/basedefinitions.h"
50 #include "gromacs/utility/real.h"
52 /* Abstract type for PME that is defined only in the routine that use them. */
54 struct nonbonded_verlet_t;
55 struct bonded_threading_t;
57 class DispersionCorrection;
65 class DeviceStreamManager;
68 class StatePropagatorDataGpu;
70 class WholeMoleculeTransform;
73 /* macros for the cginfo data in forcerec
75 * Since the tpx format support max 256 energy groups, we do the same here.
76 * Note that we thus have bits 8-14 still unused.
78 * The maximum cg size in cginfo is 63
79 * because we only have space for 6 bits in cginfo,
80 * this cg size entry is actually only read with domain decomposition.
82 #define SET_CGINFO_GID(cgi, gid) (cgi) = (((cgi) & ~255) | (gid))
83 #define GET_CGINFO_GID(cgi) ((cgi)&255)
84 #define SET_CGINFO_FEP(cgi) (cgi) = ((cgi) | (1 << 15))
85 #define GET_CGINFO_FEP(cgi) ((cgi) & (1 << 15))
86 #define SET_CGINFO_EXCL_INTER(cgi) (cgi) = ((cgi) | (1 << 17))
87 #define GET_CGINFO_EXCL_INTER(cgi) ((cgi) & (1 << 17))
88 #define SET_CGINFO_CONSTR(cgi) (cgi) = ((cgi) | (1 << 20))
89 #define GET_CGINFO_CONSTR(cgi) ((cgi) & (1 << 20))
90 #define SET_CGINFO_SETTLE(cgi) (cgi) = ((cgi) | (1 << 21))
91 #define GET_CGINFO_SETTLE(cgi) ((cgi) & (1 << 21))
92 /* This bit is only used with bBondComm in the domain decomposition */
93 #define SET_CGINFO_BOND_INTER(cgi) (cgi) = ((cgi) | (1 << 22))
94 #define GET_CGINFO_BOND_INTER(cgi) ((cgi) & (1 << 22))
95 #define SET_CGINFO_HAS_VDW(cgi) (cgi) = ((cgi) | (1 << 23))
96 #define GET_CGINFO_HAS_VDW(cgi) ((cgi) & (1 << 23))
97 #define SET_CGINFO_HAS_Q(cgi) (cgi) = ((cgi) | (1 << 24))
98 #define GET_CGINFO_HAS_Q(cgi) ((cgi) & (1 << 24))
101 /* Value to be used in mdrun for an infinite cut-off.
102 * Since we need to compare with the cut-off squared,
103 * this value should be slighlty smaller than sqrt(GMX_FLOAT_MAX).
105 #define GMX_CUTOFF_INF 1E+18
107 /* enums for the neighborlist type */
122 std::vector<int> cginfo;
126 /* Forward declaration of type for managing Ewald tables */
127 struct gmx_ewald_tab_t;
129 struct ewald_corr_thread_t;
131 /*! \brief Helper force buffers for ForceOutputs
133 * This class stores intermediate force buffers that are used
134 * internally in the force calculation and which are reduced into
135 * the output force buffer passed to the force calculation.
137 class ForceHelperBuffers
140 /*! \brief Constructs helper buffers
142 * When the forces that will be accumulated with help of these buffers
143 * have direct virial contributions, set the parameter to true, so
144 * an extra force buffer is available for these forces to enable
145 * correct virial computation.
147 ForceHelperBuffers(bool haveDirectVirialContributions);
149 //! Returns whether we have a direct virial contribution force buffer
150 bool haveDirectVirialContributions() const { return haveDirectVirialContributions_; }
152 //! Returns the buffer for direct virial contributions
153 gmx::ArrayRef<gmx::RVec> forceBufferForDirectVirialContributions()
155 GMX_ASSERT(haveDirectVirialContributions_, "Buffer can only be requested when present");
156 return forceBufferForDirectVirialContributions_;
159 //! Returns the buffer for shift forces, size SHIFTS
160 gmx::ArrayRef<gmx::RVec> shiftForces() { return shiftForces_; }
162 //! Resizes the direct virial contribution buffer, when present
163 void resize(int numAtoms);
166 //! True when we have contributions that are directly added to the virial
167 bool haveDirectVirialContributions_ = false;
168 //! Force buffer for force computation with direct virial contributions
169 std::vector<gmx::RVec> forceBufferForDirectVirialContributions_;
170 //! Shift force array for computing the virial, size SHIFTS
171 std::vector<gmx::RVec> shiftForces_;
175 { // NOLINT (clang-analyzer-optin.performance.Padding)
176 // Declare an explicit constructor and destructor, so they can be
177 // implemented in a single source file, so that not every source
178 // file that includes this one needs to understand how to find the
179 // destructors of the objects pointed to by unique_ptr members.
183 struct interaction_const_t* ic = nullptr;
186 PbcType pbcType = PbcType::Xyz;
187 //! Tells whether atoms inside a molecule can be in different periodic images,
188 // i.e. whether we need to take into account PBC when computing distances inside molecules.
189 // This determines whether PBC must be considered for e.g. bonded interactions.
190 gmx_bool bMolPBC = FALSE;
192 rvec posres_com = { 0 };
193 rvec posres_comB = { 0 };
195 gmx_bool use_simd_kernels = FALSE;
197 /* Interaction for calculated in kernels. In many cases this is similar to
198 * the electrostatics settings in the inputrecord, but the difference is that
199 * these variables always specify the actual interaction in the kernel - if
200 * we are tabulating reaction-field the inputrec will say reaction-field, but
201 * the kernel interaction will say cubic-spline-table. To be safe we also
202 * have a kernel-specific setting for the modifiers - if the interaction is
203 * tabulated we already included the inputrec modification there, so the kernel
204 * modification setting will say 'none' in that case.
206 int nbkernel_elec_interaction = 0;
207 int nbkernel_vdw_interaction = 0;
208 int nbkernel_elec_modifier = 0;
209 int nbkernel_vdw_modifier = 0;
212 * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
216 /* Charge sum for topology A/B ([0]/[1]) for Ewald corrections */
217 double qsum[2] = { 0 };
218 double q2sum[2] = { 0 };
219 double c6sum[2] = { 0 };
221 /* Dispersion correction stuff */
222 std::unique_ptr<DispersionCorrection> dispersionCorrection;
228 gmx_bool bcoultab = FALSE;
229 gmx_bool bvdwtab = FALSE;
231 t_forcetable* pairsTable = nullptr; /* for 1-4 interactions, [pairs] and [pairs_nb] */
236 /* Information about atom properties for the molecule blocks in the system */
237 std::vector<cginfo_mb_t> cginfo_mb;
238 /* Information about atom properties for local and non-local atoms */
239 std::vector<int> cginfo;
241 rvec* shift_vec = nullptr;
243 std::unique_ptr<gmx::WholeMoleculeTransform> wholeMoleculeTransform;
245 /* The Nbnxm Verlet non-bonded machinery */
246 std::unique_ptr<nonbonded_verlet_t> nbv;
248 /* The wall tables (if used) */
250 t_forcetable*** wall_tab = nullptr;
252 /* The number of atoms participating in do_force_lowlevel */
253 int natoms_force = 0;
254 /* The number of atoms participating in force calculation and constraints */
255 int natoms_force_constr = 0;
257 /* Helper buffer for ForceOutputs */
258 std::unique_ptr<ForceHelperBuffers> forceHelperBuffers;
260 /* Data for PPPM/PME/Ewald */
261 struct gmx_pme_t* pmedata = nullptr;
262 int ljpme_combination_rule = 0;
264 /* PME/Ewald stuff */
265 struct gmx_ewald_tab_t* ewald_table = nullptr;
267 /* Non bonded Parameter lists */
268 int ntype = 0; /* Number of atom types */
269 gmx_bool bBHAM = FALSE;
270 std::vector<real> nbfp;
271 real* ljpme_c6grid = nullptr; /* C6-values used on grid in LJPME */
273 /* Energy group pair flags */
274 int* egp_flags = nullptr;
276 /* Shell molecular dynamics flexible constraints */
277 real fc_stepsize = 0;
279 /* If > 0 signals Test Particle Insertion,
280 * the value is the number of atoms of the molecule to insert
281 * Only the energy difference due to the addition of the last molecule
282 * should be calculated.
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 /* Data for special listed force calculations */
300 std::unique_ptr<t_fcdata> fcdata;
302 // The listed forces calculation data, 1 entry or multiple entries with multiple time stepping
303 std::vector<ListedForces> listedForces;
305 static constexpr size_t c_listedForcesFast = 0;
307 /* TODO: Replace the pointer by an object once we got rid of C */
308 gmx::GpuBonded* gpuBonded = nullptr;
310 /* Ewald correction thread local virial and energy data */
312 struct ewald_corr_thread_t* ewc_t = nullptr;
314 gmx::ForceProviders* forceProviders = nullptr;
316 // The stateGpu object is created in runner, forcerec just keeps the copy of the pointer.
317 // TODO: This is not supposed to be here. StatePropagatorDataGpu should be a part of
318 // general StatePropagatorData object that is passed around
319 gmx::StatePropagatorDataGpu* stateGpu = nullptr;
320 // TODO: Should not be here. This is here only to pass the pointer around.
321 gmx::DeviceStreamManager* deviceStreamManager = nullptr;
323 //! GPU device context
324 DeviceContext* deviceContext = nullptr;
326 /* For PME-PP GPU communication */
327 std::unique_ptr<gmx::PmePpCommGpu> pmePpCommGpu;
330 /* Important: Starting with Gromacs-4.6, the values of c6 and c12 in the nbfp array have
331 * been scaled by 6.0 or 12.0 to save flops in the kernels. We have corrected this everywhere
332 * in the code, but beware if you are using these macros externally.
334 #define C6(nbfp, ntp, ai, aj) (nbfp)[2 * ((ntp) * (ai) + (aj))]
335 #define C12(nbfp, ntp, ai, aj) (nbfp)[2 * ((ntp) * (ai) + (aj)) + 1]
336 #define BHAMC(nbfp, ntp, ai, aj) (nbfp)[3 * ((ntp) * (ai) + (aj))]
337 #define BHAMA(nbfp, ntp, ai, aj) (nbfp)[3 * ((ntp) * (ai) + (aj)) + 1]
338 #define BHAMB(nbfp, ntp, ai, aj) (nbfp)[3 * ((ntp) * (ai) + (aj)) + 2]