2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 /* This file is completely threadsafe - keep it that way! */
46 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
47 #include "gromacs/gmxpreprocess/topio.h"
48 #include "gromacs/gmxpreprocess/toputil.h"
49 #include "gromacs/math/functions.h"
50 #include "gromacs/math/units.h"
51 #include "gromacs/math/utilities.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/mdtypes/md_enums.h"
54 #include "gromacs/topology/ifunc.h"
55 #include "gromacs/topology/topology.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/smalloc.h"
59 static int round_check(real r, int limit, int ftype, const char *name)
61 const int i = gmx::roundToInt(r);
63 if (r-i > 0.01 || r-i < -0.01)
65 gmx_fatal(FARGS, "A non-integer value (%f) was supplied for '%s' in %s",
66 r, name, interaction_function[ftype].longname);
71 gmx_fatal(FARGS, "Value of '%s' in %s is %d, which is smaller than the minimum of %d",
72 name, interaction_function[ftype].longname, i, limit);
78 static void set_ljparams(int comb, double reppow, double v, double w,
81 if (comb == eCOMB_ARITHMETIC || comb == eCOMB_GEOM_SIG_EPS)
85 *c6 = 4*w*gmx::power6(v);
86 *c12 = 4*w*std::pow(v, reppow);
90 /* Interpret negative sigma as c6=0 and c12 with -sigma */
92 *c12 = 4*w*std::pow(-v, reppow);
102 /* A return value of 0 means parameters were assigned successfully,
103 * returning -1 means this is an all-zero interaction that should not be added.
106 assign_param(t_functype ftype, t_iparams *newparam,
107 real old[MAXFORCEPARAM], int comb, double reppow)
110 bool all_param_zero = TRUE;
113 for (j = 0; (j < MAXFORCEPARAM); j++)
115 newparam->generic.buf[j] = 0.0;
116 /* If all parameters are zero we might not add some interaction types (selected below).
117 * We cannot apply this to ALL interactions, since many have valid reasons for having
118 * zero parameters (e.g. an index to a Cmap interaction, or LJ parameters), but
119 * we use it for angles and torsions that are typically generated automatically.
121 all_param_zero = all_param_zero && fabs(old[j]) < GMX_REAL_MIN;
126 if (IS_ANGLE(ftype) || IS_RESTRAINT_TYPE(ftype) || ftype == F_IDIHS ||
127 ftype == F_PDIHS || ftype == F_PIDIHS || ftype == F_RBDIHS || ftype == F_FOURDIHS)
136 /* Post processing of input data: store cosine iso angle itself */
137 newparam->harmonic.rA = cos(old[0]*DEG2RAD);
138 newparam->harmonic.krA = old[1];
139 newparam->harmonic.rB = cos(old[2]*DEG2RAD);
140 newparam->harmonic.krB = old[3];
143 /* Post processing of input data: store square of length itself */
144 newparam->harmonic.rA = gmx::square(old[0]);
145 newparam->harmonic.krA = old[1];
146 newparam->harmonic.rB = gmx::square(old[2]);
147 newparam->harmonic.krB = old[3];
150 newparam->fene.bm = old[0];
151 newparam->fene.kb = old[1];
154 newparam->restraint.lowA = old[0];
155 newparam->restraint.up1A = old[1];
156 newparam->restraint.up2A = old[2];
157 newparam->restraint.kA = old[3];
158 newparam->restraint.lowB = old[4];
159 newparam->restraint.up1B = old[5];
160 newparam->restraint.up2B = old[6];
161 newparam->restraint.kB = old[7];
167 newparam->tab.table = round_check(old[0], 0, ftype, "table index");
168 newparam->tab.kA = old[1];
169 newparam->tab.kB = old[3];
171 case F_CROSS_BOND_BONDS:
172 newparam->cross_bb.r1e = old[0];
173 newparam->cross_bb.r2e = old[1];
174 newparam->cross_bb.krr = old[2];
176 case F_CROSS_BOND_ANGLES:
177 newparam->cross_ba.r1e = old[0];
178 newparam->cross_ba.r2e = old[1];
179 newparam->cross_ba.r3e = old[2];
180 newparam->cross_ba.krt = old[3];
183 newparam->u_b.thetaA = old[0];
184 newparam->u_b.kthetaA = old[1];
185 newparam->u_b.r13A = old[2];
186 newparam->u_b.kUBA = old[3];
187 newparam->u_b.thetaB = old[4];
188 newparam->u_b.kthetaB = old[5];
189 newparam->u_b.r13B = old[6];
190 newparam->u_b.kUBB = old[7];
192 case F_QUARTIC_ANGLES:
193 newparam->qangle.theta = old[0];
194 for (i = 0; i < 5; i++)
196 newparam->qangle.c[i] = old[i+1];
199 case F_LINEAR_ANGLES:
200 newparam->linangle.aA = old[0];
201 newparam->linangle.klinA = old[1];
202 newparam->linangle.aB = old[2];
203 newparam->linangle.klinB = old[3];
209 newparam->harmonic.rA = old[0];
210 newparam->harmonic.krA = old[1];
211 newparam->harmonic.rB = old[2];
212 newparam->harmonic.krB = old[3];
215 newparam->harmonic.rA = old[0];
216 newparam->harmonic.krA = old[1];
219 newparam->morse.b0A = old[0];
220 newparam->morse.cbA = old[1];
221 newparam->morse.betaA = old[2];
222 newparam->morse.b0B = old[3];
223 newparam->morse.cbB = old[4];
224 newparam->morse.betaB = old[5];
227 newparam->cubic.b0 = old[0];
228 newparam->cubic.kb = old[1];
229 newparam->cubic.kcub = old[2];
234 newparam->polarize.alpha = old[0];
237 newparam->anharm_polarize.alpha = old[0];
238 newparam->anharm_polarize.drcut = old[1];
239 newparam->anharm_polarize.khyp = old[2];
242 newparam->wpol.al_x = old[0];
243 newparam->wpol.al_y = old[1];
244 newparam->wpol.al_z = old[2];
245 newparam->wpol.rOH = old[3];
246 newparam->wpol.rHH = old[4];
247 newparam->wpol.rOD = old[5];
250 newparam->thole.a = old[0];
251 newparam->thole.alpha1 = old[1];
252 newparam->thole.alpha2 = old[2];
253 if ((old[1] > 0) && (old[2] > 0))
255 newparam->thole.rfac = old[0]*gmx::invsixthroot(old[1]*old[2]);
259 newparam->thole.rfac = 1;
263 newparam->bham.a = old[0];
264 newparam->bham.b = old[1];
265 newparam->bham.c = old[2];
268 set_ljparams(comb, reppow, old[0], old[1], &newparam->lj14.c6A, &newparam->lj14.c12A);
269 set_ljparams(comb, reppow, old[2], old[3], &newparam->lj14.c6B, &newparam->lj14.c12B);
272 newparam->ljc14.fqq = old[0];
273 newparam->ljc14.qi = old[1];
274 newparam->ljc14.qj = old[2];
275 set_ljparams(comb, reppow, old[3], old[4], &newparam->ljc14.c6, &newparam->ljc14.c12);
278 newparam->ljcnb.qi = old[0];
279 newparam->ljcnb.qj = old[1];
280 set_ljparams(comb, reppow, old[2], old[3], &newparam->ljcnb.c6, &newparam->ljcnb.c12);
283 set_ljparams(comb, reppow, old[0], old[1], &newparam->lj.c6, &newparam->lj.c12);
289 newparam->pdihs.phiA = old[0];
290 newparam->pdihs.cpA = old[1];
292 /* Change 20100720: Amber occasionally uses negative multiplicities (mathematically OK),
293 * so I have changed the lower limit to -99 /EL
295 newparam->pdihs.phiB = old[3];
296 newparam->pdihs.cpB = old[4];
297 /* If both force constants are zero there is no interaction. Return -1 to signal
298 * this entry should NOT be added.
300 if (fabs(newparam->pdihs.cpA) < GMX_REAL_MIN && fabs(newparam->pdihs.cpB) < GMX_REAL_MIN)
305 newparam->pdihs.mult = round_check(old[2], -99, ftype, "multiplicity");
309 newparam->pdihs.phiA = old[0];
310 newparam->pdihs.cpA = old[1];
313 newparam->posres.fcA[XX] = old[0];
314 newparam->posres.fcA[YY] = old[1];
315 newparam->posres.fcA[ZZ] = old[2];
316 newparam->posres.fcB[XX] = old[3];
317 newparam->posres.fcB[YY] = old[4];
318 newparam->posres.fcB[ZZ] = old[5];
319 newparam->posres.pos0A[XX] = old[6];
320 newparam->posres.pos0A[YY] = old[7];
321 newparam->posres.pos0A[ZZ] = old[8];
322 newparam->posres.pos0B[XX] = old[9];
323 newparam->posres.pos0B[YY] = old[10];
324 newparam->posres.pos0B[ZZ] = old[11];
327 newparam->fbposres.geom = round_check(old[0], 0, ftype, "geometry");
328 if (!(newparam->fbposres.geom > efbposresZERO && newparam->fbposres.geom < efbposresNR))
330 gmx_fatal(FARGS, "Invalid geometry for flat-bottomed position restraint.\n"
331 "Expected number between 1 and %d. Found %d\n", efbposresNR-1,
332 newparam->fbposres.geom);
334 newparam->fbposres.r = old[1];
335 newparam->fbposres.k = old[2];
336 newparam->fbposres.pos0[XX] = old[3];
337 newparam->fbposres.pos0[YY] = old[4];
338 newparam->fbposres.pos0[ZZ] = old[5];
341 newparam->disres.label = round_check(old[0], 0, ftype, "label");
342 newparam->disres.type = round_check(old[1], 1, ftype, "type'");
343 newparam->disres.low = old[2];
344 newparam->disres.up1 = old[3];
345 newparam->disres.up2 = old[4];
346 newparam->disres.kfac = old[5];
349 newparam->orires.ex = round_check(old[0], 1, ftype, "experiment") - 1;
350 newparam->orires.label = round_check(old[1], 1, ftype, "label");
351 newparam->orires.power = round_check(old[2], 0, ftype, "power");
352 newparam->orires.c = old[3];
353 newparam->orires.obs = old[4];
354 newparam->orires.kfac = old[5];
357 newparam->dihres.phiA = old[0];
358 newparam->dihres.dphiA = old[1];
359 newparam->dihres.kfacA = old[2];
360 newparam->dihres.phiB = old[3];
361 newparam->dihres.dphiB = old[4];
362 newparam->dihres.kfacB = old[5];
365 for (i = 0; (i < NR_RBDIHS); i++)
367 newparam->rbdihs.rbcA[i] = old[i];
368 newparam->rbdihs.rbcB[i] = old[NR_RBDIHS+i];
372 for (i = 0; (i < NR_CBTDIHS); i++)
374 newparam->cbtdihs.cbtcA[i] = old[i];
378 /* Read the dihedral parameters to temporary arrays,
379 * and convert them to the computationally faster
380 * Ryckaert-Bellemans form.
382 /* Use conversion formula for OPLS to Ryckaert-Bellemans: */
383 newparam->rbdihs.rbcA[0] = old[1]+0.5*(old[0]+old[2]);
384 newparam->rbdihs.rbcA[1] = 0.5*(3.0*old[2]-old[0]);
385 newparam->rbdihs.rbcA[2] = 4.0*old[3]-old[1];
386 newparam->rbdihs.rbcA[3] = -2.0*old[2];
387 newparam->rbdihs.rbcA[4] = -4.0*old[3];
388 newparam->rbdihs.rbcA[5] = 0.0;
390 newparam->rbdihs.rbcB[0] = old[NR_FOURDIHS+1]+0.5*(old[NR_FOURDIHS+0]+old[NR_FOURDIHS+2]);
391 newparam->rbdihs.rbcB[1] = 0.5*(3.0*old[NR_FOURDIHS+2]-old[NR_FOURDIHS+0]);
392 newparam->rbdihs.rbcB[2] = 4.0*old[NR_FOURDIHS+3]-old[NR_FOURDIHS+1];
393 newparam->rbdihs.rbcB[3] = -2.0*old[NR_FOURDIHS+2];
394 newparam->rbdihs.rbcB[4] = -4.0*old[NR_FOURDIHS+3];
395 newparam->rbdihs.rbcB[5] = 0.0;
399 newparam->constr.dA = old[0];
400 newparam->constr.dB = old[1];
403 newparam->settle.doh = old[0];
404 newparam->settle.dhh = old[1];
412 newparam->vsite.a = old[0];
413 newparam->vsite.b = old[1];
414 newparam->vsite.c = old[2];
415 newparam->vsite.d = old[3];
416 newparam->vsite.e = old[4];
417 newparam->vsite.f = old[5];
420 newparam->vsite.a = old[1] * cos(DEG2RAD * old[0]);
421 newparam->vsite.b = old[1] * sin(DEG2RAD * old[0]);
422 newparam->vsite.c = old[2];
423 newparam->vsite.d = old[3];
424 newparam->vsite.e = old[4];
425 newparam->vsite.f = old[5];
428 newparam->vsiten.n = round_check(old[0], 1, ftype, "number of atoms");
429 newparam->vsiten.a = old[1];
432 newparam->cmap.cmapA = static_cast<int>(old[0]);
433 newparam->cmap.cmapB = static_cast<int>(old[1]);
435 case F_GB12_NOLONGERUSED:
436 case F_GB13_NOLONGERUSED:
437 case F_GB14_NOLONGERUSED:
440 gmx_fatal(FARGS, "unknown function type %d in %s line %d",
441 ftype, __FILE__, __LINE__);
446 static int enter_params(gmx_ffparams_t *ffparams, t_functype ftype,
447 real forceparams[MAXFORCEPARAM], int comb, real reppow,
448 int start, bool bAppend)
454 if ( (rc = assign_param(ftype, &newparam, forceparams, comb, reppow)) < 0)
456 /* -1 means this interaction is all-zero and should not be added */
462 for (type = start; (type < ffparams->ntypes); type++)
464 if (ffparams->functype[type] == ftype)
466 if (memcmp(&newparam, &ffparams->iparams[type], static_cast<size_t>(sizeof(newparam))) == 0)
475 type = ffparams->ntypes;
477 memcpy(&ffparams->iparams[type], &newparam, static_cast<size_t>(sizeof(newparam)));
480 ffparams->functype[type] = ftype;
485 static void append_interaction(t_ilist *ilist,
486 int type, int nral, const int a[MAXATOMLIST])
493 ilist->iatoms[where1++] = type;
494 for (i = 0; (i < nral); i++)
496 ilist->iatoms[where1++] = a[i];
500 static void enter_function(t_params *p, t_functype ftype, int comb, real reppow,
501 gmx_ffparams_t *ffparams, t_ilist *il,
503 bool bNB, bool bAppend)
505 int k, type, nr, nral, delta, start;
507 start = ffparams->ntypes;
510 for (k = 0; k < nr; k++)
512 if (*maxtypes <= ffparams->ntypes)
515 srenew(ffparams->functype, *maxtypes);
516 srenew(ffparams->iparams, *maxtypes);
518 type = enter_params(ffparams, ftype, p->param[k].c, comb, reppow, start, bAppend);
519 /* Type==-1 is used as a signal that this interaction is all-zero and should not be added. */
520 if (!bNB && type >= 0)
525 srenew(il->iatoms, il->nr+delta);
526 append_interaction(il, type, nral, p->param[k].a);
531 void convert_params(int atnr, t_params nbtypes[],
532 t_molinfo *mi, t_molinfo *intermolecular_interactions,
533 int comb, double reppow, real fudgeQQ,
544 ffp = &mtop->ffparams;
547 ffp->functype = nullptr;
548 ffp->iparams = nullptr;
549 ffp->reppow = reppow;
551 enter_function(&(nbtypes[F_LJ]), static_cast<t_functype>(F_LJ), comb, reppow, ffp, nullptr,
552 &maxtypes, TRUE, TRUE);
553 enter_function(&(nbtypes[F_BHAM]), static_cast<t_functype>(F_BHAM), comb, reppow, ffp, nullptr,
554 &maxtypes, TRUE, TRUE);
556 for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
558 molt = &mtop->moltype[mt];
559 for (i = 0; (i < F_NRE); i++)
561 molt->ilist[i].nr = 0;
562 molt->ilist[i].iatoms = nullptr;
564 plist = mi[mt].plist;
566 flags = interaction_function[i].flags;
567 if ((i != F_LJ) && (i != F_BHAM) && ((flags & IF_BOND) ||
568 (flags & IF_VSITE) ||
569 (flags & IF_CONSTRAINT)))
571 enter_function(&(plist[i]), static_cast<t_functype>(i), comb, reppow,
572 ffp, &molt->ilist[i],
573 &maxtypes, FALSE, (i == F_POSRES || i == F_FBPOSRES));
578 mtop->bIntermolecularInteractions = FALSE;
579 if (intermolecular_interactions != nullptr)
581 /* Process the intermolecular interaction list */
582 snew(mtop->intermolecular_ilist, F_NRE);
584 for (i = 0; (i < F_NRE); i++)
586 mtop->intermolecular_ilist[i].nr = 0;
587 mtop->intermolecular_ilist[i].iatoms = nullptr;
589 plist = intermolecular_interactions->plist;
593 flags = interaction_function[i].flags;
594 /* For intermolecular interactions we (currently)
595 * only support potentials.
596 * Constraints and virtual sites would be possible,
597 * but require a lot of extra (bug-prone) code.
599 if (!(flags & IF_BOND))
601 gmx_fatal(FARGS, "The intermolecular_interaction section may only contain bonded potentials");
603 else if (NRAL(i) == 1) /* e.g. position restraints */
605 gmx_fatal(FARGS, "Single atom interactions don't make sense in the intermolecular_interaction section, you can put them in the moleculetype section");
607 else if (flags & IF_CHEMBOND)
609 gmx_fatal(FARGS, "The intermolecular_interaction can not contain chemically bonding interactions");
613 enter_function(&(plist[i]), static_cast<t_functype>(i), comb, reppow,
614 ffp, &mtop->intermolecular_ilist[i],
615 &maxtypes, FALSE, FALSE);
617 mtop->bIntermolecularInteractions = TRUE;
622 if (!mtop->bIntermolecularInteractions)
624 sfree(mtop->intermolecular_ilist);
628 ffp->fudgeQQ = fudgeQQ;