6682ef7d5b3ea55747ccb6a60d1c3b29083dc5ba
[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, 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 "gromacs/utility/real.h"
39
40 #ifdef __cplusplus
41 extern "C" {
42 #endif
43
44 /* Used with force switching or a constant potential shift:
45  * rsw       = max(r - r_switch, 0)
46  * force/p   = r^-(p+1) + c2*rsw^2 + c3*rsw^3
47  * potential = r^-p + c2/3*rsw^3 + c3/4*rsw^4 + cpot
48  * With a constant potential shift c2 and c3 are both 0.
49  */
50 struct shift_consts_t
51 {
52     real c2;
53     real c3;
54     real cpot;
55 };
56
57 /* Used with potential switching:
58  * rsw        = max(r - r_switch, 0)
59  * sw         = 1 + c3*rsw^3 + c4*rsw^4 + c5*rsw^5
60  * dsw        = 3*c3*rsw^2 + 4*c4*rsw^3 + 5*c5*rsw^4
61  * force      = force*dsw - potential*sw
62  * potential *= sw
63  */
64 struct switch_consts_t
65 {
66     real c3;
67     real c4;
68     real c5;
69 };
70
71 struct interaction_const_t
72 {
73     int             cutoff_scheme;
74
75     /* VdW */
76     int                    vdwtype;
77     int                    vdw_modifier;
78     real                   rvdw;
79     real                   rvdw_switch;
80     struct shift_consts_t  dispersion_shift;
81     struct shift_consts_t  repulsion_shift;
82     struct switch_consts_t vdw_switch;
83     /* TODO: remove this variable, used for not modyfing the group kernels,
84      * it is equal to -dispersion_shift->cpot
85      */
86     real sh_invrc6;
87
88     /* type of electrostatics (defined in enums.h) */
89     int  eeltype;
90     int  coulomb_modifier;
91
92     /* Coulomb */
93     real rcoulomb;
94
95     /* Cut-off */
96     real rlist;
97
98     /* PME/Ewald */
99     real ewaldcoeff_q;
100     real ewaldcoeff_lj;
101     int  ljpme_comb_rule; /* LJ combination rule for the LJ PME mesh part */
102     real sh_ewald;        /* -sh_ewald is added to the direct space potential */
103     real sh_lj_ewald;     /* sh_lj_ewald is added to the correction potential */
104
105     /* Dielectric constant resp. multiplication factor for charges */
106     real epsilon_r;
107     real epsfac;
108
109     /* Constants for reaction-field or plain cut-off */
110     real epsilon_rf;
111     real k_rf;
112     real c_rf;
113
114     /* Force/energy interpolation tables, linear in force, quadratic in V */
115     real  tabq_scale;
116     int   tabq_size;
117     /* Coulomb force table, size of array is tabq_size (when used) */
118     real *tabq_coul_F;
119     /* Coulomb energy table, size of array is tabq_size (when used) */
120     real *tabq_coul_V;
121     /* Coulomb force+energy table, size of array is tabq_size*4,
122        entry quadruplets are: F[i], F[i+1]-F[i], V[i], 0,
123        this is used with single precision x86 SIMD for aligned loads */
124     real *tabq_coul_FDV0;
125
126     /* Vdw force table for LJ-PME, size of array is tabq_size (when used) */
127     real *tabq_vdw_F;
128     /* Vdw energy table for LJ-PME, size of array is tabq_size (when used) */
129     real *tabq_vdw_V;
130     /* Vdw force+energy table for LJ-PME, size of array is tabq_size*4, entry
131        quadruplets are: F[i], F[i+1]-F[i], V[i], 0, this is used with
132        single precision x86 SIMD for aligned loads */
133     real *tabq_vdw_FDV0;
134 };
135
136 #ifdef __cplusplus
137 }
138 #endif
139
140 #endif