Split lines with many copyright years
[gromacs.git] / src / gromacs / mdtypes / interaction_const.h
blobbc09f47c92a6e129702e5af842f3197cdf4a3acf
1 /*
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, 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
39 #include <memory>
40 #include <vector>
42 #include "gromacs/mdtypes/md_enums.h"
43 #include "gromacs/utility/alignedallocator.h"
44 #include "gromacs/utility/basedefinitions.h"
45 #include "gromacs/utility/real.h"
47 /* Used with force switching or a constant potential shift:
48 * rsw = max(r - r_switch, 0)
49 * force/p = r^-(p+1) + c2*rsw^2 + c3*rsw^3
50 * potential = r^-p + c2/3*rsw^3 + c3/4*rsw^4 + cpot
51 * With a constant potential shift c2 and c3 are both 0.
53 struct shift_consts_t
55 real c2;
56 real c3;
57 real cpot;
60 /* Used with potential switching:
61 * rsw = max(r - r_switch, 0)
62 * sw = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
63 * dsw = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
64 * force = force*dsw - potential*sw
65 * potential *= sw
67 struct switch_consts_t
69 real c3;
70 real c4;
71 real c5;
74 /* Convenience type for vector with aligned memory */
75 template<typename T>
76 using AlignedVector = std::vector<T, gmx::AlignedAllocator<T>>;
78 /* Force/energy interpolation tables for Ewald long-range corrections
80 * Interpolation is linear for the force, quadratic for the potential.
82 struct EwaldCorrectionTables
84 // 1/table_spacing, units 1/nm
85 real scale = 0;
86 // Force table
87 AlignedVector<real> tableF;
88 // Energy table
89 AlignedVector<real> tableV;
90 // Coulomb force+energy table, size of array is tabq_size*4,
91 // entry quadruplets are: F[i], F[i+1]-F[i], V[i], 0,
92 // this is used with 4-wide SIMD for aligned loads
93 AlignedVector<real> tableFDV0;
96 /* The physical interaction parameters for non-bonded interaction calculations
98 * This struct contains copies of the physical interaction parameters
99 * from the user input as well as processed values that are need in
100 * non-bonded interaction kernels.
102 * The default constructor gives plain Coulomb and LJ interactions cut off
103 * a 1 nm without potential shifting and a Coulomb pre-factor of 1.
105 struct interaction_const_t
107 int cutoff_scheme = ecutsVERLET;
109 /* VdW */
110 int vdwtype = evdwCUT;
111 int vdw_modifier = eintmodNONE;
112 double reppow = 12;
113 real rvdw = 1;
114 real rvdw_switch = 0;
115 struct shift_consts_t dispersion_shift = { 0, 0, 0 };
116 struct shift_consts_t repulsion_shift = { 0, 0, 0 };
117 struct switch_consts_t vdw_switch = { 0, 0, 0 };
118 gmx_bool useBuckingham = false;
119 real buckinghamBMax = 0;
121 /* type of electrostatics */
122 int eeltype = eelCUT;
123 int coulomb_modifier = eintmodNONE;
125 /* Coulomb */
126 real rcoulomb = 1;
127 real rcoulomb_switch = 0;
129 /* PME/Ewald */
130 real ewaldcoeff_q = 0;
131 real ewaldcoeff_lj = 0;
132 int ljpme_comb_rule = eljpmeGEOM; /* LJ combination rule for the LJ PME mesh part */
133 real sh_ewald = 0; /* -sh_ewald is added to the direct space potential */
134 real sh_lj_ewald = 0; /* sh_lj_ewald is added to the correction potential */
136 /* Dielectric constant resp. multiplication factor for charges */
137 real epsilon_r = 1;
138 real epsfac = 1;
140 /* Constants for reaction-field or plain cut-off */
141 real epsilon_rf = 1;
142 real k_rf = 0;
143 real c_rf = 0;
145 // Coulomb Ewald correction table
146 std::unique_ptr<EwaldCorrectionTables> coulombEwaldTables;
147 /* Note that a Van der Waals Ewald correction table
148 * of type EwaldCorrectionTables can be added here if wanted.
152 #endif