Move ifunc stuff from idef.h to ifunc.h
[alexxy/gromacs.git] / src / gromacs / topology / ifunc.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2012,2014,2015,2016,2018, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #ifndef GMX_TOPOLOGY_IFUNC_H
38 #define GMX_TOPOLOGY_IFUNC_H
39
40 #include "gromacs/math/vectypes.h"
41
42 struct t_fcdata;
43 struct t_graph;
44 union t_iparams;
45 struct t_mdatoms;
46 struct t_pbc;
47
48 /* TODO: Remove this typedef when t_ilist is removed */
49 typedef int t_iatom;
50
51 /* Real vector type with an additional, unused 4th element.
52  * This type is used to allow aligned 4-wide SIMD loads and stores.
53  */
54 typedef real rvec4[4];
55
56 typedef real t_ifunc (int nbonds, const t_iatom iatoms[],
57                       const t_iparams iparams[],
58                       const rvec x[], rvec4 f[], rvec fshift[],
59                       const struct t_pbc *pbc, const struct t_graph *g,
60                       real lambda, real *dvdlambda,
61                       const struct t_mdatoms *md, struct t_fcdata *fcd,
62                       int *ddgatindex);
63
64 /*
65  * The function type t_ifunc() calculates one interaction, using iatoms[]
66  * and iparams. Within the function the number of atoms to be used is
67  * known. Within the function only the atomid part of the iatoms[] array
68  * is supplied, not the type field (see also t_ilist). The function
69  * returns the potential energy. If pbc==NULL the coordinates in x are
70  * assumed to be such that no calculation of PBC is necessary,
71  * If pbc!=NULL a full PBC calculation is performed.
72  * If g!=NULL it is used for determining the shift forces.
73  * With domain decomposition ddgatindex can be used for getting global
74  * atom numbers for warnings and error messages.
75  * ddgatindex is NULL when domain decomposition is not used.
76  */
77
78 constexpr unsigned int IF_NULL       = 0;
79 constexpr unsigned int IF_BOND       = 1 << 0;
80 constexpr unsigned int IF_VSITE      = 1 << 1;
81 constexpr unsigned int IF_CONSTRAINT = 1 << 2;
82 constexpr unsigned int IF_CHEMBOND   = 1 << 3;
83 constexpr unsigned int IF_BTYPE      = 1 << 4;
84 constexpr unsigned int IF_ATYPE      = 1 << 5;
85 constexpr unsigned int IF_TABULATED  = 1 << 6;
86 constexpr unsigned int IF_LIMZERO    = 1 << 7;
87 /* These flags tell to some of the routines what can be done with this
88  * item in the list.
89  * With IF_BOND a bonded interaction will be calculated.
90  * With IF_BTYPE grompp can convert the bond to a Morse potential.
91  * With IF_BTYPE or IF_ATYPE the bond/angle can be converted to
92  * a constraint or used for vsite parameter determination by grompp.
93  * IF_LIMZERO indicates that for a bonded interaction the potential
94  * does goes to zero for large distances, thus if such an interaction
95  * it not assigned to any node by the domain decompostion, the simulation
96  * still continue, if mdrun has been told so.
97  */
98
99 struct t_interaction_function // NOLINT (clang-analyzer-optin.performance.Padding)
100 {
101     const char *name;         /* the name of this function                      */
102     const char *longname;     /* The name for printing etc.                   */
103     int         nratoms;      /* nr of atoms needed for this function           */
104     int         nrfpA, nrfpB; /* number of parameters for this function.      */
105                               /* this corresponds to the number of params in  */
106                               /* iparams struct! (see idef.h)                 */
107     /* A and B are for normal and free energy components respectively.    */
108     unsigned int    flags;    /* Flags (see above)                            */
109     int             nrnb_ind; /* index for nrnb (-1 if unknown)               */
110     t_ifunc        *ifunc;    /* the function it self                           */
111 };
112
113 #define NRFPA(ftype) (interaction_function[(ftype)].nrfpA)
114 #define NRFPB(ftype) (interaction_function[(ftype)].nrfpB)
115 #define NRFP(ftype)  (NRFPA(ftype)+NRFPB(ftype))
116 #define NRAL(ftype) (interaction_function[(ftype)].nratoms)
117
118 #define IS_CHEMBOND(ftype) (interaction_function[(ftype)].nratoms == 2 && (interaction_function[(ftype)].flags & IF_CHEMBOND))
119 /* IS_CHEMBOND tells if function type ftype represents a chemical bond */
120
121 /* IS_ANGLE tells if a function type ftype represents an angle
122  * Per Larsson, 2007-11-06
123  */
124 #define IS_ANGLE(ftype) (interaction_function[(ftype)].nratoms == 3 && (interaction_function[(ftype)].flags & IF_ATYPE))
125 #define IS_VSITE(ftype) (interaction_function[(ftype)].flags & IF_VSITE)
126
127 #define IS_TABULATED(ftype) (interaction_function[(ftype)].flags & IF_TABULATED)
128
129 /* this MUST correspond to the
130    t_interaction_function[F_NRE] in gmxlib/ifunc.c */
131 enum
132 {
133     F_BONDS,
134     F_G96BONDS,
135     F_MORSE,
136     F_CUBICBONDS,
137     F_CONNBONDS,
138     F_HARMONIC,
139     F_FENEBONDS,
140     F_TABBONDS,
141     F_TABBONDSNC,
142     F_RESTRBONDS,
143     F_ANGLES,
144     F_G96ANGLES,
145     F_RESTRANGLES,
146     F_LINEAR_ANGLES,
147     F_CROSS_BOND_BONDS,
148     F_CROSS_BOND_ANGLES,
149     F_UREY_BRADLEY,
150     F_QUARTIC_ANGLES,
151     F_TABANGLES,
152     F_PDIHS,
153     F_RBDIHS,
154     F_RESTRDIHS,
155     F_CBTDIHS,
156     F_FOURDIHS,
157     F_IDIHS,
158     F_PIDIHS,
159     F_TABDIHS,
160     F_CMAP,
161     F_GB12_NOLONGERUSED,
162     F_GB13_NOLONGERUSED,
163     F_GB14_NOLONGERUSED,
164     F_GBPOL_NOLONGERUSED,
165     F_NPSOLVATION_NOLONGERUSED,
166     F_LJ14,
167     F_COUL14,
168     F_LJC14_Q,
169     F_LJC_PAIRS_NB,
170     F_LJ,
171     F_BHAM,
172     F_LJ_LR_NOLONGERUSED,
173     F_BHAM_LR_NOLONGERUSED,
174     F_DISPCORR,
175     F_COUL_SR,
176     F_COUL_LR_NOLONGERUSED,
177     F_RF_EXCL,
178     F_COUL_RECIP,
179     F_LJ_RECIP,
180     F_DPD,
181     F_POLARIZATION,
182     F_WATER_POL,
183     F_THOLE_POL,
184     F_ANHARM_POL,
185     F_POSRES,
186     F_FBPOSRES,
187     F_DISRES,
188     F_DISRESVIOL,
189     F_ORIRES,
190     F_ORIRESDEV,
191     F_ANGRES,
192     F_ANGRESZ,
193     F_DIHRES,
194     F_DIHRESVIOL,
195     F_CONSTR,
196     F_CONSTRNC,
197     F_SETTLE,
198     F_VSITE2,
199     F_VSITE3,
200     F_VSITE3FD,
201     F_VSITE3FAD,
202     F_VSITE3OUT,
203     F_VSITE4FD,
204     F_VSITE4FDN,
205     F_VSITEN,
206     F_COM_PULL,
207     F_EQM,
208     F_EPOT,
209     F_EKIN,
210     F_ETOT,
211     F_ECONSERVED,
212     F_TEMP,
213     F_VTEMP_NOLONGERUSED,
214     F_PDISPCORR,
215     F_PRES,
216     F_DVDL_CONSTR,
217     F_DVDL,
218     F_DKDL,
219     F_DVDL_COUL,
220     F_DVDL_VDW,
221     F_DVDL_BONDED,
222     F_DVDL_RESTRAINT,
223     F_DVDL_TEMPERATURE, /* not calculated for now, but should just be the energy (NVT) or enthalpy (NPT), or 0 (NVE) */
224     F_NRE               /* This number is for the total number of energies      */
225 };
226
227 static inline bool IS_RESTRAINT_TYPE(int ifunc)
228 {
229     return
230         ifunc == F_POSRES || ifunc == F_FBPOSRES ||
231         ifunc == F_DISRES || ifunc == F_RESTRBONDS || ifunc == F_DISRESVIOL ||
232         ifunc == F_ORIRES || ifunc == F_ORIRESDEV ||
233         ifunc == F_ANGRES || ifunc == F_ANGRESZ || ifunc == F_DIHRES;
234 }
235
236 /* Maximum allowed number of atoms, parameters and terms in interaction_function.
237  * Check kernel/toppush.c when you change these numbers.
238  */
239 constexpr int MAXATOMLIST   = 6;
240 constexpr int MAXFORCEPARAM = 12;
241 constexpr int NR_RBDIHS     = 6;
242 constexpr int NR_CBTDIHS    = 6;
243 constexpr int NR_FOURDIHS   = 4;
244
245 extern const t_interaction_function interaction_function[F_NRE];
246 /* initialised interaction functions descriptor                         */
247
248 #endif