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