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,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.
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! */
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"
62 static int round_check(real r, int limit, int ftype, const char *name)
64 const int i = gmx::roundToInt(r);
66 if (r-i > 0.01 || r-i < -0.01)
68 gmx_fatal(FARGS, "A non-integer value (%f) was supplied for '%s' in %s",
69 r, name, interaction_function[ftype].longname);
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);
81 static void set_ljparams(int comb, double reppow, double v, double w,
84 if (comb == eCOMB_ARITHMETIC || comb == eCOMB_GEOM_SIG_EPS)
88 *c6 = 4*w*gmx::power6(v);
89 *c12 = 4*w*std::pow(v, reppow);
93 /* Interpret negative sigma as c6=0 and c12 with -sigma */
95 *c12 = 4*w*std::pow(-v, reppow);
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.
109 assign_param(t_functype ftype, t_iparams *newparam,
110 real old[MAXFORCEPARAM], int comb, double reppow)
113 bool all_param_zero = TRUE;
116 for (j = 0; (j < MAXFORCEPARAM); j++)
118 newparam->generic.buf[j] = 0.0;
119 /* If all parameters are zero we might not add some interaction types (selected below).
120 * We cannot apply this to ALL interactions, since many have valid reasons for having
121 * zero parameters (e.g. an index to a Cmap interaction, or LJ parameters), but
122 * we use it for angles and torsions that are typically generated automatically.
124 all_param_zero = all_param_zero && fabs(old[j]) < GMX_REAL_MIN;
129 if (IS_ANGLE(ftype) || IS_RESTRAINT_TYPE(ftype) || ftype == F_IDIHS ||
130 ftype == F_PDIHS || ftype == F_PIDIHS || ftype == F_RBDIHS || ftype == F_FOURDIHS)
139 /* Post processing of input data: store cosine iso angle itself */
140 newparam->harmonic.rA = cos(old[0]*DEG2RAD);
141 newparam->harmonic.krA = old[1];
142 newparam->harmonic.rB = cos(old[2]*DEG2RAD);
143 newparam->harmonic.krB = old[3];
146 /* Post processing of input data: store square of length itself */
147 newparam->harmonic.rA = gmx::square(old[0]);
148 newparam->harmonic.krA = old[1];
149 newparam->harmonic.rB = gmx::square(old[2]);
150 newparam->harmonic.krB = old[3];
153 newparam->fene.bm = old[0];
154 newparam->fene.kb = old[1];
157 newparam->restraint.lowA = old[0];
158 newparam->restraint.up1A = old[1];
159 newparam->restraint.up2A = old[2];
160 newparam->restraint.kA = old[3];
161 newparam->restraint.lowB = old[4];
162 newparam->restraint.up1B = old[5];
163 newparam->restraint.up2B = old[6];
164 newparam->restraint.kB = old[7];
170 newparam->tab.table = round_check(old[0], 0, ftype, "table index");
171 newparam->tab.kA = old[1];
172 newparam->tab.kB = old[3];
174 case F_CROSS_BOND_BONDS:
175 newparam->cross_bb.r1e = old[0];
176 newparam->cross_bb.r2e = old[1];
177 newparam->cross_bb.krr = old[2];
179 case F_CROSS_BOND_ANGLES:
180 newparam->cross_ba.r1e = old[0];
181 newparam->cross_ba.r2e = old[1];
182 newparam->cross_ba.r3e = old[2];
183 newparam->cross_ba.krt = old[3];
186 newparam->u_b.thetaA = old[0];
187 newparam->u_b.kthetaA = old[1];
188 newparam->u_b.r13A = old[2];
189 newparam->u_b.kUBA = old[3];
190 newparam->u_b.thetaB = old[4];
191 newparam->u_b.kthetaB = old[5];
192 newparam->u_b.r13B = old[6];
193 newparam->u_b.kUBB = old[7];
195 case F_QUARTIC_ANGLES:
196 newparam->qangle.theta = old[0];
197 for (i = 0; i < 5; i++)
199 newparam->qangle.c[i] = old[i+1];
202 case F_LINEAR_ANGLES:
203 newparam->linangle.aA = old[0];
204 newparam->linangle.klinA = old[1];
205 newparam->linangle.aB = old[2];
206 newparam->linangle.klinB = old[3];
212 newparam->harmonic.rA = old[0];
213 newparam->harmonic.krA = old[1];
214 newparam->harmonic.rB = old[2];
215 newparam->harmonic.krB = old[3];
218 newparam->harmonic.rA = old[0];
219 newparam->harmonic.krA = old[1];
222 newparam->morse.b0A = old[0];
223 newparam->morse.cbA = old[1];
224 newparam->morse.betaA = old[2];
225 newparam->morse.b0B = old[3];
226 newparam->morse.cbB = old[4];
227 newparam->morse.betaB = old[5];
230 newparam->cubic.b0 = old[0];
231 newparam->cubic.kb = old[1];
232 newparam->cubic.kcub = old[2];
237 newparam->polarize.alpha = old[0];
240 newparam->anharm_polarize.alpha = old[0];
241 newparam->anharm_polarize.drcut = old[1];
242 newparam->anharm_polarize.khyp = old[2];
245 newparam->wpol.al_x = old[0];
246 newparam->wpol.al_y = old[1];
247 newparam->wpol.al_z = old[2];
248 newparam->wpol.rOH = old[3];
249 newparam->wpol.rHH = old[4];
250 newparam->wpol.rOD = old[5];
253 newparam->thole.a = old[0];
254 newparam->thole.alpha1 = old[1];
255 newparam->thole.alpha2 = old[2];
256 if ((old[1] > 0) && (old[2] > 0))
258 newparam->thole.rfac = old[0]*gmx::invsixthroot(old[1]*old[2]);
262 newparam->thole.rfac = 1;
266 newparam->bham.a = old[0];
267 newparam->bham.b = old[1];
268 newparam->bham.c = old[2];
271 set_ljparams(comb, reppow, old[0], old[1], &newparam->lj14.c6A, &newparam->lj14.c12A);
272 set_ljparams(comb, reppow, old[2], old[3], &newparam->lj14.c6B, &newparam->lj14.c12B);
275 newparam->ljc14.fqq = old[0];
276 newparam->ljc14.qi = old[1];
277 newparam->ljc14.qj = old[2];
278 set_ljparams(comb, reppow, old[3], old[4], &newparam->ljc14.c6, &newparam->ljc14.c12);
281 newparam->ljcnb.qi = old[0];
282 newparam->ljcnb.qj = old[1];
283 set_ljparams(comb, reppow, old[2], old[3], &newparam->ljcnb.c6, &newparam->ljcnb.c12);
286 set_ljparams(comb, reppow, old[0], old[1], &newparam->lj.c6, &newparam->lj.c12);
292 newparam->pdihs.phiA = old[0];
293 newparam->pdihs.cpA = old[1];
295 /* Change 20100720: Amber occasionally uses negative multiplicities (mathematically OK),
296 * so I have changed the lower limit to -99 /EL
298 newparam->pdihs.phiB = old[3];
299 newparam->pdihs.cpB = old[4];
300 /* If both force constants are zero there is no interaction. Return -1 to signal
301 * this entry should NOT be added.
303 if (fabs(newparam->pdihs.cpA) < GMX_REAL_MIN && fabs(newparam->pdihs.cpB) < GMX_REAL_MIN)
308 newparam->pdihs.mult = round_check(old[2], -99, ftype, "multiplicity");
312 newparam->pdihs.phiA = old[0];
313 newparam->pdihs.cpA = old[1];
316 newparam->posres.fcA[XX] = old[0];
317 newparam->posres.fcA[YY] = old[1];
318 newparam->posres.fcA[ZZ] = old[2];
319 newparam->posres.fcB[XX] = old[3];
320 newparam->posres.fcB[YY] = old[4];
321 newparam->posres.fcB[ZZ] = old[5];
322 newparam->posres.pos0A[XX] = old[6];
323 newparam->posres.pos0A[YY] = old[7];
324 newparam->posres.pos0A[ZZ] = old[8];
325 newparam->posres.pos0B[XX] = old[9];
326 newparam->posres.pos0B[YY] = old[10];
327 newparam->posres.pos0B[ZZ] = old[11];
330 newparam->fbposres.geom = round_check(old[0], 0, ftype, "geometry");
331 if (!(newparam->fbposres.geom > efbposresZERO && newparam->fbposres.geom < efbposresNR))
333 gmx_fatal(FARGS, "Invalid geometry for flat-bottomed position restraint.\n"
334 "Expected number between 1 and %d. Found %d\n", efbposresNR-1,
335 newparam->fbposres.geom);
337 newparam->fbposres.r = old[1];
338 newparam->fbposres.k = old[2];
339 newparam->fbposres.pos0[XX] = old[3];
340 newparam->fbposres.pos0[YY] = old[4];
341 newparam->fbposres.pos0[ZZ] = old[5];
344 newparam->disres.label = round_check(old[0], 0, ftype, "label");
345 newparam->disres.type = round_check(old[1], 1, ftype, "type'");
346 newparam->disres.low = old[2];
347 newparam->disres.up1 = old[3];
348 newparam->disres.up2 = old[4];
349 newparam->disres.kfac = old[5];
352 newparam->orires.ex = round_check(old[0], 1, ftype, "experiment") - 1;
353 newparam->orires.label = round_check(old[1], 1, ftype, "label");
354 newparam->orires.power = round_check(old[2], 0, ftype, "power");
355 newparam->orires.c = old[3];
356 newparam->orires.obs = old[4];
357 newparam->orires.kfac = old[5];
360 newparam->dihres.phiA = old[0];
361 newparam->dihres.dphiA = old[1];
362 newparam->dihres.kfacA = old[2];
363 newparam->dihres.phiB = old[3];
364 newparam->dihres.dphiB = old[4];
365 newparam->dihres.kfacB = old[5];
368 for (i = 0; (i < NR_RBDIHS); i++)
370 newparam->rbdihs.rbcA[i] = old[i];
371 newparam->rbdihs.rbcB[i] = old[NR_RBDIHS+i];
375 for (i = 0; (i < NR_CBTDIHS); i++)
377 newparam->cbtdihs.cbtcA[i] = old[i];
381 /* Read the dihedral parameters to temporary arrays,
382 * and convert them to the computationally faster
383 * Ryckaert-Bellemans form.
385 /* Use conversion formula for OPLS to Ryckaert-Bellemans: */
386 newparam->rbdihs.rbcA[0] = old[1]+0.5*(old[0]+old[2]);
387 newparam->rbdihs.rbcA[1] = 0.5*(3.0*old[2]-old[0]);
388 newparam->rbdihs.rbcA[2] = 4.0*old[3]-old[1];
389 newparam->rbdihs.rbcA[3] = -2.0*old[2];
390 newparam->rbdihs.rbcA[4] = -4.0*old[3];
391 newparam->rbdihs.rbcA[5] = 0.0;
393 newparam->rbdihs.rbcB[0] = old[NR_FOURDIHS+1]+0.5*(old[NR_FOURDIHS+0]+old[NR_FOURDIHS+2]);
394 newparam->rbdihs.rbcB[1] = 0.5*(3.0*old[NR_FOURDIHS+2]-old[NR_FOURDIHS+0]);
395 newparam->rbdihs.rbcB[2] = 4.0*old[NR_FOURDIHS+3]-old[NR_FOURDIHS+1];
396 newparam->rbdihs.rbcB[3] = -2.0*old[NR_FOURDIHS+2];
397 newparam->rbdihs.rbcB[4] = -4.0*old[NR_FOURDIHS+3];
398 newparam->rbdihs.rbcB[5] = 0.0;
402 newparam->constr.dA = old[0];
403 newparam->constr.dB = old[1];
406 newparam->settle.doh = old[0];
407 newparam->settle.dhh = old[1];
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];
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];
431 newparam->vsiten.n = round_check(old[0], 1, ftype, "number of atoms");
432 newparam->vsiten.a = old[1];
435 newparam->cmap.cmapA = static_cast<int>(old[0]);
436 newparam->cmap.cmapB = static_cast<int>(old[1]);
438 case F_GB12_NOLONGERUSED:
439 case F_GB13_NOLONGERUSED:
440 case F_GB14_NOLONGERUSED:
443 gmx_fatal(FARGS, "unknown function type %d in %s line %d",
444 ftype, __FILE__, __LINE__);
449 static int enter_params(gmx_ffparams_t *ffparams, t_functype ftype,
450 real forceparams[MAXFORCEPARAM], int comb, real reppow,
451 int start, bool bAppend)
457 if ( (rc = assign_param(ftype, &newparam, forceparams, comb, reppow)) < 0)
459 /* -1 means this interaction is all-zero and should not be added */
465 for (type = start; (type < ffparams->numTypes()); type++)
467 if (ffparams->functype[type] == ftype)
469 if (memcmp(&newparam, &ffparams->iparams[type], static_cast<size_t>(sizeof(newparam))) == 0)
478 type = ffparams->numTypes();
481 ffparams->iparams.push_back(newparam);
482 ffparams->functype.push_back(ftype);
484 GMX_ASSERT(ffparams->iparams.size() == ffparams->functype.size(),
485 "sizes should match");
490 static void append_interaction(InteractionList *ilist,
491 int type, int nral, const int a[MAXATOMLIST])
493 ilist->iatoms.push_back(type);
494 for (int i = 0; (i < nral); i++)
496 ilist->iatoms.push_back(a[i]);
500 static void enter_function(t_params *p, t_functype ftype, int comb, real reppow,
501 gmx_ffparams_t *ffparams, InteractionList *il,
502 bool bNB, bool bAppend)
504 int k, type, nr, nral, start;
506 start = ffparams->numTypes();
509 for (k = 0; k < nr; k++)
511 type = enter_params(ffparams, ftype, p->param[k].c, comb, reppow, start, bAppend);
512 /* Type==-1 is used as a signal that this interaction is all-zero and should not be added. */
513 if (!bNB && type >= 0)
517 append_interaction(il, type, nral, p->param[k].a);
522 void convert_params(int atnr, t_params nbtypes[],
523 t_molinfo *mi, t_molinfo *intermolecular_interactions,
524 int comb, double reppow, real fudgeQQ,
533 ffp = &mtop->ffparams;
535 ffp->functype.clear();
536 ffp->iparams.clear();
537 ffp->reppow = reppow;
539 enter_function(&(nbtypes[F_LJ]), static_cast<t_functype>(F_LJ), comb, reppow, ffp, nullptr,
541 enter_function(&(nbtypes[F_BHAM]), static_cast<t_functype>(F_BHAM), comb, reppow, ffp, nullptr,
544 for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
546 molt = &mtop->moltype[mt];
547 for (i = 0; (i < F_NRE); i++)
549 molt->ilist[i].iatoms.clear();
551 plist = mi[mt].plist;
553 flags = interaction_function[i].flags;
554 if ((i != F_LJ) && (i != F_BHAM) && ((flags & IF_BOND) ||
555 (flags & IF_VSITE) ||
556 (flags & IF_CONSTRAINT)))
558 enter_function(&(plist[i]), static_cast<t_functype>(i), comb, reppow,
559 ffp, &molt->ilist[i],
560 FALSE, (i == F_POSRES || i == F_FBPOSRES));
565 mtop->bIntermolecularInteractions = FALSE;
566 if (intermolecular_interactions != nullptr)
568 /* Process the intermolecular interaction list */
569 mtop->intermolecular_ilist = std::make_unique<InteractionLists>();
571 for (i = 0; (i < F_NRE); i++)
573 (*mtop->intermolecular_ilist)[i].iatoms.clear();
575 plist = intermolecular_interactions->plist;
579 flags = interaction_function[i].flags;
580 /* For intermolecular interactions we (currently)
581 * only support potentials.
582 * Constraints and virtual sites would be possible,
583 * but require a lot of extra (bug-prone) code.
585 if (!(flags & IF_BOND))
587 gmx_fatal(FARGS, "The intermolecular_interaction section may only contain bonded potentials");
589 else if (NRAL(i) == 1) /* e.g. position restraints */
591 gmx_fatal(FARGS, "Single atom interactions don't make sense in the intermolecular_interaction section, you can put them in the moleculetype section");
593 else if (flags & IF_CHEMBOND)
595 gmx_fatal(FARGS, "The intermolecular_interaction can not contain chemically bonding interactions");
599 enter_function(&(plist[i]), static_cast<t_functype>(i), comb, reppow,
600 ffp, &(*mtop->intermolecular_ilist)[i],
603 mtop->bIntermolecularInteractions = TRUE;
608 if (!mtop->bIntermolecularInteractions)
610 mtop->intermolecular_ilist.reset(nullptr);
614 ffp->fudgeQQ = fudgeQQ;