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