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,2019, 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/real.h"
50 typedef union t_iparams
{
51 /* Some parameters have A and B values for free energy calculations.
52 * The B values are not used for regular simulations of course.
53 * Free Energy for nonbondeds can be computed by changing the atom type.
54 * The harmonic type is used for all harmonic potentials:
55 * bonds, angles and improper dihedrals
63 real rA
, krA
, rB
, krB
;
67 real klinA
, aA
, klinB
, aB
;
71 real lowA
, up1A
, up2A
, kA
, lowB
, up1B
, up2B
, kB
;
73 /* No free energy supported for cubic bonds, FENE, WPOL or cross terms */
88 real r1e
, r2e
, r3e
, krt
;
92 real thetaA
, kthetaA
, r13A
, kUBA
, thetaB
, kthetaB
, r13B
, kUBB
;
104 real alpha
, drcut
, khyp
;
108 real al_x
, al_y
, al_z
, rOH
, rHH
, rOD
;
112 real a
, alpha1
, alpha2
, rfac
;
120 real c6A
, c12A
, c6B
, c12B
;
124 real fqq
, qi
, qj
, c6
, c12
;
128 real qi
, qj
, c6
, c12
;
130 /* Proper dihedrals can not have different multiplicity when
131 * doing free energy calculations, because the potential would not
132 * be periodic anymore.
144 /* Settle can not be used for Free energy calculations of water bond geometry.
145 * Use shake (or lincs) instead if you have to change the water bonds.
153 real b0A
, cbA
, betaA
, b0B
, cbB
, betaB
;
157 real pos0A
[DIM
], fcA
[DIM
], pos0B
[DIM
], fcB
[DIM
];
161 real pos0
[DIM
], r
, k
;
166 real rbcA
[NR_RBDIHS
], rbcB
[NR_RBDIHS
];
170 real cbtcA
[NR_CBTDIHS
], cbtcB
[NR_CBTDIHS
];
174 real a
, b
, c
, d
, e
, f
;
181 /* NOTE: npair is only set after reading the tpx file */
184 real low
, up1
, up2
, kfac
;
185 int type
, label
, npair
;
189 real phiA
, dphiA
, kfacA
, phiB
, dphiB
, kfacB
;
193 int ex
, power
, label
;
208 real buf
[MAXFORCEPARAM
];
209 } generic
; /* Conversion */
212 typedef int t_functype
;
214 /* List of listed interactions, see description further down.
216 * TODO: Consider storing the function type as well.
217 * TODO: Consider providing per interaction access.
219 struct InteractionList
221 /* Returns the total number of elements in iatoms */
222 int size() const { return gmx::ssize(iatoms
); }
224 /* List of interactions, see explanation further down */
225 std::vector
<int> iatoms
;
228 /* List of interaction lists, one list for each interaction type
230 * TODO: Consider only including entries in use instead of all F_NRE
232 typedef std::array
<InteractionList
, F_NRE
> InteractionLists
;
234 /* Deprecated list of listed interactions.
236 * The nonperturbed/perturbed interactions are now separated (sorted) in the
237 * ilist, such that the first 0..(nr_nonperturbed-1) ones are exactly that, and
238 * the remaining ones from nr_nonperturbed..(nr-1) are perturbed bonded
243 /* Returns the total number of elements in iatoms */
244 int size() const { return nr
; }
251 /* TODO: Replace t_ilist in gmx_localtop_t by InteractionList.
252 * The nr_nonperturbed functionality needs to be ported.
254 * Remove t_ilist and remove templating on list type
255 * in mshift.cpp, constr.cpp, vsite.cpp and domdec_topology.cpp.
259 * The structs InteractionList and t_ilist defines a list of atoms with their interactions.
260 * General field description:
262 * the size (nr elements) of the interactions array (iatoms[]).
264 * specifies which atoms are involved in an interaction of a certain
265 * type. The layout of this array is as follows:
267 * +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
268 * |type1|at1|at2|at3|type2|at1|at2|type1|at1|at2|at3|type3|at1|at2|
269 * +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
271 * So for interaction type type1 3 atoms are needed, and for type2 and
272 * type3 only 2. The type identifier is used to select the function to
273 * calculate the interaction and its actual parameters. This type
274 * identifier is an index in a params[] and functype[] array.
277 /*! \brief Type for returning a list of InteractionList references
279 * TODO: Remove when the function type is made part of InteractionList
281 struct InteractionListHandle
283 const int functionType
; //!< The function type
284 const std::vector
<int>& iatoms
; //!< Reference to interaction list
287 /*! \brief Returns a list of all non-empty InteractionList entries with any of the interaction flags in \p flags set
289 * \param[in] ilists Set of interaction lists
290 * \param[in] flags Bit mask with one or more IF_... bits set
292 static inline std::vector
<InteractionListHandle
> extractILists(const InteractionLists
& ilists
, int flags
)
294 std::vector
<InteractionListHandle
> handles
;
295 for (size_t ftype
= 0; ftype
< ilists
.size(); ftype
++)
297 if ((interaction_function
[ftype
].flags
& flags
) && ilists
[ftype
].size() > 0)
299 handles
.push_back({ static_cast<int>(ftype
), ilists
[ftype
].iatoms
});
305 /*! \brief Returns the stride for the iatoms array in \p ilistHandle
307 * \param[in] ilistHandle The ilist to return the stride for
309 static inline int ilistStride(const InteractionListHandle
& ilistHandle
)
311 return 1 + NRAL(ilistHandle
.functionType
);
314 struct gmx_cmapdata_t
316 std::vector
<real
> cmap
; /* Has length 4*grid_spacing*grid_spacing, */
317 /* there are 4 entries for each cmap type (V,dVdx,dVdy,d2dVdxdy) */
322 int grid_spacing
= 0; /* Grid spacing */
323 std::vector
<gmx_cmapdata_t
> cmapdata
; /* Lists of grids with actual, pre-interpolated data */
335 typedef struct t_idef
339 t_functype
* functype
;
342 gmx_cmap_t
* cmap_grid
;
343 t_iparams
* iparams_posres
, *iparams_fbposres
;
344 int iparams_posres_nalloc
, iparams_fbposres_nalloc
;
347 /* The number of non-perturbed interactions at the start of each entry in il */
348 int numNonperturbedInteractions
[F_NRE
];
353 * The struct t_idef defines all the interactions for the complete
354 * simulation. The structure is setup in such a way that the multinode
355 * version of the program can use it as easy as the single node version.
356 * General field description:
358 * defines the number of elements in functype[] and param[].
360 * the node id (if parallel machines)
362 * the number of atomtypes
363 * t_functype *functype
364 * array of length ntypes, defines for every force type what type of
365 * function to use. Every "bond" with the same function but different
366 * force parameters is a different force type. The type identifier in the
367 * forceatoms[] array is an index in this array.
369 * array of length ntypes, defines the parameters for every interaction
370 * type. The type identifier in the actual interaction list
371 * (ilist[ftype].iatoms[]) is an index in this array.
372 * gmx_cmap_t cmap_grid
373 * the grid for the dihedral pair correction maps.
374 * t_iparams *iparams_posres, *iparams_fbposres
375 * defines the parameters for position restraints only.
376 * Position restraints are the only interactions that have different
377 * parameters (reference positions) for different molecules
378 * of the same type. ilist[F_POSRES].iatoms[] is an index in this array.
380 * The list of interactions for each type. Note that some,
381 * such as LJ and COUL will have 0 entries.
383 * The state of the sorting of il, values are provided above.
391 void printInteractionParameters(gmx::TextWriter
* writer
, t_functype ftype
, const t_iparams
& iparams
);
392 void pr_iparams(FILE* fp
, t_functype ftype
, const t_iparams
& iparams
);
393 void pr_ilist(FILE* fp
,
396 const t_functype
* functype
,
397 const InteractionList
& ilist
,
398 gmx_bool bShowNumbers
,
399 gmx_bool bShowParameters
,
400 const t_iparams
* iparams
);
401 void pr_idef(FILE* fp
, int indent
, const char* title
, const t_idef
* idef
, gmx_bool bShowNumbers
, gmx_bool bShowParameters
);
404 * Properly initialize idef struct.
406 * \param[in] idef Pointer to idef struct to initialize.
408 void init_idef(t_idef
* idef
);
411 * Properly clean up idef struct.
413 * \param[in] idef Pointer to idef struct to clean up.
415 void done_idef(t_idef
* idef
);
417 void copy_ilist(const t_ilist
* src
, t_ilist
* dst
);