Remove obsolete mdp option ns-type
[gromacs.git] / src / gromacs / mdtypes / forcerec.h
blob78e826a4901dc47e393ccdd5f3fb890c8ee59b50
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,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
40 #include <array>
41 #include <memory>
42 #include <vector>
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. */
53 struct gmx_ns_t;
54 struct gmx_pme_t;
55 struct nonbonded_verlet_t;
56 struct bonded_threading_t;
57 class DispersionCorrection;
58 struct t_forcetable;
59 struct t_QMMMrec;
61 namespace gmx
63 class GpuBonded;
66 /* macros for the cginfo data in forcerec
68 * Since the tpx format support max 256 energy groups, we do the same here.
69 * Note that we thus have bits 8-14 still unused.
71 * The maximum cg size in cginfo is 63
72 * because we only have space for 6 bits in cginfo,
73 * this cg size entry is actually only read with domain decomposition.
75 #define SET_CGINFO_GID(cgi, gid) (cgi) = (((cgi) & ~255) | (gid))
76 #define GET_CGINFO_GID(cgi) ( (cgi) & 255)
77 #define SET_CGINFO_FEP(cgi) (cgi) = ((cgi) | (1<<15))
78 #define GET_CGINFO_FEP(cgi) ( (cgi) & (1<<15))
79 #define SET_CGINFO_EXCL_INTRA(cgi) (cgi) = ((cgi) | (1<<16))
80 #define GET_CGINFO_EXCL_INTRA(cgi) ( (cgi) & (1<<16))
81 #define SET_CGINFO_EXCL_INTER(cgi) (cgi) = ((cgi) | (1<<17))
82 #define GET_CGINFO_EXCL_INTER(cgi) ( (cgi) & (1<<17))
83 #define SET_CGINFO_SOLOPT(cgi, opt) (cgi) = (((cgi) & ~(3<<18)) | ((opt)<<18))
84 #define GET_CGINFO_SOLOPT(cgi) (((cgi)>>18) & 3)
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))
96 #define SET_CGINFO_NATOMS(cgi, opt) (cgi) = (((cgi) & ~(63<<25)) | ((opt)<<25))
97 #define GET_CGINFO_NATOMS(cgi) (((cgi)>>25) & 63)
100 /* Value to be used in mdrun for an infinite cut-off.
101 * Since we need to compare with the cut-off squared,
102 * this value should be slighlty smaller than sqrt(GMX_FLOAT_MAX).
104 #define GMX_CUTOFF_INF 1E+18
106 /* enums for the neighborlist type */
107 enum {
108 enbvdwNONE, enbvdwLJ, enbvdwBHAM, enbvdwTAB, enbvdwNR
111 struct cginfo_mb_t
113 int cg_start;
114 int cg_end;
115 int cg_mod;
116 int *cginfo;
120 /* Forward declaration of type for managing Ewald tables */
121 struct gmx_ewald_tab_t;
123 struct ewald_corr_thread_t;
125 struct t_forcerec { // NOLINT (clang-analyzer-optin.performance.Padding)
126 struct interaction_const_t *ic = nullptr;
128 /* PBC stuff */
129 int ePBC = 0;
130 //! Tells whether atoms inside a molecule can be in different periodic images,
131 // i.e. whether we need to take into account PBC when computing distances inside molecules.
132 // This determines whether PBC must be considered for e.g. bonded interactions.
133 gmx_bool bMolPBC = FALSE;
134 int rc_scaling = 0;
135 rvec posres_com = { 0 };
136 rvec posres_comB = { 0 };
138 gmx_bool use_simd_kernels = FALSE;
140 /* Interaction for calculated in kernels. In many cases this is similar to
141 * the electrostatics settings in the inputrecord, but the difference is that
142 * these variables always specify the actual interaction in the kernel - if
143 * we are tabulating reaction-field the inputrec will say reaction-field, but
144 * the kernel interaction will say cubic-spline-table. To be safe we also
145 * have a kernel-specific setting for the modifiers - if the interaction is
146 * tabulated we already included the inputrec modification there, so the kernel
147 * modification setting will say 'none' in that case.
149 int nbkernel_elec_interaction = 0;
150 int nbkernel_vdw_interaction = 0;
151 int nbkernel_elec_modifier = 0;
152 int nbkernel_vdw_modifier = 0;
154 /* Cut-Off stuff.
155 * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
157 real rlist = 0;
159 /* Charge sum and dipole for topology A/B ([0]/[1]) for Ewald corrections */
160 double qsum[2] = { 0 };
161 double q2sum[2] = { 0 };
162 double c6sum[2] = { 0 };
163 rvec mu_tot[2] = { { 0 } };
165 /* Dispersion correction stuff */
166 std::unique_ptr<DispersionCorrection> dispersionCorrection;
168 /* Fudge factors */
169 real fudgeQQ = 0;
171 /* Table stuff */
172 gmx_bool bcoultab = FALSE;
173 gmx_bool bvdwtab = FALSE;
175 t_forcetable *pairsTable = nullptr; /* for 1-4 interactions, [pairs] and [pairs_nb] */
177 /* Free energy */
178 int efep = 0;
179 real sc_alphavdw = 0;
180 real sc_alphacoul = 0;
181 int sc_power = 0;
182 real sc_r_power = 0;
183 real sc_sigma6_def = 0;
184 real sc_sigma6_min = 0;
186 /* NS Stuff */
187 int cg0 = 0;
188 int hcg = 0;
189 /* solvent_opt contains the enum for the most common solvent
190 * in the system, which will be optimized.
191 * It can be set to esolNO to disable all water optimization */
192 int solvent_opt = 0;
193 int nWatMol = 0;
194 gmx_bool bExcl_IntraCGAll_InterCGNone = FALSE;
195 struct cginfo_mb_t *cginfo_mb = nullptr;
196 std::vector<int> cginfo;
197 rvec *shift_vec = nullptr;
199 int cutoff_scheme = 0; /* group- or Verlet-style cutoff */
200 gmx_bool bNonbonded = FALSE; /* true if nonbonded calculations are *not* turned off */
202 /* The Nbnxm Verlet non-bonded machinery */
203 std::unique_ptr<nonbonded_verlet_t> nbv;
205 /* The wall tables (if used) */
206 int nwall = 0;
207 t_forcetable ***wall_tab = nullptr;
209 /* The number of charge groups participating in do_force_lowlevel */
210 int ncg_force = 0;
211 /* The number of atoms participating in do_force_lowlevel */
212 int natoms_force = 0;
213 /* The number of atoms participating in force and constraints */
214 int natoms_force_constr = 0;
215 /* The allocation size of vectors of size natoms_force */
216 int nalloc_force = 0;
218 /* Forces that should not enter into the coord x force virial summation:
219 * PPPM/PME/Ewald/posres/ForceProviders
221 /* True when we have contributions that are directly added to the virial */
222 bool haveDirectVirialContributions = false;
223 /* Force buffer for force computation with direct virial contributions */
224 std::vector<gmx::RVec> forceBufferForDirectVirialContributions;
226 /* Data for PPPM/PME/Ewald */
227 struct gmx_pme_t *pmedata = nullptr;
228 int ljpme_combination_rule = 0;
230 /* PME/Ewald stuff */
231 struct gmx_ewald_tab_t *ewald_table = nullptr;
233 /* Shift force array for computing the virial, size SHIFTS */
234 std::vector<gmx::RVec> shiftForces;
236 /* Non bonded Parameter lists */
237 int ntype = 0; /* Number of atom types */
238 gmx_bool bBHAM = FALSE;
239 real *nbfp = nullptr;
240 real *ljpme_c6grid = nullptr; /* C6-values used on grid in LJPME */
242 /* Energy group pair flags */
243 int *egp_flags = nullptr;
245 /* Shell molecular dynamics flexible constraints */
246 real fc_stepsize = 0;
248 /* If > 0 signals Test Particle Insertion,
249 * the value is the number of atoms of the molecule to insert
250 * Only the energy difference due to the addition of the last molecule
251 * should be calculated.
253 int n_tpi = 0;
255 /* QMMM stuff */
256 gmx_bool bQMMM = FALSE;
257 struct t_QMMMrec *qr = nullptr;
259 /* QM-MM neighborlists */
260 struct t_nblist *QMMMlist = nullptr;
262 /* Limit for printing large forces, negative is don't print */
263 real print_force = 0;
265 /* User determined parameters, copied from the inputrec */
266 int userint1 = 0;
267 int userint2 = 0;
268 int userint3 = 0;
269 int userint4 = 0;
270 real userreal1 = 0;
271 real userreal2 = 0;
272 real userreal3 = 0;
273 real userreal4 = 0;
275 /* Pointer to struct for managing threading of bonded force calculation */
276 struct bonded_threading_t *bondedThreading = nullptr;
278 /* TODO: Replace the pointer by an object once we got rid of C */
279 gmx::GpuBonded *gpuBonded = nullptr;
281 /* Ewald correction thread local virial and energy data */
282 int nthread_ewc = 0;
283 struct ewald_corr_thread_t *ewc_t = nullptr;
285 struct ForceProviders *forceProviders = nullptr;
288 /* Important: Starting with Gromacs-4.6, the values of c6 and c12 in the nbfp array have
289 * been scaled by 6.0 or 12.0 to save flops in the kernels. We have corrected this everywhere
290 * in the code, but beware if you are using these macros externally.
292 #define C6(nbfp, ntp, ai, aj) (nbfp)[2*((ntp)*(ai)+(aj))]
293 #define C12(nbfp, ntp, ai, aj) (nbfp)[2*((ntp)*(ai)+(aj))+1]
294 #define BHAMC(nbfp, ntp, ai, aj) (nbfp)[3*((ntp)*(ai)+(aj))]
295 #define BHAMA(nbfp, ntp, ai, aj) (nbfp)[3*((ntp)*(ai)+(aj))+1]
296 #define BHAMB(nbfp, ntp, ai, aj) (nbfp)[3*((ntp)*(ai)+(aj))+2]
298 #endif