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