2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2017 by the GROMACS development team.
5 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
36 #ifndef GMX_MDTYPES_INTERACTION_CONST_H
37 #define GMX_MDTYPES_INTERACTION_CONST_H
44 #include "gromacs/mdtypes/md_enums.h"
45 #include "gromacs/utility/alignedallocator.h"
46 #include "gromacs/utility/real.h"
52 /* Used with force switching or a constant potential shift:
53 * rsw = max(r - r_switch, 0)
54 * force/p = r^-(p+1) + c2*rsw^2 + c3*rsw^3
55 * potential = r^-p + c2/3*rsw^3 + c3/4*rsw^4 + cpot
56 * With a constant potential shift c2 and c3 are both 0.
65 /* Used with potential switching:
66 * rsw = max(r - r_switch, 0)
67 * sw = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
68 * dsw = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
69 * force = force*dsw - potential*sw
72 struct switch_consts_t
79 /* Convenience type for vector with aligned memory */
81 using AlignedVector = std::vector<T, gmx::AlignedAllocator<T>>;
83 /* Force/energy interpolation tables for Ewald long-range corrections
85 * Interpolation is linear for the force, quadratic for the potential.
87 struct EwaldCorrectionTables
89 // 1/table_spacing, units 1/nm
92 AlignedVector<real> tableF;
94 AlignedVector<real> tableV;
95 // Coulomb force+energy table, size of array is tabq_size*4,
96 // entry quadruplets are: F[i], F[i+1]-F[i], V[i], 0,
97 // this is used with 4-wide SIMD for aligned loads
98 AlignedVector<real> tableFDV0;
101 /* The physical interaction parameters for non-bonded interaction calculations
103 * This struct contains copies of the physical interaction parameters
104 * from the user input as well as processed values that are need in
105 * non-bonded interaction kernels.
107 * The default constructor gives plain Coulomb and LJ interactions cut off
108 * a 1 nm without potential shifting and a Coulomb pre-factor of 1.
110 struct interaction_const_t
112 /* This struct contains the soft-core parameters from t_lambda,
113 * but processed for direct use in the kernels.
115 struct SoftCoreParameters
118 SoftCoreParameters(const t_lambda& fepvals);
120 // Alpha parameter for Van der Waals interactions
122 // Alpha parameter for Coulomb interactions
124 // Exponent for the dependence of the soft-core on lambda
126 // Value for sigma^6 for LJ interaction with C6<=0 and/or C12<=0
127 real sigma6WithInvalidSigma;
128 // Minimum value for sigma^6, used when soft-core is applied to Coulomb interactions
130 // soft-core function
131 SoftcoreType softcoreType;
132 // (gapsys sc) linearization point scaling for vdW interactions
133 real gapsysScaleLinpointVdW;
134 // (gapsys sc) linearization point scaling for Coulomb interactions
135 real gapsysScaleLinpointCoul;
136 // (gapsys sc) lower bound/replacement for c12/c6 in vdw interactions
137 real gapsysSigma6VdW;
141 VanDerWaalsType vdwtype = VanDerWaalsType::Cut;
142 InteractionModifiers vdw_modifier = InteractionModifiers::None;
145 real rvdw_switch = 0;
146 struct shift_consts_t dispersion_shift = { 0, 0, 0 };
147 struct shift_consts_t repulsion_shift = { 0, 0, 0 };
148 struct switch_consts_t vdw_switch = { 0, 0, 0 };
149 bool useBuckingham = false;
150 real buckinghamBMax = 0;
152 /* type of electrostatics */
153 CoulombInteractionType eeltype = CoulombInteractionType::Cut;
154 InteractionModifiers coulomb_modifier = InteractionModifiers::None;
158 real rcoulomb_switch = 0;
161 real ewaldcoeff_q = 0;
162 real ewaldcoeff_lj = 0;
163 LongRangeVdW ljpme_comb_rule = LongRangeVdW::Geom; /* LJ combination rule for the LJ PME mesh part */
164 real sh_ewald = 0; /* -sh_ewald is added to the direct space potential */
165 real sh_lj_ewald = 0; /* sh_lj_ewald is added to the correction potential */
167 /* Dielectric constant resp. multiplication factor for charges */
171 /* Constants for reaction-field or plain cut-off */
172 //! Dielectric constant for reaction field beyond the cutoff distance
173 real reactionFieldPermitivity = 1;
174 //! Coefficient for reaction field; scales relation between epsilon_r and reactionFieldPermitivity
175 real reactionFieldCoefficient = 0;
176 //! Constant shift to reaction field Coulomb interaction to make potential an integral of force
177 real reactionFieldShift = 0;
179 // Coulomb Ewald correction table
180 std::unique_ptr<EwaldCorrectionTables> coulombEwaldTables;
181 // Van der Waals Ewald correction table
182 std::unique_ptr<EwaldCorrectionTables> vdwEwaldTables;
184 // Free-energy parameters, only present when free-energy calculations are requested
185 std::unique_ptr<SoftCoreParameters> softCoreParameters;
188 /*! \brief Construct interaction constants
190 * This data is used (particularly) by search and force code for
191 * short-range interactions. Many of these are constant for the whole
192 * simulation; some are constant only after PME tuning completes.
194 interaction_const_t init_interaction_const(FILE* fp,
195 const t_inputrec& ir,
196 const gmx_mtop_t& mtop,
197 bool systemHasNetCharge);