Merge branch release-2018 into release-2019
[alexxy/gromacs.git] / src / gromacs / topology / idef.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #ifndef GMX_TOPOLOGY_IDEF_H
38 #define GMX_TOPOLOGY_IDEF_H
39
40 #include <cstdio>
41
42 #include <array>
43 #include <vector>
44
45 #include "gromacs/math/vectypes.h"
46 #include "gromacs/topology/ifunc.h"
47 #include "gromacs/utility/basedefinitions.h"
48 #include "gromacs/utility/real.h"
49
50 typedef union t_iparams
51 {
52     /* Some parameters have A and B values for free energy calculations.
53      * The B values are not used for regular simulations of course.
54      * Free Energy for nonbondeds can be computed by changing the atom type.
55      * The harmonic type is used for all harmonic potentials:
56      * bonds, angles and improper dihedrals
57      */
58     struct {
59         real a, b, c;
60     } bham;
61     struct {
62         real rA, krA, rB, krB;
63     } harmonic;
64     struct {
65         real klinA, aA, klinB, aB;
66     } linangle;
67     struct {
68         real lowA, up1A, up2A, kA, lowB, up1B, up2B, kB;
69     } restraint;
70     /* No free energy supported for cubic bonds, FENE, WPOL or cross terms */
71     struct {
72         real b0, kb, kcub;
73     } cubic;
74     struct {
75         real bm, kb;
76     } fene;
77     struct {
78         real r1e, r2e, krr;
79     } cross_bb;
80     struct {
81         real r1e, r2e, r3e, krt;
82     } cross_ba;
83     struct {
84         real thetaA, kthetaA, r13A, kUBA, thetaB, kthetaB, r13B, kUBB;
85     } u_b;
86     struct {
87         real theta, c[5];
88     } qangle;
89     struct {
90         real alpha;
91     } polarize;
92     struct {
93         real alpha, drcut, khyp;
94     } anharm_polarize;
95     struct {
96         real al_x, al_y, al_z, rOH, rHH, rOD;
97     } wpol;
98     struct {
99         real a, alpha1, alpha2, rfac;
100     } thole;
101     struct {
102         real c6, c12;
103     } lj;
104     struct {
105         real c6A, c12A, c6B, c12B;
106     } lj14;
107     struct {
108         real fqq, qi, qj, c6, c12;
109     } ljc14;
110     struct {
111         real qi, qj, c6, c12;
112     } ljcnb;
113     /* Proper dihedrals can not have different multiplicity when
114      * doing free energy calculations, because the potential would not
115      * be periodic anymore.
116      */
117     struct {
118         real phiA, cpA; int mult; real phiB, cpB;
119     } pdihs;
120     struct {
121         real dA, dB;
122     } constr;
123     /* Settle can not be used for Free energy calculations of water bond geometry.
124      * Use shake (or lincs) instead if you have to change the water bonds.
125      */
126     struct {
127         real doh, dhh;
128     } settle;
129     struct {
130         real b0A, cbA, betaA, b0B, cbB, betaB;
131     } morse;
132     struct {
133         real pos0A[DIM], fcA[DIM], pos0B[DIM], fcB[DIM];
134     } posres;
135     struct {
136         real pos0[DIM], r, k; int geom;
137     } fbposres;
138     struct {
139         real rbcA[NR_RBDIHS], rbcB[NR_RBDIHS];
140     } rbdihs;
141     struct {
142         real cbtcA[NR_CBTDIHS], cbtcB[NR_CBTDIHS];
143     } cbtdihs;
144     struct {
145         real a, b, c, d, e, f;
146     } vsite;
147     struct {
148         int  n; real a;
149     } vsiten;
150     /* NOTE: npair is only set after reading the tpx file */
151     struct {
152         real low, up1, up2, kfac; int type, label, npair;
153     } disres;
154     struct {
155         real phiA, dphiA, kfacA, phiB, dphiB, kfacB;
156     } dihres;
157     struct {
158         int  ex, power, label; real c, obs, kfac;
159     } orires;
160     struct {
161         int  table; real kA; real kB;
162     } tab;
163     struct {
164         int cmapA, cmapB;
165     } cmap;
166     struct {
167         real buf[MAXFORCEPARAM];
168     } generic;                                               /* Conversion */
169 } t_iparams;
170
171 typedef int t_functype;
172
173 /* List of listed interactions, see description further down.
174  *
175  * TODO: Consider storing the function type as well.
176  * TODO: Consider providing per interaction access.
177  */
178 struct InteractionList
179 {
180     /* Returns the total number of elements in iatoms */
181     int size() const
182     {
183         return iatoms.size();
184     }
185
186     /* List of interactions, see explanation further down */
187     std::vector<int> iatoms;
188 };
189
190 /* List of interaction lists, one list for each interaction type
191  *
192  * TODO: Consider only including entries in use instead of all F_NRE
193  */
194 typedef std::array<InteractionList, F_NRE> InteractionLists;
195
196 /* Deprecated list of listed interactions.
197  *
198  * The nonperturbed/perturbed interactions are now separated (sorted) in the
199  * ilist, such that the first 0..(nr_nonperturbed-1) ones are exactly that, and
200  * the remaining ones from nr_nonperturbed..(nr-1) are perturbed bonded
201  * interactions.
202  */
203 struct t_ilist
204 {
205     /* Returns the total number of elements in iatoms */
206     int size() const
207     {
208         return nr;
209     }
210
211     int      nr;
212     int      nr_nonperturbed;
213     t_iatom *iatoms;
214     int      nalloc;
215 };
216
217 /* TODO: Replace t_ilist in gmx_localtop_t by InteractionList.
218  *       The nr_nonperturbed functionality needs to be ported.
219  *       Remove t_topology.
220  *       Remove t_ilist and remove templating on list type
221  *       in mshift.cpp, constr.cpp, vsite.cpp and domdec_topology.cpp.
222  */
223
224 /*
225  * The structs InteractionList and t_ilist defines a list of atoms with their interactions.
226  * General field description:
227  *   int nr
228  *      the size (nr elements) of the interactions array (iatoms[]).
229  *   t_iatom *iatoms
230  *  specifies which atoms are involved in an interaction of a certain
231  *       type. The layout of this array is as follows:
232  *
233  *        +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
234  *        |type1|at1|at2|at3|type2|at1|at2|type1|at1|at2|at3|type3|at1|at2|
235  *        +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
236  *
237  *  So for interaction type type1 3 atoms are needed, and for type2 and
238  *      type3 only 2. The type identifier is used to select the function to
239  *      calculate the interaction and its actual parameters. This type
240  *      identifier is an index in a params[] and functype[] array.
241  */
242
243 /*! \brief Type for returning a list of InteractionList references
244  *
245  * TODO: Remove when the function type is made part of InteractionList
246  */
247 struct InteractionListHandle
248 {
249     const int               functionType; //!< The function type
250     const std::vector<int> &iatoms;       //!< Reference to interaction list
251 };
252
253 /*! \brief Returns a list of all non-empty InteractionList entries with any of the interaction flags in \p flags set
254  *
255  * \param[in] ilists  Set of interaction lists
256  * \param[in] flags   Bit mask with one or more IF_... bits set
257  */
258 static inline const std::vector<InteractionListHandle>
259 extractILists(const InteractionLists &ilists,
260               int                     flags)
261 {
262     std::vector<InteractionListHandle> handles;
263     for (size_t ftype = 0; ftype < ilists.size(); ftype++)
264     {
265         if ((interaction_function[ftype].flags & flags) && ilists[ftype].size() > 0)
266         {
267             handles.push_back({ static_cast<int>(ftype), ilists[ftype].iatoms });
268         }
269     }
270     return handles;
271 }
272
273 /*! \brief Returns the stride for the iatoms array in \p ilistHandle
274  *
275  * \param[in] ilistHandle  The ilist to return the stride for
276  */
277 static inline int ilistStride(const InteractionListHandle &ilistHandle)
278 {
279     return 1 + NRAL(ilistHandle.functionType);
280 }
281
282 struct gmx_cmapdata_t
283 {
284     std::vector<real> cmap; /* Has length 4*grid_spacing*grid_spacing, */
285     /* there are 4 entries for each cmap type (V,dVdx,dVdy,d2dVdxdy) */
286 };
287
288 struct gmx_cmap_t
289 {
290     int                         grid_spacing = 0; /* Grid spacing */
291     std::vector<gmx_cmapdata_t> cmapdata;         /* Lists of grids with actual, pre-interpolated data */
292 };
293
294
295 enum {
296     ilsortUNKNOWN, ilsortNO_FE, ilsortFE_UNSORTED, ilsortFE_SORTED
297 };
298
299 typedef struct t_idef
300 {
301     int         ntypes;
302     int         atnr;
303     t_functype *functype;
304     t_iparams  *iparams;
305     real        fudgeQQ;
306     gmx_cmap_t *cmap_grid;
307     t_iparams  *iparams_posres, *iparams_fbposres;
308     int         iparams_posres_nalloc, iparams_fbposres_nalloc;
309
310     t_ilist     il[F_NRE];
311     int         ilsort;
312 } t_idef;
313
314 /*
315  * The struct t_idef defines all the interactions for the complete
316  * simulation. The structure is setup in such a way that the multinode
317  * version of the program  can use it as easy as the single node version.
318  * General field description:
319  *   int ntypes
320  *      defines the number of elements in functype[] and param[].
321  *   int nodeid
322  *      the node id (if parallel machines)
323  *   int atnr
324  *      the number of atomtypes
325  *   t_functype *functype
326  *      array of length ntypes, defines for every force type what type of
327  *      function to use. Every "bond" with the same function but different
328  *      force parameters is a different force type. The type identifier in the
329  *      forceatoms[] array is an index in this array.
330  *   t_iparams *iparams
331  *      array of length ntypes, defines the parameters for every interaction
332  *      type. The type identifier in the actual interaction list
333  *      (ilist[ftype].iatoms[]) is an index in this array.
334  *   gmx_cmap_t cmap_grid
335  *      the grid for the dihedral pair correction maps.
336  *   t_iparams *iparams_posres, *iparams_fbposres
337  *      defines the parameters for position restraints only.
338  *      Position restraints are the only interactions that have different
339  *      parameters (reference positions) for different molecules
340  *      of the same type. ilist[F_POSRES].iatoms[] is an index in this array.
341  *   t_ilist il[F_NRE]
342  *      The list of interactions for each type. Note that some,
343  *      such as LJ and COUL will have 0 entries.
344  *   int ilsort
345  *      The state of the sorting of il, values are provided above.
346  */
347
348 void pr_iparams(FILE *fp, t_functype ftype, const t_iparams *iparams);
349 void pr_ilist(FILE *fp, int indent, const char *title,
350               const t_functype *functype, const InteractionList &ilist,
351               gmx_bool bShowNumbers,
352               gmx_bool bShowParameters, const t_iparams *iparams);
353 void pr_idef(FILE *fp, int indent, const char *title, const t_idef *idef,
354              gmx_bool bShowNumbers, gmx_bool bShowParameters);
355
356 /*! \brief
357  * Properly initialize idef struct.
358  *
359  * \param[in] idef Pointer to idef struct to initialize.
360  */
361 void init_idef(t_idef *idef);
362
363 /*! \brief
364  * Properly clean up idef struct.
365  *
366  * \param[in] idef Pointer to idef struct to clean up.
367  */
368 void done_idef(t_idef *idef);
369
370 void copy_ilist(const t_ilist *src, t_ilist *dst);
371
372 #endif