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