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) 2013,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_IDEF_H
38 #define GMX_TOPOLOGY_IDEF_H
45 #include "gromacs/math/vectypes.h"
46 #include "gromacs/topology/ifunc.h"
47 #include "gromacs/utility/basedefinitions.h"
48 #include "gromacs/utility/gmxassert.h"
49 #include "gromacs/utility/real.h"
51 typedef union t_iparams
53 /* Some parameters have A and B values for free energy calculations.
54 * The B values are not used for regular simulations of course.
55 * Free Energy for nonbondeds can be computed by changing the atom type.
56 * The harmonic type is used for all harmonic potentials:
57 * bonds, angles and improper dihedrals
63 real rA, krA, rB, krB;
66 real klinA, aA, klinB, aB;
69 real lowA, up1A, up2A, kA, lowB, up1B, up2B, kB;
71 /* No free energy supported for cubic bonds, FENE, WPOL or cross terms */
82 real r1e, r2e, r3e, krt;
85 real thetaA, kthetaA, r13A, kUBA, thetaB, kthetaB, r13B, kUBB;
94 real alpha, drcut, khyp;
97 real al_x, al_y, al_z, rOH, rHH, rOD;
100 real a, alpha1, alpha2, rfac;
106 real c6A, c12A, c6B, c12B;
109 real fqq, qi, qj, c6, c12;
112 real qi, qj, c6, c12;
114 /* Proper dihedrals can not have different multiplicity when
115 * doing free energy calculations, because the potential would not
116 * be periodic anymore.
119 real phiA, cpA; int mult; real phiB, cpB;
124 /* Settle can not be used for Free energy calculations of water bond geometry.
125 * Use shake (or lincs) instead if you have to change the water bonds.
131 real b0A, cbA, betaA, b0B, cbB, betaB;
134 real pos0A[DIM], fcA[DIM], pos0B[DIM], fcB[DIM];
137 real pos0[DIM], r, k; int geom;
140 real rbcA[NR_RBDIHS], rbcB[NR_RBDIHS];
143 real cbtcA[NR_CBTDIHS], cbtcB[NR_CBTDIHS];
146 real a, b, c, d, e, f;
151 /* NOTE: npair is only set after reading the tpx file */
153 real low, up1, up2, kfac; int type, label, npair;
156 real phiA, dphiA, kfacA, phiB, dphiB, kfacB;
159 int ex, power, label; real c, obs, kfac;
162 int table; real kA; real kB;
168 real buf[MAXFORCEPARAM];
169 } generic; /* Conversion */
172 typedef int t_functype;
174 /* List of listed interactions, see description further down.
176 * TODO: Consider storing the function type as well.
177 * TODO: Consider providing per interaction access.
179 struct InteractionList
181 /* Returns the total number of elements in iatoms */
184 return iatoms.size();
187 /* List of interactions, see explanation further down */
188 std::vector<int> iatoms;
191 /* List of interaction lists, one list for each interaction type
193 * TODO: Consider only including entries in use instead of all F_NRE
195 typedef std::array<InteractionList, F_NRE> InteractionLists;
197 /* Deprecated list of listed interactions.
199 * The nonperturbed/perturbed interactions are now separated (sorted) in the
200 * ilist, such that the first 0..(nr_nonperturbed-1) ones are exactly that, and
201 * the remaining ones from nr_nonperturbed..(nr-1) are perturbed bonded
206 /* Returns the total number of elements in iatoms */
218 /* TODO: Replace t_ilist in gmx_localtop_t by InteractionList.
219 * The nr_nonperturbed functionality needs to be ported.
221 * Remove t_ilist and remove templating on list type
222 * in mshift.cpp, constr.cpp, vsite.cpp and domdec_topology.cpp.
226 * The structs InteractionList and t_ilist defines a list of atoms with their interactions.
227 * General field description:
229 * the size (nr elements) of the interactions array (iatoms[]).
231 * specifies which atoms are involved in an interaction of a certain
232 * type. The layout of this array is as follows:
234 * +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
235 * |type1|at1|at2|at3|type2|at1|at2|type1|at1|at2|at3|type3|at1|at2|
236 * +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
238 * So for interaction type type1 3 atoms are needed, and for type2 and
239 * type3 only 2. The type identifier is used to select the function to
240 * calculate the interaction and its actual parameters. This type
241 * identifier is an index in a params[] and functype[] array.
244 /*! \brief Type for returning a list of InteractionList references
246 * TODO: Remove when the function type is made part of InteractionList
248 struct InteractionListHandle
250 const int functionType; //!< The function type
251 const std::vector<int> &iatoms; //!< Reference to interaction list
254 /*! \brief Returns a list of all non-empty InteractionList entries with any of the interaction flags in \p flags set
256 * \param[in] ilists Set of interaction lists
257 * \param[in] flags Bit mask with one or more IF_... bits set
259 static inline const std::vector<InteractionListHandle>
260 extractILists(const InteractionLists &ilists,
263 std::vector<InteractionListHandle> handles;
264 for (size_t ftype = 0; ftype < ilists.size(); ftype++)
266 if ((interaction_function[ftype].flags & flags) && ilists[ftype].size() > 0)
268 handles.push_back({ static_cast<int>(ftype), ilists[ftype].iatoms });
276 real *cmap; /* Has length 4*grid_spacing*grid_spacing, */
277 /* there are 4 entries for each cmap type (V,dVdx,dVdy,d2dVdxdy) */
280 typedef struct gmx_cmap_t
282 int ngrid; /* Number of allocated cmap (cmapdata_t ) grids */
283 int grid_spacing; /* Grid spacing */
284 gmx_cmapdata_t *cmapdata; /* Pointer to grid with actual, pre-interpolated data */
288 /* Struct that holds all force field parameters for the simulated system */
289 struct gmx_ffparams_t
291 /* Returns the number of function types, which matches the number of elements in iparams */
294 GMX_ASSERT(iparams.size() == functype.size(), "Parameters and function types go together");
296 return functype.size();
299 /* TODO: Consider merging functype and iparams, either by storing
300 * the functype in t_iparams or by putting both in a single class.
302 int atnr; /* The number of non-bonded atom types */
303 std::vector<t_functype> functype; /* The function type per type */
304 std::vector<t_iparams> iparams; /* Force field parameters per type */
305 double reppow; /* The repulsion power for VdW: C12*r^-reppow */
306 real fudgeQQ; /* The scaling factor for Coulomb 1-4: f*q1*q2 */
307 gmx_cmap_t cmap_grid; /* The dihedral correction maps */
311 ilsortUNKNOWN, ilsortNO_FE, ilsortFE_UNSORTED, ilsortFE_SORTED
314 typedef struct t_idef
318 t_functype *functype;
321 gmx_cmap_t cmap_grid;
322 t_iparams *iparams_posres, *iparams_fbposres;
323 int iparams_posres_nalloc, iparams_fbposres_nalloc;
330 * The struct t_idef defines all the interactions for the complete
331 * simulation. The structure is setup in such a way that the multinode
332 * version of the program can use it as easy as the single node version.
333 * General field description:
335 * defines the number of elements in functype[] and param[].
337 * the node id (if parallel machines)
339 * the number of atomtypes
340 * t_functype *functype
341 * array of length ntypes, defines for every force type what type of
342 * function to use. Every "bond" with the same function but different
343 * force parameters is a different force type. The type identifier in the
344 * forceatoms[] array is an index in this array.
346 * array of length ntypes, defines the parameters for every interaction
347 * type. The type identifier in the actual interaction list
348 * (ilist[ftype].iatoms[]) is an index in this array.
349 * gmx_cmap_t cmap_grid
350 * the grid for the dihedral pair correction maps.
351 * t_iparams *iparams_posres, *iparams_fbposres
352 * defines the parameters for position restraints only.
353 * Position restraints are the only interactions that have different
354 * parameters (reference positions) for different molecules
355 * of the same type. ilist[F_POSRES].iatoms[] is an index in this array.
357 * The list of interactions for each type. Note that some,
358 * such as LJ and COUL will have 0 entries.
360 * The state of the sorting of il, values are provided above.
363 void pr_iparams(FILE *fp, t_functype ftype, const t_iparams *iparams);
364 void pr_ilist(FILE *fp, int indent, const char *title,
365 const t_functype *functype, const InteractionList &ilist,
366 gmx_bool bShowNumbers,
367 gmx_bool bShowParameters, const t_iparams *iparams);
368 void pr_ffparams(FILE *fp, int indent, const char *title,
369 const gmx_ffparams_t *ffparams, gmx_bool bShowNumbers);
370 void pr_idef(FILE *fp, int indent, const char *title, const t_idef *idef,
371 gmx_bool bShowNumbers, gmx_bool bShowParameters);
374 * Properly initialize idef struct.
376 * \param[in] idef Pointer to idef struct to initialize.
378 void init_idef(t_idef *idef);
381 * Properly clean up idef struct.
383 * \param[in] idef Pointer to idef struct to clean up.
385 void done_idef(t_idef *idef);