Update instructions in containers.rst
[gromacs.git] / src / gromacs / topology / idef.h
blob40a4913cfbb16b15652112801468ec393f35c78f
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) 2013,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_IDEF_H
39 #define GMX_TOPOLOGY_IDEF_H
41 #include <cstdio>
43 #include <array>
44 #include <vector>
46 #include "gromacs/math/vectypes.h"
47 #include "gromacs/topology/ifunc.h"
48 #include "gromacs/utility/basedefinitions.h"
49 #include "gromacs/utility/real.h"
51 struct gmx_ffparams_t;
53 typedef union t_iparams {
54 /* Some parameters have A and B values for free energy calculations.
55 * The B values are not used for regular simulations of course.
56 * Free Energy for nonbondeds can be computed by changing the atom type.
57 * The harmonic type is used for all harmonic potentials:
58 * bonds, angles and improper dihedrals
60 struct
62 real a, b, c;
63 } bham;
64 struct
66 real rA, krA, rB, krB;
67 } harmonic;
68 struct
70 real klinA, aA, klinB, aB;
71 } linangle;
72 struct
74 real lowA, up1A, up2A, kA, lowB, up1B, up2B, kB;
75 } restraint;
76 /* No free energy supported for cubic bonds, FENE, WPOL or cross terms */
77 struct
79 real b0, kb, kcub;
80 } cubic;
81 struct
83 real bm, kb;
84 } fene;
85 struct
87 real r1e, r2e, krr;
88 } cross_bb;
89 struct
91 real r1e, r2e, r3e, krt;
92 } cross_ba;
93 struct
95 real thetaA, kthetaA, r13A, kUBA, thetaB, kthetaB, r13B, kUBB;
96 } u_b;
97 struct
99 real theta, c[5];
100 } qangle;
101 struct
103 real alpha;
104 } polarize;
105 struct
107 real alpha, drcut, khyp;
108 } anharm_polarize;
109 struct
111 real al_x, al_y, al_z, rOH, rHH, rOD;
112 } wpol;
113 struct
115 real a, alpha1, alpha2, rfac;
116 } thole;
117 struct
119 real c6, c12;
120 } lj;
121 struct
123 real c6A, c12A, c6B, c12B;
124 } lj14;
125 struct
127 real fqq, qi, qj, c6, c12;
128 } ljc14;
129 struct
131 real qi, qj, c6, c12;
132 } ljcnb;
133 /* Proper dihedrals can not have different multiplicity when
134 * doing free energy calculations, because the potential would not
135 * be periodic anymore.
137 struct
139 real phiA, cpA;
140 int mult;
141 real phiB, cpB;
142 } pdihs;
143 struct
145 real dA, dB;
146 } constr;
147 /* Settle can not be used for Free energy calculations of water bond geometry.
148 * Use shake (or lincs) instead if you have to change the water bonds.
150 struct
152 real doh, dhh;
153 } settle;
154 struct
156 real b0A, cbA, betaA, b0B, cbB, betaB;
157 } morse;
158 struct
160 real pos0A[DIM], fcA[DIM], pos0B[DIM], fcB[DIM];
161 } posres;
162 struct
164 real pos0[DIM], r, k;
165 int geom;
166 } fbposres;
167 struct
169 real rbcA[NR_RBDIHS], rbcB[NR_RBDIHS];
170 } rbdihs;
171 struct
173 real cbtcA[NR_CBTDIHS], cbtcB[NR_CBTDIHS];
174 } cbtdihs;
175 struct
177 real a, b, c, d, e, f;
178 } vsite;
179 struct
181 int n;
182 real a;
183 } vsiten;
184 /* NOTE: npair is only set after reading the tpx file */
185 struct
187 real low, up1, up2, kfac;
188 int type, label, npair;
189 } disres;
190 struct
192 real phiA, dphiA, kfacA, phiB, dphiB, kfacB;
193 } dihres;
194 struct
196 int ex, power, label;
197 real c, obs, kfac;
198 } orires;
199 struct
201 int table;
202 real kA;
203 real kB;
204 } tab;
205 struct
207 int cmapA, cmapB;
208 } cmap;
209 struct
211 real buf[MAXFORCEPARAM];
212 } generic; /* Conversion */
213 } t_iparams;
215 typedef int t_functype;
217 /* List of listed interactions, see description further down.
219 * TODO: Consider storing the function type as well.
220 * TODO: Consider providing per interaction access.
222 struct InteractionList
224 /* Returns the total number of elements in iatoms */
225 int size() const { return static_cast<int>(iatoms.size()); }
227 /* Returns whether the list is empty */
228 bool empty() const { return iatoms.empty(); }
230 /* Adds one interaction to the list */
231 template<std::size_t numAtoms>
232 void push_back(const int parameterType, const std::array<int, numAtoms>& atoms)
234 const std::size_t oldSize = iatoms.size();
235 iatoms.resize(iatoms.size() + 1 + numAtoms);
236 iatoms[oldSize] = parameterType;
237 for (std::size_t i = 0; i < numAtoms; i++)
239 iatoms[oldSize + 1 + i] = atoms[i];
243 /* Adds one interaction to the list */
244 void push_back(const int parameterType, const int numAtoms, const int* atoms)
246 const std::size_t oldSize = iatoms.size();
247 iatoms.resize(iatoms.size() + 1 + numAtoms);
248 iatoms[oldSize] = parameterType;
249 for (int i = 0; i < numAtoms; i++)
251 iatoms[oldSize + 1 + i] = atoms[i];
255 /* Appends \p ilist at the back of the list */
256 void append(const InteractionList& ilist)
258 iatoms.insert(iatoms.end(), ilist.iatoms.begin(), ilist.iatoms.end());
261 /* Clears the list */
262 void clear() { iatoms.clear(); }
264 /* List of interactions, see explanation further down */
265 std::vector<int> iatoms;
268 /* List of interaction lists, one list for each interaction type
270 * TODO: Consider only including entries in use instead of all F_NRE
272 using InteractionLists = std::array<InteractionList, F_NRE>;
274 /* Deprecated list of listed interactions */
275 struct t_ilist
277 /* Returns the total number of elements in iatoms */
278 int size() const { return nr; }
280 /* Returns whether the list is empty */
281 bool empty() const { return nr == 0; }
283 int nr;
284 t_iatom* iatoms;
285 int nalloc;
288 /* TODO: Remove t_ilist and remove templating on list type in mshift.cpp */
291 * The structs InteractionList and t_ilist defines a list of atoms with their interactions.
292 * General field description:
293 * int nr
294 * the size (nr elements) of the interactions array (iatoms[]).
295 * t_iatom *iatoms
296 * specifies which atoms are involved in an interaction of a certain
297 * type. The layout of this array is as follows:
299 * +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
300 * |type1|at1|at2|at3|type2|at1|at2|type1|at1|at2|at3|type3|at1|at2|
301 * +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
303 * So for interaction type type1 3 atoms are needed, and for type2 and
304 * type3 only 2. The type identifier is used to select the function to
305 * calculate the interaction and its actual parameters. This type
306 * identifier is an index in a params[] and functype[] array.
309 /*! \brief Type for returning a list of InteractionList references
311 * TODO: Remove when the function type is made part of InteractionList
313 struct InteractionListHandle
315 const int functionType; //!< The function type
316 const std::vector<int>& iatoms; //!< Reference to interaction list
319 /*! \brief Returns a list of all non-empty InteractionList entries with any of the interaction flags in \p flags set
321 * \param[in] ilists Set of interaction lists
322 * \param[in] flags Bit mask with one or more IF_... bits set
324 static inline std::vector<InteractionListHandle> extractILists(const InteractionLists& ilists, int flags)
326 std::vector<InteractionListHandle> handles;
327 for (size_t ftype = 0; ftype < ilists.size(); ftype++)
329 if ((interaction_function[ftype].flags & flags) && !ilists[ftype].empty())
331 handles.push_back({ static_cast<int>(ftype), ilists[ftype].iatoms });
334 return handles;
337 /*! \brief Returns the stride for the iatoms array in \p ilistHandle
339 * \param[in] ilistHandle The ilist to return the stride for
341 static inline int ilistStride(const InteractionListHandle& ilistHandle)
343 return 1 + NRAL(ilistHandle.functionType);
346 struct gmx_cmapdata_t
348 std::vector<real> cmap; /* Has length 4*grid_spacing*grid_spacing, */
349 /* there are 4 entries for each cmap type (V,dVdx,dVdy,d2dVdxdy) */
352 struct gmx_cmap_t
354 int grid_spacing = 0; /* Grid spacing */
355 std::vector<gmx_cmapdata_t> cmapdata; /* Lists of grids with actual, pre-interpolated data */
359 enum
361 ilsortUNKNOWN,
362 ilsortNO_FE,
363 ilsortFE_UNSORTED,
364 ilsortFE_SORTED
367 /* Struct with list of interaction parameters and lists of interactions
369 * TODO: Convert to a proper class with private data members so we can
370 * ensure that the free-energy sorting and sorting setting is consistent.
372 class InteractionDefinitions
374 public:
375 /* Constructor
377 * \param[in] ffparams The interaction parameters, the lifetime of the created object should not exceed the lifetime of the passed parameters
379 InteractionDefinitions(const gmx_ffparams_t& ffparams);
381 // Clears data not read in from ffparams
382 void clear();
384 // The interaction parameters
385 const std::vector<t_iparams>& iparams;
386 // The function type per type
387 const std::vector<int>& functype;
388 // Position restraint interaction parameters
389 std::vector<t_iparams> iparams_posres;
390 // Flat-bottomed position restraint parameters
391 std::vector<t_iparams> iparams_fbposres;
392 // The list of interactions for each type. Note that some, such as LJ and COUL will have 0 entries.
393 std::array<InteractionList, F_NRE> il;
394 /* The number of non-perturbed interactions at the start of each entry in il */
395 std::array<int, F_NRE> numNonperturbedInteractions;
396 // The sorting state of interaction in il
397 int ilsort = ilsortUNKNOWN;
398 // The dihedral correction maps
399 gmx_cmap_t cmap_grid;
402 /* Deprecated interation definitions, used in t_topology */
403 struct t_idef
405 int ntypes;
406 int atnr;
407 t_functype* functype;
408 t_iparams* iparams;
409 real fudgeQQ;
410 t_iparams * iparams_posres, *iparams_fbposres;
412 t_ilist il[F_NRE];
413 int ilsort;
417 * The struct t_idef defines all the interactions for the complete
418 * simulation. The structure is setup in such a way that the multinode
419 * version of the program can use it as easy as the single node version.
420 * General field description:
421 * int ntypes
422 * defines the number of elements in functype[] and param[].
423 * int nodeid
424 * the node id (if parallel machines)
425 * int atnr
426 * the number of atomtypes
427 * t_functype *functype
428 * array of length ntypes, defines for every force type what type of
429 * function to use. Every "bond" with the same function but different
430 * force parameters is a different force type. The type identifier in the
431 * forceatoms[] array is an index in this array.
432 * t_iparams *iparams
433 * array of length ntypes, defines the parameters for every interaction
434 * type. The type identifier in the actual interaction list
435 * (ilist[ftype].iatoms[]) is an index in this array.
436 * gmx_cmap_t cmap_grid
437 * the grid for the dihedral pair correction maps.
438 * t_iparams *iparams_posres, *iparams_fbposres
439 * defines the parameters for position restraints only.
440 * Position restraints are the only interactions that have different
441 * parameters (reference positions) for different molecules
442 * of the same type. ilist[F_POSRES].iatoms[] is an index in this array.
443 * t_ilist il[F_NRE]
444 * The list of interactions for each type. Note that some,
445 * such as LJ and COUL will have 0 entries.
446 * int ilsort
447 * The state of the sorting of il, values are provided above.
450 namespace gmx
452 class TextWriter;
453 } // namespace gmx
455 void printInteractionParameters(gmx::TextWriter* writer, t_functype ftype, const t_iparams& iparams);
456 void pr_iparams(FILE* fp, t_functype ftype, const t_iparams& iparams);
457 void pr_ilist(FILE* fp,
458 int indent,
459 const char* title,
460 const t_functype* functype,
461 const InteractionList& ilist,
462 gmx_bool bShowNumbers,
463 gmx_bool bShowParameters,
464 const t_iparams* iparams);
465 void pr_idef(FILE* fp, int indent, const char* title, const t_idef* idef, gmx_bool bShowNumbers, gmx_bool bShowParameters);
467 /*! \brief
468 * Properly initialize idef struct.
470 * \param[in] idef Pointer to idef struct to initialize.
472 void init_idef(t_idef* idef);
474 /*! \brief
475 * Properly clean up idef struct.
477 * \param[in] idef Pointer to idef struct to clean up.
479 void done_idef(t_idef* idef);
481 #endif