Add C++ version of t_ilist
[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/utility/basedefinitions.h"
47 #include "gromacs/utility/real.h"
48
49 /* check kernel/toppush.c when you change these numbers */
50 #define MAXATOMLIST 6
51 #define MAXFORCEPARAM   12
52 #define NR_RBDIHS   6
53 #define NR_CBTDIHS   6
54 #define NR_FOURDIHS     4
55
56 typedef int t_iatom;
57
58 /* this MUST correspond to the
59    t_interaction_function[F_NRE] in gmxlib/ifunc.c */
60 enum {
61     F_BONDS,
62     F_G96BONDS,
63     F_MORSE,
64     F_CUBICBONDS,
65     F_CONNBONDS,
66     F_HARMONIC,
67     F_FENEBONDS,
68     F_TABBONDS,
69     F_TABBONDSNC,
70     F_RESTRBONDS,
71     F_ANGLES,
72     F_G96ANGLES,
73     F_RESTRANGLES,
74     F_LINEAR_ANGLES,
75     F_CROSS_BOND_BONDS,
76     F_CROSS_BOND_ANGLES,
77     F_UREY_BRADLEY,
78     F_QUARTIC_ANGLES,
79     F_TABANGLES,
80     F_PDIHS,
81     F_RBDIHS,
82     F_RESTRDIHS,
83     F_CBTDIHS,
84     F_FOURDIHS,
85     F_IDIHS,
86     F_PIDIHS,
87     F_TABDIHS,
88     F_CMAP,
89     F_GB12_NOLONGERUSED,
90     F_GB13_NOLONGERUSED,
91     F_GB14_NOLONGERUSED,
92     F_GBPOL_NOLONGERUSED,
93     F_NPSOLVATION_NOLONGERUSED,
94     F_LJ14,
95     F_COUL14,
96     F_LJC14_Q,
97     F_LJC_PAIRS_NB,
98     F_LJ,
99     F_BHAM,
100     F_LJ_LR_NOLONGERUSED,
101     F_BHAM_LR_NOLONGERUSED,
102     F_DISPCORR,
103     F_COUL_SR,
104     F_COUL_LR_NOLONGERUSED,
105     F_RF_EXCL,
106     F_COUL_RECIP,
107     F_LJ_RECIP,
108     F_DPD,
109     F_POLARIZATION,
110     F_WATER_POL,
111     F_THOLE_POL,
112     F_ANHARM_POL,
113     F_POSRES,
114     F_FBPOSRES,
115     F_DISRES,
116     F_DISRESVIOL,
117     F_ORIRES,
118     F_ORIRESDEV,
119     F_ANGRES,
120     F_ANGRESZ,
121     F_DIHRES,
122     F_DIHRESVIOL,
123     F_CONSTR,
124     F_CONSTRNC,
125     F_SETTLE,
126     F_VSITE2,
127     F_VSITE3,
128     F_VSITE3FD,
129     F_VSITE3FAD,
130     F_VSITE3OUT,
131     F_VSITE4FD,
132     F_VSITE4FDN,
133     F_VSITEN,
134     F_COM_PULL,
135     F_EQM,
136     F_EPOT,
137     F_EKIN,
138     F_ETOT,
139     F_ECONSERVED,
140     F_TEMP,
141     F_VTEMP_NOLONGERUSED,
142     F_PDISPCORR,
143     F_PRES,
144     F_DVDL_CONSTR,
145     F_DVDL,
146     F_DKDL,
147     F_DVDL_COUL,
148     F_DVDL_VDW,
149     F_DVDL_BONDED,
150     F_DVDL_RESTRAINT,
151     F_DVDL_TEMPERATURE, /* not calculated for now, but should just be the energy (NVT) or enthalpy (NPT), or 0 (NVE) */
152     F_NRE               /* This number is for the total number of energies      */
153 };
154
155 #define IS_RESTRAINT_TYPE(ifunc) ((((ifunc) == F_POSRES) || ((ifunc) == F_FBPOSRES) || ((ifunc) == F_DISRES) || ((ifunc) == F_RESTRBONDS) || ((ifunc) == F_DISRESVIOL) || ((ifunc) == F_ORIRES) || ((ifunc) == F_ORIRESDEV) || ((ifunc) == F_ANGRES) || ((ifunc) == F_ANGRESZ) || ((ifunc) == F_DIHRES)))
156
157 typedef union t_iparams
158 {
159     /* Some parameters have A and B values for free energy calculations.
160      * The B values are not used for regular simulations of course.
161      * Free Energy for nonbondeds can be computed by changing the atom type.
162      * The harmonic type is used for all harmonic potentials:
163      * bonds, angles and improper dihedrals
164      */
165     struct {
166         real a, b, c;
167     } bham;
168     struct {
169         real rA, krA, rB, krB;
170     } harmonic;
171     struct {
172         real klinA, aA, klinB, aB;
173     } linangle;
174     struct {
175         real lowA, up1A, up2A, kA, lowB, up1B, up2B, kB;
176     } restraint;
177     /* No free energy supported for cubic bonds, FENE, WPOL or cross terms */
178     struct {
179         real b0, kb, kcub;
180     } cubic;
181     struct {
182         real bm, kb;
183     } fene;
184     struct {
185         real r1e, r2e, krr;
186     } cross_bb;
187     struct {
188         real r1e, r2e, r3e, krt;
189     } cross_ba;
190     struct {
191         real thetaA, kthetaA, r13A, kUBA, thetaB, kthetaB, r13B, kUBB;
192     } u_b;
193     struct {
194         real theta, c[5];
195     } qangle;
196     struct {
197         real alpha;
198     } polarize;
199     struct {
200         real alpha, drcut, khyp;
201     } anharm_polarize;
202     struct {
203         real al_x, al_y, al_z, rOH, rHH, rOD;
204     } wpol;
205     struct {
206         real a, alpha1, alpha2, rfac;
207     } thole;
208     struct {
209         real c6, c12;
210     } lj;
211     struct {
212         real c6A, c12A, c6B, c12B;
213     } lj14;
214     struct {
215         real fqq, qi, qj, c6, c12;
216     } ljc14;
217     struct {
218         real qi, qj, c6, c12;
219     } ljcnb;
220     /* Proper dihedrals can not have different multiplicity when
221      * doing free energy calculations, because the potential would not
222      * be periodic anymore.
223      */
224     struct {
225         real phiA, cpA; int mult; real phiB, cpB;
226     } pdihs;
227     struct {
228         real dA, dB;
229     } constr;
230     /* Settle can not be used for Free energy calculations of water bond geometry.
231      * Use shake (or lincs) instead if you have to change the water bonds.
232      */
233     struct {
234         real doh, dhh;
235     } settle;
236     struct {
237         real b0A, cbA, betaA, b0B, cbB, betaB;
238     } morse;
239     struct {
240         real pos0A[DIM], fcA[DIM], pos0B[DIM], fcB[DIM];
241     } posres;
242     struct {
243         real pos0[DIM], r, k; int geom;
244     } fbposres;
245     struct {
246         real rbcA[NR_RBDIHS], rbcB[NR_RBDIHS];
247     } rbdihs;
248     struct {
249         real cbtcA[NR_CBTDIHS], cbtcB[NR_CBTDIHS];
250     } cbtdihs;
251     struct {
252         real a, b, c, d, e, f;
253     } vsite;
254     struct {
255         int  n; real a;
256     } vsiten;
257     /* NOTE: npair is only set after reading the tpx file */
258     struct {
259         real low, up1, up2, kfac; int type, label, npair;
260     } disres;
261     struct {
262         real phiA, dphiA, kfacA, phiB, dphiB, kfacB;
263     } dihres;
264     struct {
265         int  ex, power, label; real c, obs, kfac;
266     } orires;
267     struct {
268         int  table; real kA; real kB;
269     } tab;
270     struct {
271         int cmapA, cmapB;
272     } cmap;
273     struct {
274         real buf[MAXFORCEPARAM];
275     } generic;                                               /* Conversion */
276 } t_iparams;
277
278 typedef int t_functype;
279
280 /* List of listed interactions, see description further down.
281  *
282  * TODO: Consider storing the function type as well.
283  * TODO: Consider providing per interaction access.
284  */
285 struct InteractionList
286 {
287     /* Returns the total number of elements in iatoms */
288     int size() const
289     {
290         return iatoms.size();
291     }
292
293     /* List of interactions, see explanation further down */
294     std::vector<int> iatoms;
295 };
296
297 /* List of interaction lists, one list for each interaction type
298  *
299  * TODO: Consider only including entries in use instead of all F_NRE
300  */
301 typedef std::array<InteractionList, F_NRE> InteractionLists;
302
303 /* Deprecated list of listed interactions.
304  *
305  * The nonperturbed/perturbed interactions are now separated (sorted) in the
306  * ilist, such that the first 0..(nr_nonperturbed-1) ones are exactly that, and
307  * the remaining ones from nr_nonperturbed..(nr-1) are perturbed bonded
308  * interactions.
309  */
310 struct t_ilist
311 {
312     /* Returns the total number of elements in iatoms */
313     int size() const
314     {
315         return nr;
316     }
317
318     int      nr;
319     int      nr_nonperturbed;
320     t_iatom *iatoms;
321     int      nalloc;
322 };
323
324 /* TODO: Replace t_ilist in gmx_localtop_t by InteractionList.
325  *       The nr_nonperturbed functionality needs to be ported.
326  *       Remove t_topology.
327  *       Remove t_ilist and remove templating on list type
328  *       in mshift.cpp, constr.cpp, vsite.cpp and domdec_topology.cpp.
329  */
330
331 /*
332  * The structs InteractionList and t_ilist defines a list of atoms with their interactions.
333  * General field description:
334  *   int nr
335  *      the size (nr elements) of the interactions array (iatoms[]).
336  *   t_iatom *iatoms
337  *  specifies which atoms are involved in an interaction of a certain
338  *       type. The layout of this array is as follows:
339  *
340  *        +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
341  *        |type1|at1|at2|at3|type2|at1|at2|type1|at1|at2|at3|type3|at1|at2|
342  *        +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
343  *
344  *  So for interaction type type1 3 atoms are needed, and for type2 and
345  *      type3 only 2. The type identifier is used to select the function to
346  *      calculate the interaction and its actual parameters. This type
347  *      identifier is an index in a params[] and functype[] array.
348  */
349
350 typedef struct
351 {
352     real *cmap; /* Has length 4*grid_spacing*grid_spacing, */
353     /* there are 4 entries for each cmap type (V,dVdx,dVdy,d2dVdxdy) */
354 } gmx_cmapdata_t;
355
356 typedef struct gmx_cmap_t
357 {
358     int             ngrid;        /* Number of allocated cmap (cmapdata_t ) grids */
359     int             grid_spacing; /* Grid spacing */
360     gmx_cmapdata_t *cmapdata;     /* Pointer to grid with actual, pre-interpolated data */
361 } gmx_cmap_t;
362
363
364 typedef struct gmx_ffparams_t
365 {
366     int         ntypes;
367     int         atnr;
368     t_functype *functype;
369     t_iparams  *iparams;
370     double      reppow;    /* The repulsion power for VdW: C12*r^-reppow   */
371     real        fudgeQQ;   /* The scaling factor for Coulomb 1-4: f*q1*q2  */
372     gmx_cmap_t  cmap_grid; /* The dihedral correction maps                 */
373 } gmx_ffparams_t;
374
375 enum {
376     ilsortUNKNOWN, ilsortNO_FE, ilsortFE_UNSORTED, ilsortFE_SORTED
377 };
378
379 typedef struct t_idef
380 {
381     int         ntypes;
382     int         atnr;
383     t_functype *functype;
384     t_iparams  *iparams;
385     real        fudgeQQ;
386     gmx_cmap_t  cmap_grid;
387     t_iparams  *iparams_posres, *iparams_fbposres;
388     int         iparams_posres_nalloc, iparams_fbposres_nalloc;
389
390     t_ilist     il[F_NRE];
391     int         ilsort;
392 } t_idef;
393
394 /*
395  * The struct t_idef defines all the interactions for the complete
396  * simulation. The structure is setup in such a way that the multinode
397  * version of the program  can use it as easy as the single node version.
398  * General field description:
399  *   int ntypes
400  *      defines the number of elements in functype[] and param[].
401  *   int nodeid
402  *      the node id (if parallel machines)
403  *   int atnr
404  *      the number of atomtypes
405  *   t_functype *functype
406  *      array of length ntypes, defines for every force type what type of
407  *      function to use. Every "bond" with the same function but different
408  *      force parameters is a different force type. The type identifier in the
409  *      forceatoms[] array is an index in this array.
410  *   t_iparams *iparams
411  *      array of length ntypes, defines the parameters for every interaction
412  *      type. The type identifier in the actual interaction list
413  *      (ilist[ftype].iatoms[]) is an index in this array.
414  *   gmx_cmap_t cmap_grid
415  *      the grid for the dihedral pair correction maps.
416  *   t_iparams *iparams_posres, *iparams_fbposres
417  *      defines the parameters for position restraints only.
418  *      Position restraints are the only interactions that have different
419  *      parameters (reference positions) for different molecules
420  *      of the same type. ilist[F_POSRES].iatoms[] is an index in this array.
421  *   t_ilist il[F_NRE]
422  *      The list of interactions for each type. Note that some,
423  *      such as LJ and COUL will have 0 entries.
424  *   int ilsort
425  *      The state of the sorting of il, values are provided above.
426  */
427
428 void pr_iparams(FILE *fp, t_functype ftype, const t_iparams *iparams);
429 void pr_ilist(FILE *fp, int indent, const char *title,
430               const t_functype *functype, const InteractionList *ilist,
431               gmx_bool bShowNumbers,
432               gmx_bool bShowParameters, const t_iparams *iparams);
433 void pr_ffparams(FILE *fp, int indent, const char *title,
434                  const gmx_ffparams_t *ffparams, gmx_bool bShowNumbers);
435 void pr_idef(FILE *fp, int indent, const char *title, const t_idef *idef,
436              gmx_bool bShowNumbers, gmx_bool bShowParameters);
437
438 #endif