bafd389beef3e8f175c292fffb3efab499b72c20
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / convparm.cpp
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,2017 by the GROMACS development team.
7  * Copyright (c) 2018,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 /* This file is completely threadsafe - keep it that way! */
39 #include "gmxpre.h"
40
41 #include "convparm.h"
42
43 #include <cassert>
44 #include <cmath>
45 #include <cstring>
46
47 #include <memory>
48
49 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
50 #include "gromacs/gmxpreprocess/grompp_impl.h"
51 #include "gromacs/gmxpreprocess/topio.h"
52 #include "gromacs/gmxpreprocess/toputil.h"
53 #include "gromacs/math/functions.h"
54 #include "gromacs/math/units.h"
55 #include "gromacs/math/utilities.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/mdtypes/md_enums.h"
58 #include "gromacs/topology/ifunc.h"
59 #include "gromacs/topology/topology.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/smalloc.h"
62
63 static int round_check(real r, int limit, int ftype, const char* name)
64 {
65     const int i = gmx::roundToInt(r);
66
67     if (r - i > 0.01 || r - i < -0.01)
68     {
69         gmx_fatal(FARGS,
70                   "A non-integer value (%f) was supplied for '%s' in %s",
71                   r,
72                   name,
73                   interaction_function[ftype].longname);
74     }
75
76     if (i < limit)
77     {
78         gmx_fatal(FARGS,
79                   "Value of '%s' in %s is %d, which is smaller than the minimum of %d",
80                   name,
81                   interaction_function[ftype].longname,
82                   i,
83                   limit);
84     }
85
86     return i;
87 }
88
89 static void set_ljparams(int comb, double reppow, double v, double w, real* c6, real* c12)
90 {
91     if (comb == eCOMB_ARITHMETIC || comb == eCOMB_GEOM_SIG_EPS)
92     {
93         if (v >= 0)
94         {
95             *c6  = 4 * w * gmx::power6(v);
96             *c12 = 4 * w * std::pow(v, reppow);
97         }
98         else
99         {
100             /* Interpret negative sigma as c6=0 and c12 with -sigma */
101             *c6  = 0;
102             *c12 = 4 * w * std::pow(-v, reppow);
103         }
104     }
105     else
106     {
107         *c6  = v;
108         *c12 = w;
109     }
110 }
111
112 /* A return value of 0 means parameters were assigned successfully,
113  * returning -1 means this is an all-zero interaction that should not be added.
114  */
115 static int assign_param(t_functype ftype, t_iparams* newparam, gmx::ArrayRef<const real> old, int comb, double reppow)
116 {
117     bool all_param_zero = true;
118
119     /* Set to zero */
120     for (int j = 0; (j < MAXFORCEPARAM); j++)
121     {
122         newparam->generic.buf[j] = 0.0;
123         /* If all parameters are zero we might not add some interaction types (selected below).
124          * We cannot apply this to ALL interactions, since many have valid reasons for having
125          * zero parameters (e.g. an index to a Cmap interaction, or LJ parameters), but
126          * we use it for angles and torsions that are typically generated automatically.
127          */
128         all_param_zero = all_param_zero && fabs(old[j]) < GMX_REAL_MIN;
129     }
130
131     if (all_param_zero)
132     {
133         if (IS_ANGLE(ftype) || IS_RESTRAINT_TYPE(ftype) || ftype == F_IDIHS || ftype == F_PDIHS
134             || ftype == F_PIDIHS || ftype == F_RBDIHS || ftype == F_FOURDIHS)
135         {
136             return -1;
137         }
138     }
139
140     switch (ftype)
141     {
142         case F_G96ANGLES:
143             /* Post processing of input data: store cosine iso angle itself */
144             newparam->harmonic.rA  = cos(old[0] * DEG2RAD);
145             newparam->harmonic.krA = old[1];
146             newparam->harmonic.rB  = cos(old[2] * DEG2RAD);
147             newparam->harmonic.krB = old[3];
148             break;
149         case F_G96BONDS:
150             /* Post processing of input data: store square of length itself */
151             newparam->harmonic.rA  = gmx::square(old[0]);
152             newparam->harmonic.krA = old[1];
153             newparam->harmonic.rB  = gmx::square(old[2]);
154             newparam->harmonic.krB = old[3];
155             break;
156         case F_FENEBONDS:
157             newparam->fene.bm = old[0];
158             newparam->fene.kb = old[1];
159             break;
160         case F_RESTRBONDS:
161             newparam->restraint.lowA = old[0];
162             newparam->restraint.up1A = old[1];
163             newparam->restraint.up2A = old[2];
164             newparam->restraint.kA   = old[3];
165             newparam->restraint.lowB = old[4];
166             newparam->restraint.up1B = old[5];
167             newparam->restraint.up2B = old[6];
168             newparam->restraint.kB   = old[7];
169             break;
170         case F_TABBONDS:
171         case F_TABBONDSNC:
172         case F_TABANGLES:
173         case F_TABDIHS:
174             newparam->tab.table = round_check(old[0], 0, ftype, "table index");
175             newparam->tab.kA    = old[1];
176             newparam->tab.kB    = old[3];
177             break;
178         case F_CROSS_BOND_BONDS:
179             newparam->cross_bb.r1e = old[0];
180             newparam->cross_bb.r2e = old[1];
181             newparam->cross_bb.krr = old[2];
182             break;
183         case F_CROSS_BOND_ANGLES:
184             newparam->cross_ba.r1e = old[0];
185             newparam->cross_ba.r2e = old[1];
186             newparam->cross_ba.r3e = old[2];
187             newparam->cross_ba.krt = old[3];
188             break;
189         case F_UREY_BRADLEY:
190             newparam->u_b.thetaA  = old[0];
191             newparam->u_b.kthetaA = old[1];
192             newparam->u_b.r13A    = old[2];
193             newparam->u_b.kUBA    = old[3];
194             newparam->u_b.thetaB  = old[4];
195             newparam->u_b.kthetaB = old[5];
196             newparam->u_b.r13B    = old[6];
197             newparam->u_b.kUBB    = old[7];
198             break;
199         case F_QUARTIC_ANGLES:
200             newparam->qangle.theta = old[0];
201             for (int i = 0; i < 5; i++)
202             {
203                 newparam->qangle.c[i] = old[i + 1];
204             }
205             break;
206         case F_LINEAR_ANGLES:
207             newparam->linangle.aA    = old[0];
208             newparam->linangle.klinA = old[1];
209             newparam->linangle.aB    = old[2];
210             newparam->linangle.klinB = old[3];
211             break;
212         case F_BONDS:
213         case F_ANGLES:
214         case F_HARMONIC:
215         case F_IDIHS:
216             newparam->harmonic.rA  = old[0];
217             newparam->harmonic.krA = old[1];
218             newparam->harmonic.rB  = old[2];
219             newparam->harmonic.krB = old[3];
220             break;
221         case F_RESTRANGLES:
222             newparam->harmonic.rA  = old[0];
223             newparam->harmonic.krA = old[1];
224             break;
225         case F_MORSE:
226             newparam->morse.b0A   = old[0];
227             newparam->morse.cbA   = old[1];
228             newparam->morse.betaA = old[2];
229             newparam->morse.b0B   = old[3];
230             newparam->morse.cbB   = old[4];
231             newparam->morse.betaB = old[5];
232             break;
233         case F_CUBICBONDS:
234             newparam->cubic.b0   = old[0];
235             newparam->cubic.kb   = old[1];
236             newparam->cubic.kcub = old[2];
237             break;
238         case F_CONNBONDS: break;
239         case F_POLARIZATION: newparam->polarize.alpha = old[0]; break;
240         case F_ANHARM_POL:
241             newparam->anharm_polarize.alpha = old[0];
242             newparam->anharm_polarize.drcut = old[1];
243             newparam->anharm_polarize.khyp  = old[2];
244             break;
245         case F_WATER_POL:
246             newparam->wpol.al_x = old[0];
247             newparam->wpol.al_y = old[1];
248             newparam->wpol.al_z = old[2];
249             newparam->wpol.rOH  = old[3];
250             newparam->wpol.rHH  = old[4];
251             newparam->wpol.rOD  = old[5];
252             break;
253         case F_THOLE_POL:
254             newparam->thole.a      = old[0];
255             newparam->thole.alpha1 = old[1];
256             newparam->thole.alpha2 = old[2];
257             if ((old[1] > 0) && (old[2] > 0))
258             {
259                 newparam->thole.rfac = old[0] * gmx::invsixthroot(old[1] * old[2]);
260             }
261             else
262             {
263                 newparam->thole.rfac = 1;
264             }
265             break;
266         case F_BHAM:
267             newparam->bham.a = old[0];
268             newparam->bham.b = old[1];
269             newparam->bham.c = old[2];
270             break;
271         case F_LJ14:
272             set_ljparams(comb, reppow, old[0], old[1], &newparam->lj14.c6A, &newparam->lj14.c12A);
273             set_ljparams(comb, reppow, old[2], old[3], &newparam->lj14.c6B, &newparam->lj14.c12B);
274             break;
275         case F_LJC14_Q:
276             newparam->ljc14.fqq = old[0];
277             newparam->ljc14.qi  = old[1];
278             newparam->ljc14.qj  = old[2];
279             set_ljparams(comb, reppow, old[3], old[4], &newparam->ljc14.c6, &newparam->ljc14.c12);
280             break;
281         case F_LJC_PAIRS_NB:
282             newparam->ljcnb.qi = old[0];
283             newparam->ljcnb.qj = old[1];
284             set_ljparams(comb, reppow, old[2], old[3], &newparam->ljcnb.c6, &newparam->ljcnb.c12);
285             break;
286         case F_LJ:
287             set_ljparams(comb, reppow, old[0], old[1], &newparam->lj.c6, &newparam->lj.c12);
288             break;
289         case F_PDIHS:
290         case F_PIDIHS:
291         case F_ANGRES:
292         case F_ANGRESZ:
293             newparam->pdihs.phiA = old[0];
294             newparam->pdihs.cpA  = old[1];
295
296             /* Change 20100720: Amber occasionally uses negative multiplicities (mathematically OK),
297              * so I have changed the lower limit to -99 /EL
298              */
299             newparam->pdihs.phiB = old[3];
300             newparam->pdihs.cpB  = old[4];
301             /* If both force constants are zero there is no interaction. Return -1 to signal
302              * this entry should NOT be added.
303              */
304             if (fabs(newparam->pdihs.cpA) < GMX_REAL_MIN && fabs(newparam->pdihs.cpB) < GMX_REAL_MIN)
305             {
306                 return -1;
307             }
308
309             newparam->pdihs.mult = round_check(old[2], -99, ftype, "multiplicity");
310
311             break;
312         case F_RESTRDIHS:
313             newparam->pdihs.phiA = old[0];
314             newparam->pdihs.cpA  = old[1];
315             break;
316         case F_POSRES:
317             newparam->posres.fcA[XX]   = old[0];
318             newparam->posres.fcA[YY]   = old[1];
319             newparam->posres.fcA[ZZ]   = old[2];
320             newparam->posres.fcB[XX]   = old[3];
321             newparam->posres.fcB[YY]   = old[4];
322             newparam->posres.fcB[ZZ]   = old[5];
323             newparam->posres.pos0A[XX] = old[6];
324             newparam->posres.pos0A[YY] = old[7];
325             newparam->posres.pos0A[ZZ] = old[8];
326             newparam->posres.pos0B[XX] = old[9];
327             newparam->posres.pos0B[YY] = old[10];
328             newparam->posres.pos0B[ZZ] = old[11];
329             break;
330         case F_FBPOSRES:
331             newparam->fbposres.geom = round_check(old[0], 0, ftype, "geometry");
332             if (!(newparam->fbposres.geom > efbposresZERO && newparam->fbposres.geom < efbposresNR))
333             {
334                 gmx_fatal(FARGS,
335                           "Invalid geometry for flat-bottomed position restraint.\n"
336                           "Expected number between 1 and %d. Found %d\n",
337                           efbposresNR - 1,
338                           newparam->fbposres.geom);
339             }
340             newparam->fbposres.r        = old[1];
341             newparam->fbposres.k        = old[2];
342             newparam->fbposres.pos0[XX] = old[3];
343             newparam->fbposres.pos0[YY] = old[4];
344             newparam->fbposres.pos0[ZZ] = old[5];
345             break;
346         case F_DISRES:
347             newparam->disres.label = round_check(old[0], 0, ftype, "label");
348             newparam->disres.type  = round_check(old[1], 1, ftype, "type'");
349             newparam->disres.low   = old[2];
350             newparam->disres.up1   = old[3];
351             newparam->disres.up2   = old[4];
352             newparam->disres.kfac  = old[5];
353             break;
354         case F_ORIRES:
355             newparam->orires.ex    = round_check(old[0], 1, ftype, "experiment") - 1;
356             newparam->orires.label = round_check(old[1], 1, ftype, "label");
357             newparam->orires.power = round_check(old[2], 0, ftype, "power");
358             newparam->orires.c     = old[3];
359             newparam->orires.obs   = old[4];
360             newparam->orires.kfac  = old[5];
361             break;
362         case F_DIHRES:
363             newparam->dihres.phiA  = old[0];
364             newparam->dihres.dphiA = old[1];
365             newparam->dihres.kfacA = old[2];
366             newparam->dihres.phiB  = old[3];
367             newparam->dihres.dphiB = old[4];
368             newparam->dihres.kfacB = old[5];
369             break;
370         case F_RBDIHS:
371             for (int i = 0; (i < NR_RBDIHS); i++)
372             {
373                 newparam->rbdihs.rbcA[i] = old[i];
374                 newparam->rbdihs.rbcB[i] = old[NR_RBDIHS + i];
375             }
376             break;
377         case F_CBTDIHS:
378             for (int i = 0; (i < NR_CBTDIHS); i++)
379             {
380                 newparam->cbtdihs.cbtcA[i] = old[i];
381             }
382             break;
383         case F_FOURDIHS:
384             /* Read the dihedral parameters to temporary arrays,
385              * and convert them to the computationally faster
386              * Ryckaert-Bellemans form.
387              */
388             /* Use conversion formula for OPLS to Ryckaert-Bellemans: */
389             newparam->rbdihs.rbcA[0] = old[1] + 0.5 * (old[0] + old[2]);
390             newparam->rbdihs.rbcA[1] = 0.5 * (3.0 * old[2] - old[0]);
391             newparam->rbdihs.rbcA[2] = 4.0 * old[3] - old[1];
392             newparam->rbdihs.rbcA[3] = -2.0 * old[2];
393             newparam->rbdihs.rbcA[4] = -4.0 * old[3];
394             newparam->rbdihs.rbcA[5] = 0.0;
395
396             newparam->rbdihs.rbcB[0] =
397                     old[NR_FOURDIHS + 1] + 0.5 * (old[NR_FOURDIHS + 0] + old[NR_FOURDIHS + 2]);
398             newparam->rbdihs.rbcB[1] = 0.5 * (3.0 * old[NR_FOURDIHS + 2] - old[NR_FOURDIHS + 0]);
399             newparam->rbdihs.rbcB[2] = 4.0 * old[NR_FOURDIHS + 3] - old[NR_FOURDIHS + 1];
400             newparam->rbdihs.rbcB[3] = -2.0 * old[NR_FOURDIHS + 2];
401             newparam->rbdihs.rbcB[4] = -4.0 * old[NR_FOURDIHS + 3];
402             newparam->rbdihs.rbcB[5] = 0.0;
403             break;
404         case F_CONSTR:
405         case F_CONSTRNC:
406             newparam->constr.dA = old[0];
407             newparam->constr.dB = old[1];
408             break;
409         case F_SETTLE:
410             newparam->settle.doh = old[0];
411             newparam->settle.dhh = old[1];
412             break;
413         case F_VSITE1:
414         case F_VSITE2:
415         case F_VSITE2FD:
416         case F_VSITE3:
417         case F_VSITE3FD:
418         case F_VSITE3OUT:
419         case F_VSITE4FD:
420         case F_VSITE4FDN:
421             newparam->vsite.a = old[0];
422             newparam->vsite.b = old[1];
423             newparam->vsite.c = old[2];
424             newparam->vsite.d = old[3];
425             newparam->vsite.e = old[4];
426             newparam->vsite.f = old[5];
427             break;
428         case F_VSITE3FAD:
429             newparam->vsite.a = old[1] * cos(DEG2RAD * old[0]);
430             newparam->vsite.b = old[1] * sin(DEG2RAD * old[0]);
431             newparam->vsite.c = old[2];
432             newparam->vsite.d = old[3];
433             newparam->vsite.e = old[4];
434             newparam->vsite.f = old[5];
435             break;
436         case F_VSITEN:
437             newparam->vsiten.n = round_check(old[0], 1, ftype, "number of atoms");
438             newparam->vsiten.a = old[1];
439             break;
440         case F_CMAP:
441             newparam->cmap.cmapA = static_cast<int>(old[0]);
442             newparam->cmap.cmapB = static_cast<int>(old[1]);
443             break;
444         case F_GB12_NOLONGERUSED:
445         case F_GB13_NOLONGERUSED:
446         case F_GB14_NOLONGERUSED: break;
447         default:
448             gmx_fatal(FARGS, "unknown function type %d in %s line %d", ftype, __FILE__, __LINE__);
449     }
450     return 0;
451 }
452
453 static int enter_params(gmx_ffparams_t*           ffparams,
454                         t_functype                ftype,
455                         gmx::ArrayRef<const real> forceparams,
456                         int                       comb,
457                         real                      reppow,
458                         int                       start,
459                         bool                      bAppend)
460 {
461     t_iparams newparam;
462     int       rc;
463
464     if ((rc = assign_param(ftype, &newparam, forceparams, comb, reppow)) < 0)
465     {
466         /* -1 means this interaction is all-zero and should not be added */
467         return rc;
468     }
469
470     if (!bAppend)
471     {
472         if (ftype != F_DISRES)
473         {
474             for (int type = start; type < ffparams->numTypes(); type++)
475             {
476                 // Note that the first condition is always met by starting the loop at start
477                 if (ffparams->functype[type] == ftype
478                     && memcmp(&newparam, &ffparams->iparams[type], static_cast<size_t>(sizeof(newparam))) == 0)
479                 {
480                     return type;
481                 }
482             }
483         }
484         else
485         {
486             // Distance restraints should have unique labels and pairs with the same label
487             // should be consecutive, so we here we only need to check the last type in the list.
488             // This changes the complexity from quadratic to linear in the number of restraints.
489             const int type = ffparams->numTypes() - 1;
490             if (type >= 0 && ffparams->functype[type] == ftype
491                 && memcmp(&newparam, &ffparams->iparams[type], static_cast<size_t>(sizeof(newparam))) == 0)
492             {
493                 return type;
494             }
495         }
496     }
497
498     const int type = ffparams->numTypes();
499
500     ffparams->iparams.push_back(newparam);
501     ffparams->functype.push_back(ftype);
502
503     GMX_ASSERT(ffparams->iparams.size() == ffparams->functype.size(), "sizes should match");
504
505     return type;
506 }
507
508 static void append_interaction(InteractionList* ilist, int type, gmx::ArrayRef<const int> a)
509 {
510     ilist->iatoms.push_back(type);
511     for (const auto& atom : a)
512     {
513         ilist->iatoms.push_back(atom);
514     }
515 }
516
517 static void enter_function(const InteractionsOfType* p,
518                            t_functype                ftype,
519                            int                       comb,
520                            real                      reppow,
521                            gmx_ffparams_t*           ffparams,
522                            InteractionList*          il,
523                            bool                      bNB,
524                            bool                      bAppend)
525 {
526     int start = ffparams->numTypes();
527
528     for (auto& parm : p->interactionTypes)
529     {
530         int type = enter_params(ffparams, ftype, parm.forceParam(), comb, reppow, start, bAppend);
531         /* Type==-1 is used as a signal that this interaction is all-zero and should not be added. */
532         if (!bNB && type >= 0)
533         {
534             GMX_RELEASE_ASSERT(il, "Need valid interaction list");
535             GMX_RELEASE_ASSERT(parm.atoms().ssize() == NRAL(ftype),
536                                "Need to have correct number of atoms for the parameter");
537             append_interaction(il, type, parm.atoms());
538         }
539     }
540 }
541
542 void convertInteractionsOfType(int                                      atnr,
543                                gmx::ArrayRef<const InteractionsOfType>  nbtypes,
544                                gmx::ArrayRef<const MoleculeInformation> mi,
545                                const MoleculeInformation*               intermolecular_interactions,
546                                int                                      comb,
547                                double                                   reppow,
548                                real                                     fudgeQQ,
549                                gmx_mtop_t*                              mtop)
550 {
551     int             i;
552     unsigned long   flags;
553     gmx_ffparams_t* ffp;
554     gmx_moltype_t*  molt;
555
556     ffp       = &mtop->ffparams;
557     ffp->atnr = atnr;
558     ffp->functype.clear();
559     ffp->iparams.clear();
560     ffp->reppow = reppow;
561
562     enter_function(&(nbtypes[F_LJ]), static_cast<t_functype>(F_LJ), comb, reppow, ffp, nullptr, TRUE, TRUE);
563     enter_function(
564             &(nbtypes[F_BHAM]), static_cast<t_functype>(F_BHAM), comb, reppow, ffp, nullptr, TRUE, TRUE);
565
566     for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
567     {
568         molt = &mtop->moltype[mt];
569         for (i = 0; (i < F_NRE); i++)
570         {
571             molt->ilist[i].iatoms.clear();
572
573             gmx::ArrayRef<const InteractionsOfType> interactions = mi[mt].interactions;
574
575             flags = interaction_function[i].flags;
576             if ((i != F_LJ) && (i != F_BHAM)
577                 && ((flags & IF_BOND) || (flags & IF_VSITE) || (flags & IF_CONSTRAINT)))
578             {
579                 enter_function(&(interactions[i]),
580                                static_cast<t_functype>(i),
581                                comb,
582                                reppow,
583                                ffp,
584                                &molt->ilist[i],
585                                FALSE,
586                                (i == F_POSRES || i == F_FBPOSRES));
587             }
588         }
589     }
590
591     mtop->bIntermolecularInteractions = FALSE;
592     if (intermolecular_interactions != nullptr)
593     {
594         /* Process the intermolecular interaction list */
595         mtop->intermolecular_ilist = std::make_unique<InteractionLists>();
596
597         for (i = 0; (i < F_NRE); i++)
598         {
599             (*mtop->intermolecular_ilist)[i].iatoms.clear();
600
601             gmx::ArrayRef<const InteractionsOfType> interactions = intermolecular_interactions->interactions;
602
603             if (!interactions[i].interactionTypes.empty())
604             {
605                 flags = interaction_function[i].flags;
606                 /* For intermolecular interactions we (currently)
607                  * only support potentials.
608                  * Constraints and virtual sites would be possible,
609                  * but require a lot of extra (bug-prone) code.
610                  */
611                 if (!(flags & IF_BOND))
612                 {
613                     gmx_fatal(FARGS,
614                               "The intermolecular_interaction section may only contain bonded "
615                               "potentials");
616                 }
617                 else if (NRAL(i) == 1) /* e.g. position restraints */
618                 {
619                     gmx_fatal(FARGS,
620                               "Single atom interactions don't make sense in the "
621                               "intermolecular_interaction section, you can put them in the "
622                               "moleculetype section");
623                 }
624                 else if (flags & IF_CHEMBOND)
625                 {
626                     gmx_fatal(FARGS,
627                               "The intermolecular_interaction can not contain chemically bonding "
628                               "interactions");
629                 }
630                 else
631                 {
632                     enter_function(&(interactions[i]),
633                                    static_cast<t_functype>(i),
634                                    comb,
635                                    reppow,
636                                    ffp,
637                                    &(*mtop->intermolecular_ilist)[i],
638                                    FALSE,
639                                    FALSE);
640
641                     mtop->bIntermolecularInteractions = TRUE;
642                 }
643             }
644         }
645
646         if (!mtop->bIntermolecularInteractions)
647         {
648             mtop->intermolecular_ilist.reset(nullptr);
649         }
650     }
651
652     ffp->fudgeQQ = fudgeQQ;
653 }