Move nr_nonperturbed out of t_ilist
[gromacs.git] / src / gromacs / topology / idef.h
blob269d531fc0380b7a7fa022e4f0ff0cf22ecee8c9
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,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
40 #include <cstdio>
42 #include <array>
43 #include <vector>
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
57 struct
59 real a, b, c;
60 } bham;
61 struct
63 real rA, krA, rB, krB;
64 } harmonic;
65 struct
67 real klinA, aA, klinB, aB;
68 } linangle;
69 struct
71 real lowA, up1A, up2A, kA, lowB, up1B, up2B, kB;
72 } restraint;
73 /* No free energy supported for cubic bonds, FENE, WPOL or cross terms */
74 struct
76 real b0, kb, kcub;
77 } cubic;
78 struct
80 real bm, kb;
81 } fene;
82 struct
84 real r1e, r2e, krr;
85 } cross_bb;
86 struct
88 real r1e, r2e, r3e, krt;
89 } cross_ba;
90 struct
92 real thetaA, kthetaA, r13A, kUBA, thetaB, kthetaB, r13B, kUBB;
93 } u_b;
94 struct
96 real theta, c[5];
97 } qangle;
98 struct
100 real alpha;
101 } polarize;
102 struct
104 real alpha, drcut, khyp;
105 } anharm_polarize;
106 struct
108 real al_x, al_y, al_z, rOH, rHH, rOD;
109 } wpol;
110 struct
112 real a, alpha1, alpha2, rfac;
113 } thole;
114 struct
116 real c6, c12;
117 } lj;
118 struct
120 real c6A, c12A, c6B, c12B;
121 } lj14;
122 struct
124 real fqq, qi, qj, c6, c12;
125 } ljc14;
126 struct
128 real qi, qj, c6, c12;
129 } ljcnb;
130 /* Proper dihedrals can not have different multiplicity when
131 * doing free energy calculations, because the potential would not
132 * be periodic anymore.
134 struct
136 real phiA, cpA;
137 int mult;
138 real phiB, cpB;
139 } pdihs;
140 struct
142 real dA, dB;
143 } constr;
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.
147 struct
149 real doh, dhh;
150 } settle;
151 struct
153 real b0A, cbA, betaA, b0B, cbB, betaB;
154 } morse;
155 struct
157 real pos0A[DIM], fcA[DIM], pos0B[DIM], fcB[DIM];
158 } posres;
159 struct
161 real pos0[DIM], r, k;
162 int geom;
163 } fbposres;
164 struct
166 real rbcA[NR_RBDIHS], rbcB[NR_RBDIHS];
167 } rbdihs;
168 struct
170 real cbtcA[NR_CBTDIHS], cbtcB[NR_CBTDIHS];
171 } cbtdihs;
172 struct
174 real a, b, c, d, e, f;
175 } vsite;
176 struct
178 int n;
179 real a;
180 } vsiten;
181 /* NOTE: npair is only set after reading the tpx file */
182 struct
184 real low, up1, up2, kfac;
185 int type, label, npair;
186 } disres;
187 struct
189 real phiA, dphiA, kfacA, phiB, dphiB, kfacB;
190 } dihres;
191 struct
193 int ex, power, label;
194 real c, obs, kfac;
195 } orires;
196 struct
198 int table;
199 real kA;
200 real kB;
201 } tab;
202 struct
204 int cmapA, cmapB;
205 } cmap;
206 struct
208 real buf[MAXFORCEPARAM];
209 } generic; /* Conversion */
210 } t_iparams;
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
239 * interactions.
241 struct t_ilist
243 /* Returns the total number of elements in iatoms */
244 int size() const { return nr; }
246 int nr;
247 t_iatom* iatoms;
248 int nalloc;
251 /* TODO: Replace t_ilist in gmx_localtop_t by InteractionList.
252 * The nr_nonperturbed functionality needs to be ported.
253 * Remove t_topology.
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:
261 * int nr
262 * the size (nr elements) of the interactions array (iatoms[]).
263 * t_iatom *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 });
302 return handles;
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) */
320 struct gmx_cmap_t
322 int grid_spacing = 0; /* Grid spacing */
323 std::vector<gmx_cmapdata_t> cmapdata; /* Lists of grids with actual, pre-interpolated data */
327 enum
329 ilsortUNKNOWN,
330 ilsortNO_FE,
331 ilsortFE_UNSORTED,
332 ilsortFE_SORTED
335 typedef struct t_idef
337 int ntypes;
338 int atnr;
339 t_functype* functype;
340 t_iparams* iparams;
341 real fudgeQQ;
342 gmx_cmap_t* cmap_grid;
343 t_iparams * iparams_posres, *iparams_fbposres;
344 int iparams_posres_nalloc, iparams_fbposres_nalloc;
346 t_ilist il[F_NRE];
347 /* The number of non-perturbed interactions at the start of each entry in il */
348 int numNonperturbedInteractions[F_NRE];
349 int ilsort;
350 } t_idef;
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:
357 * int ntypes
358 * defines the number of elements in functype[] and param[].
359 * int nodeid
360 * the node id (if parallel machines)
361 * int atnr
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.
368 * t_iparams *iparams
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.
379 * t_ilist il[F_NRE]
380 * The list of interactions for each type. Note that some,
381 * such as LJ and COUL will have 0 entries.
382 * int ilsort
383 * The state of the sorting of il, values are provided above.
386 namespace gmx
388 class TextWriter;
389 } // namespace gmx
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,
394 int indent,
395 const char* title,
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);
403 /*! \brief
404 * Properly initialize idef struct.
406 * \param[in] idef Pointer to idef struct to initialize.
408 void init_idef(t_idef* idef);
410 /*! \brief
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);
419 #endif