Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / mdtypes / interaction_const.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2017,2018,2019, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 #ifndef GMX_MDTYPES_INTERACTION_CONST_H
36 #define GMX_MDTYPES_INTERACTION_CONST_H
37
38 #include <memory>
39 #include <vector>
40
41 #include "gromacs/mdtypes/md_enums.h"
42 #include "gromacs/utility/alignedallocator.h"
43 #include "gromacs/utility/basedefinitions.h"
44 #include "gromacs/utility/real.h"
45
46 /* Used with force switching or a constant potential shift:
47  * rsw       = max(r - r_switch, 0)
48  * force/p   = r^-(p+1) + c2*rsw^2 + c3*rsw^3
49  * potential = r^-p + c2/3*rsw^3 + c3/4*rsw^4 + cpot
50  * With a constant potential shift c2 and c3 are both 0.
51  */
52 struct shift_consts_t
53 {
54     real c2;
55     real c3;
56     real cpot;
57 };
58
59 /* Used with potential switching:
60  * rsw        = max(r - r_switch, 0)
61  * sw         = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
62  * dsw        = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
63  * force      = force*dsw - potential*sw
64  * potential *= sw
65  */
66 struct switch_consts_t
67 {
68     real c3;
69     real c4;
70     real c5;
71 };
72
73 /* Convenience type for vector with aligned memory */
74 template<typename T>
75 using AlignedVector = std::vector<T, gmx::AlignedAllocator<T>>;
76
77 /* Force/energy interpolation tables for Ewald long-range corrections
78  *
79  * Interpolation is linear for the force, quadratic for the potential.
80  */
81 struct EwaldCorrectionTables
82 {
83     // 1/table_spacing, units 1/nm
84     real scale = 0;
85     // Force table
86     AlignedVector<real> tableF;
87     // Energy table
88     AlignedVector<real> tableV;
89     // Coulomb force+energy table, size of array is tabq_size*4,
90     // entry quadruplets are: F[i], F[i+1]-F[i], V[i], 0,
91     // this is used with 4-wide SIMD for aligned loads
92     AlignedVector<real> tableFDV0;
93 };
94
95 /* The physical interaction parameters for non-bonded interaction calculations
96  *
97  * This struct contains copies of the physical interaction parameters
98  * from the user input as well as processed values that are need in
99  * non-bonded interaction kernels.
100  *
101  * The default constructor gives plain Coulomb and LJ interactions cut off
102  * a 1 nm without potential shifting and a Coulomb pre-factor of 1.
103  */
104 struct interaction_const_t
105 {
106     int cutoff_scheme = ecutsVERLET;
107
108     /* VdW */
109     int                    vdwtype          = evdwCUT;
110     int                    vdw_modifier     = eintmodNONE;
111     double                 reppow           = 12;
112     real                   rvdw             = 1;
113     real                   rvdw_switch      = 0;
114     struct shift_consts_t  dispersion_shift = { 0, 0, 0 };
115     struct shift_consts_t  repulsion_shift  = { 0, 0, 0 };
116     struct switch_consts_t vdw_switch       = { 0, 0, 0 };
117     gmx_bool               useBuckingham    = false;
118     real                   buckinghamBMax   = 0;
119
120     /* type of electrostatics */
121     int eeltype          = eelCUT;
122     int coulomb_modifier = eintmodNONE;
123
124     /* Coulomb */
125     real rcoulomb        = 1;
126     real rcoulomb_switch = 0;
127
128     /* PME/Ewald */
129     real ewaldcoeff_q    = 0;
130     real ewaldcoeff_lj   = 0;
131     int  ljpme_comb_rule = eljpmeGEOM; /* LJ combination rule for the LJ PME mesh part */
132     real sh_ewald        = 0;          /* -sh_ewald is added to the direct space potential */
133     real sh_lj_ewald     = 0;          /* sh_lj_ewald is added to the correction potential */
134
135     /* Dielectric constant resp. multiplication factor for charges */
136     real epsilon_r = 1;
137     real epsfac    = 1;
138
139     /* Constants for reaction-field or plain cut-off */
140     real epsilon_rf = 1;
141     real k_rf       = 0;
142     real c_rf       = 0;
143
144     // Coulomb Ewald correction table
145     std::unique_ptr<EwaldCorrectionTables> coulombEwaldTables;
146     /* Note that a Van der Waals Ewald correction table
147      * of type EwaldCorrectionTables can be added here if wanted.
148      */
149 };
150
151 #endif