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,2021, 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/ishift.h"
48 #include "gromacs/pbcutil/pbc.h"
49 #include "gromacs/utility/arrayref.h"
50 #include "gromacs/utility/real.h"
54 /* Abstract type for PME that is defined only in the routine that use them. */
56 struct nonbonded_verlet_t;
57 struct bonded_threading_t;
59 class DispersionCorrection;
63 struct interaction_const_t;
67 class DeviceStreamManager;
69 class GpuForceReduction;
71 class StatePropagatorDataGpu;
73 class WholeMoleculeTransform;
76 /* macros for the cginfo data in forcerec
78 * Since the tpx format support max 256 energy groups, we do the same here.
79 * Note that we thus have bits 8-14 still unused.
81 * The maximum cg size in cginfo is 63
82 * because we only have space for 6 bits in cginfo,
83 * this cg size entry is actually only read with domain decomposition.
85 #define SET_CGINFO_GID(cgi, gid) (cgi) = (((cgi) & ~255) | (gid))
86 #define GET_CGINFO_GID(cgi) ((cgi)&255)
87 #define SET_CGINFO_FEP(cgi) (cgi) = ((cgi) | (1 << 15))
88 #define GET_CGINFO_FEP(cgi) ((cgi) & (1 << 15))
89 #define SET_CGINFO_EXCL_INTER(cgi) (cgi) = ((cgi) | (1 << 17))
90 #define GET_CGINFO_EXCL_INTER(cgi) ((cgi) & (1 << 17))
91 #define SET_CGINFO_CONSTR(cgi) (cgi) = ((cgi) | (1 << 20))
92 #define GET_CGINFO_CONSTR(cgi) ((cgi) & (1 << 20))
93 #define SET_CGINFO_SETTLE(cgi) (cgi) = ((cgi) | (1 << 21))
94 #define GET_CGINFO_SETTLE(cgi) ((cgi) & (1 << 21))
95 /* This bit is only used with bBondComm in the domain decomposition */
96 #define SET_CGINFO_BOND_INTER(cgi) (cgi) = ((cgi) | (1 << 22))
97 #define GET_CGINFO_BOND_INTER(cgi) ((cgi) & (1 << 22))
98 #define SET_CGINFO_HAS_VDW(cgi) (cgi) = ((cgi) | (1 << 23))
99 #define GET_CGINFO_HAS_VDW(cgi) ((cgi) & (1 << 23))
100 #define SET_CGINFO_HAS_Q(cgi) (cgi) = ((cgi) | (1 << 24))
101 #define GET_CGINFO_HAS_Q(cgi) ((cgi) & (1 << 24))
104 /* Value to be used in mdrun for an infinite cut-off.
105 * Since we need to compare with the cut-off squared,
106 * this value should be slighlty smaller than sqrt(GMX_FLOAT_MAX).
108 #define GMX_CUTOFF_INF 1E+18
115 std::vector<int> cginfo;
119 /* Forward declaration of type for managing Ewald tables */
120 struct gmx_ewald_tab_t;
122 struct ewald_corr_thread_t;
124 /*! \brief Helper force buffers for ForceOutputs
126 * This class stores intermediate force buffers that are used
127 * internally in the force calculation and which are reduced into
128 * the output force buffer passed to the force calculation.
130 class ForceHelperBuffers
133 /*! \brief Constructs helper buffers
135 * When the forces that will be accumulated with help of these buffers
136 * have direct virial contributions, set the parameter to true, so
137 * an extra force buffer is available for these forces to enable
138 * correct virial computation.
140 ForceHelperBuffers(bool haveDirectVirialContributions);
142 //! Returns whether we have a direct virial contribution force buffer
143 bool haveDirectVirialContributions() const { return haveDirectVirialContributions_; }
145 //! Returns the buffer for direct virial contributions
146 gmx::ArrayRef<gmx::RVec> forceBufferForDirectVirialContributions()
148 GMX_ASSERT(haveDirectVirialContributions_, "Buffer can only be requested when present");
149 return forceBufferForDirectVirialContributions_;
152 //! Returns the buffer for shift forces, size SHIFTS
153 gmx::ArrayRef<gmx::RVec> shiftForces() { return shiftForces_; }
155 //! Resizes the direct virial contribution buffer, when present
156 void resize(int numAtoms);
159 //! True when we have contributions that are directly added to the virial
160 bool haveDirectVirialContributions_ = false;
161 //! Force buffer for force computation with direct virial contributions
162 std::vector<gmx::RVec> forceBufferForDirectVirialContributions_;
163 //! Shift force array for computing the virial, size SHIFTS
164 std::vector<gmx::RVec> shiftForces_;
166 // NOLINTNEXTLINE (clang-analyzer-optin.performance.Padding)
169 // Declare an explicit constructor and destructor, so they can be
170 // implemented in a single source file, so that not every source
171 // file that includes this one needs to understand how to find the
172 // destructors of the objects pointed to by unique_ptr members.
176 std::unique_ptr<interaction_const_t> ic;
179 PbcType pbcType = PbcType::Xyz;
180 //! Tells whether atoms inside a molecule can be in different periodic images,
181 // i.e. whether we need to take into account PBC when computing distances inside molecules.
182 // This determines whether PBC must be considered for e.g. bonded interactions.
183 bool bMolPBC = false;
184 RefCoordScaling rc_scaling = RefCoordScaling::No;
185 gmx::RVec posres_com = { 0, 0, 0 };
186 gmx::RVec posres_comB = { 0, 0, 0 };
188 bool use_simd_kernels = false;
190 /* Interaction for calculated in kernels. In many cases this is similar to
191 * the electrostatics settings in the inputrecord, but the difference is that
192 * these variables always specify the actual interaction in the kernel - if
193 * we are tabulating reaction-field the inputrec will say reaction-field, but
194 * the kernel interaction will say cubic-spline-table. To be safe we also
195 * have a kernel-specific setting for the modifiers - if the interaction is
196 * tabulated we already included the inputrec modification there, so the kernel
197 * modification setting will say 'none' in that case.
199 NbkernelElecType nbkernel_elec_interaction = NbkernelElecType::None;
200 NbkernelVdwType nbkernel_vdw_interaction = NbkernelVdwType::None;
201 InteractionModifiers nbkernel_elec_modifier = InteractionModifiers::None;
202 InteractionModifiers nbkernel_vdw_modifier = InteractionModifiers::None;
205 * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
209 /* Charge sum for topology A/B ([0]/[1]) for Ewald corrections */
210 std::array<double, 2> qsum = { 0 };
211 std::array<double, 2> q2sum = { 0 };
212 std::array<double, 2> c6sum = { 0 };
214 /* Dispersion correction stuff */
215 std::unique_ptr<DispersionCorrection> dispersionCorrection;
220 std::unique_ptr<t_forcetable> pairsTable; /* for 1-4 interactions, [pairs] and [pairs_nb] */
223 FreeEnergyPerturbationType efep = FreeEnergyPerturbationType::No;
225 /* Information about atom properties for the molecule blocks in the system */
226 std::vector<cginfo_mb_t> cginfo_mb;
227 /* Information about atom properties for local and non-local atoms */
228 std::vector<int> cginfo;
230 std::vector<gmx::RVec> shift_vec;
232 std::unique_ptr<gmx::WholeMoleculeTransform> wholeMoleculeTransform;
234 /* The Nbnxm Verlet non-bonded machinery */
235 std::unique_ptr<nonbonded_verlet_t> nbv;
237 /* The wall tables (if used) */
239 std::vector<std::vector<std::unique_ptr<t_forcetable>>> wall_tab;
241 /* The number of atoms participating in do_force_lowlevel */
242 int natoms_force = 0;
243 /* The number of atoms participating in force calculation and constraints */
244 int natoms_force_constr = 0;
246 /* List of helper buffers for ForceOutputs, one for each time step with MTS */
247 std::vector<ForceHelperBuffers> forceHelperBuffers;
249 /* Data for PPPM/PME/Ewald */
250 gmx_pme_t* pmedata = nullptr;
251 LongRangeVdW ljpme_combination_rule = LongRangeVdW::Geom;
253 /* PME/Ewald stuff */
254 std::unique_ptr<gmx_ewald_tab_t> ewald_table;
256 /* Non bonded Parameter lists */
257 int ntype = 0; /* Number of atom types */
258 bool haveBuckingham = false;
259 std::vector<real> nbfp;
260 std::vector<real> ljpme_c6grid; /* C6-values used on grid in LJPME */
262 /* Energy group pair flags */
263 int* egp_flags = nullptr;
265 /* Shell molecular dynamics flexible constraints */
266 real fc_stepsize = 0;
268 /* If > 0 signals Test Particle Insertion,
269 * the value is the number of atoms of the molecule to insert
270 * Only the energy difference due to the addition of the last molecule
271 * should be calculated.
275 /* Limit for printing large forces, negative is don't print */
276 real print_force = 0;
278 /* User determined parameters, copied from the inputrec */
288 /* Tells whether we use multiple time stepping, computing some forces less frequently */
291 /* Data for special listed force calculations */
292 std::unique_ptr<t_fcdata> fcdata;
294 // The listed forces calculation data, 1 entry or multiple entries with multiple time stepping
295 std::vector<ListedForces> listedForces;
297 /* TODO: Replace the pointer by an object once we got rid of C */
298 gmx::GpuBonded* gpuBonded = nullptr;
300 /* Ewald correction thread local virial and energy data */
302 std::vector<ewald_corr_thread_t> ewc_t;
304 gmx::ForceProviders* forceProviders = nullptr;
306 // The stateGpu object is created in runner, forcerec just keeps the copy of the pointer.
307 // TODO: This is not supposed to be here. StatePropagatorDataGpu should be a part of
308 // general StatePropagatorData object that is passed around
309 gmx::StatePropagatorDataGpu* stateGpu = nullptr;
310 // TODO: Should not be here. This is here only to pass the pointer around.
311 gmx::DeviceStreamManager* deviceStreamManager = nullptr;
313 //! GPU device context
314 DeviceContext* deviceContext = nullptr;
316 /* For PME-PP GPU communication */
317 std::unique_ptr<gmx::PmePpCommGpu> pmePpCommGpu;
319 /* For GPU force reduction (on both local and non-local atoms) */
320 gmx::EnumerationArray<gmx::AtomLocality, std::unique_ptr<gmx::GpuForceReduction>> gpuForceReduction;
323 /* Important: Starting with Gromacs-4.6, the values of c6 and c12 in the nbfp array have
324 * been scaled by 6.0 or 12.0 to save flops in the kernels. We have corrected this everywhere
325 * in the code, but beware if you are using these macros externally.
327 #define C6(nbfp, ntp, ai, aj) (nbfp)[2 * ((ntp) * (ai) + (aj))]
328 #define C12(nbfp, ntp, ai, aj) (nbfp)[2 * ((ntp) * (ai) + (aj)) + 1]
329 #define BHAMC(nbfp, ntp, ai, aj) (nbfp)[3 * ((ntp) * (ai) + (aj))]
330 #define BHAMA(nbfp, ntp, ai, aj) (nbfp)[3 * ((ntp) * (ai) + (aj)) + 1]
331 #define BHAMB(nbfp, ntp, ai, aj) (nbfp)[3 * ((ntp) * (ai) + (aj)) + 2]