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