Use proper doxygen tags in modular simulator
[gromacs.git] / src / gromacs / topology / ifunc.h
blob584c57f0b76fb6c34276478de5798b2b9363a74e
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
7 * Copyright (c) 2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 #ifndef GMX_TOPOLOGY_IFUNC_H
39 #define GMX_TOPOLOGY_IFUNC_H
41 #include "gromacs/math/vectypes.h"
43 struct t_fcdata;
44 struct t_graph;
45 union t_iparams;
46 struct t_mdatoms;
47 struct t_pbc;
49 /* TODO: Remove this typedef when t_ilist is removed */
50 typedef int t_iatom;
52 /* Real vector type with an additional, unused 4th element.
53 * This type is used to allow aligned 4-wide SIMD loads and stores.
55 typedef real rvec4[4];
58 * The function type t_ifunc() calculates one interaction, using iatoms[]
59 * and iparams. Within the function the number of atoms to be used is
60 * known. Within the function only the atomid part of the iatoms[] array
61 * is supplied, not the type field (see also t_ilist). The function
62 * returns the potential energy. If pbc==NULL the coordinates in x are
63 * assumed to be such that no calculation of PBC is necessary,
64 * If pbc!=NULL a full PBC calculation is performed.
65 * If g!=NULL it is used for determining the shift forces.
66 * With domain decomposition ddgatindex can be used for getting global
67 * atom numbers for warnings and error messages.
68 * ddgatindex is NULL when domain decomposition is not used.
71 constexpr unsigned int IF_NULL = 0;
72 constexpr unsigned int IF_BOND = 1 << 0;
73 constexpr unsigned int IF_VSITE = 1 << 1;
74 constexpr unsigned int IF_CONSTRAINT = 1 << 2;
75 constexpr unsigned int IF_CHEMBOND = 1 << 3;
76 constexpr unsigned int IF_BTYPE = 1 << 4;
77 constexpr unsigned int IF_ATYPE = 1 << 5;
78 constexpr unsigned int IF_TABULATED = 1 << 6;
79 constexpr unsigned int IF_LIMZERO = 1 << 7;
80 /* These flags tell to some of the routines what can be done with this
81 * item in the list.
82 * With IF_BOND a bonded interaction will be calculated.
83 * With IF_BTYPE grompp can convert the bond to a Morse potential.
84 * With IF_BTYPE or IF_ATYPE the bond/angle can be converted to
85 * a constraint or used for vsite parameter determination by grompp.
86 * IF_LIMZERO indicates that for a bonded interaction the potential
87 * does goes to zero for large distances, thus if such an interaction
88 * it not assigned to any node by the domain decompostion, the simulation
89 * still continue, if mdrun has been told so.
92 struct t_interaction_function // NOLINT (clang-analyzer-optin.performance.Padding)
94 const char* name; /* the name of this function */
95 const char* longname; /* The name for printing etc. */
96 int nratoms; /* nr of atoms needed for this function */
97 int nrfpA, nrfpB; /* number of parameters for this function. */
98 /* this corresponds to the number of params in */
99 /* iparams struct! (see idef.h) */
100 /* A and B are for normal and free energy components respectively. */
101 unsigned int flags; /* Flags (see above) */
104 #define NRFPA(ftype) (interaction_function[(ftype)].nrfpA)
105 #define NRFPB(ftype) (interaction_function[(ftype)].nrfpB)
106 #define NRFP(ftype) (NRFPA(ftype) + NRFPB(ftype))
107 #define NRAL(ftype) (interaction_function[(ftype)].nratoms)
109 #define IS_CHEMBOND(ftype) \
110 (interaction_function[(ftype)].nratoms == 2 && (interaction_function[(ftype)].flags & IF_CHEMBOND))
111 /* IS_CHEMBOND tells if function type ftype represents a chemical bond */
113 /* IS_ANGLE tells if a function type ftype represents an angle
114 * Per Larsson, 2007-11-06
116 #define IS_ANGLE(ftype) \
117 (interaction_function[(ftype)].nratoms == 3 && (interaction_function[(ftype)].flags & IF_ATYPE))
118 #define IS_VSITE(ftype) (interaction_function[(ftype)].flags & IF_VSITE)
120 #define IS_TABULATED(ftype) (interaction_function[(ftype)].flags & IF_TABULATED)
122 /* this MUST correspond to the
123 t_interaction_function[F_NRE] in gmxlib/ifunc.cpp */
124 enum
126 F_BONDS,
127 F_G96BONDS,
128 F_MORSE,
129 F_CUBICBONDS,
130 F_CONNBONDS,
131 F_HARMONIC,
132 F_FENEBONDS,
133 F_TABBONDS,
134 F_TABBONDSNC,
135 F_RESTRBONDS,
136 F_ANGLES,
137 F_G96ANGLES,
138 F_RESTRANGLES,
139 F_LINEAR_ANGLES,
140 F_CROSS_BOND_BONDS,
141 F_CROSS_BOND_ANGLES,
142 F_UREY_BRADLEY,
143 F_QUARTIC_ANGLES,
144 F_TABANGLES,
145 F_PDIHS,
146 F_RBDIHS,
147 F_RESTRDIHS,
148 F_CBTDIHS,
149 F_FOURDIHS,
150 F_IDIHS,
151 F_PIDIHS,
152 F_TABDIHS,
153 F_CMAP,
154 F_GB12_NOLONGERUSED,
155 F_GB13_NOLONGERUSED,
156 F_GB14_NOLONGERUSED,
157 F_GBPOL_NOLONGERUSED,
158 F_NPSOLVATION_NOLONGERUSED,
159 F_LJ14,
160 F_COUL14,
161 F_LJC14_Q,
162 F_LJC_PAIRS_NB,
163 F_LJ,
164 F_BHAM,
165 F_LJ_LR_NOLONGERUSED,
166 F_BHAM_LR_NOLONGERUSED,
167 F_DISPCORR,
168 F_COUL_SR,
169 F_COUL_LR_NOLONGERUSED,
170 F_RF_EXCL,
171 F_COUL_RECIP,
172 F_LJ_RECIP,
173 F_DPD,
174 F_POLARIZATION,
175 F_WATER_POL,
176 F_THOLE_POL,
177 F_ANHARM_POL,
178 F_POSRES,
179 F_FBPOSRES,
180 F_DISRES,
181 F_DISRESVIOL,
182 F_ORIRES,
183 F_ORIRESDEV,
184 F_ANGRES,
185 F_ANGRESZ,
186 F_DIHRES,
187 F_DIHRESVIOL,
188 F_CONSTR,
189 F_CONSTRNC,
190 F_SETTLE,
191 F_VSITE1,
192 F_VSITE2,
193 F_VSITE2FD,
194 F_VSITE3,
195 F_VSITE3FD,
196 F_VSITE3FAD,
197 F_VSITE3OUT,
198 F_VSITE4FD,
199 F_VSITE4FDN,
200 F_VSITEN,
201 F_COM_PULL,
202 F_DENSITYFITTING,
203 F_EQM,
204 F_EPOT,
205 F_EKIN,
206 F_ETOT,
207 F_ECONSERVED,
208 F_TEMP,
209 F_VTEMP_NOLONGERUSED,
210 F_PDISPCORR,
211 F_PRES,
212 F_DVDL_CONSTR,
213 F_DVDL,
214 F_DKDL,
215 F_DVDL_COUL,
216 F_DVDL_VDW,
217 F_DVDL_BONDED,
218 F_DVDL_RESTRAINT,
219 F_DVDL_TEMPERATURE, /* not calculated for now, but should just be the energy (NVT) or enthalpy (NPT), or 0 (NVE) */
220 F_NRE /* This number is for the total number of energies */
223 static inline bool IS_RESTRAINT_TYPE(int ifunc)
225 return ifunc == F_POSRES || ifunc == F_FBPOSRES || ifunc == F_DISRES || ifunc == F_RESTRBONDS
226 || ifunc == F_DISRESVIOL || ifunc == F_ORIRES || ifunc == F_ORIRESDEV
227 || ifunc == F_ANGRES || ifunc == F_ANGRESZ || ifunc == F_DIHRES;
230 /* Maximum allowed number of atoms, parameters and terms in interaction_function.
231 * Check kernel/toppush.c when you change these numbers.
233 constexpr int MAXATOMLIST = 6;
234 constexpr int MAXFORCEPARAM = 12;
235 constexpr int NR_RBDIHS = 6;
236 constexpr int NR_CBTDIHS = 6;
237 constexpr int NR_FOURDIHS = 4;
239 extern const t_interaction_function interaction_function[F_NRE];
240 /* initialised interaction functions descriptor */
242 #endif