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