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", r, name,
69 interaction_function[ftype].longname);
74 gmx_fatal(FARGS, "Value of '%s' in %s is %d, which is smaller than the minimum of %d", name,
75 interaction_function[ftype].longname, i, limit);
81 static void set_ljparams(int comb, double reppow, double v, double w, real* c6, real* c12)
83 if (comb == eCOMB_ARITHMETIC || comb == eCOMB_GEOM_SIG_EPS)
87 *c6 = 4 * w * gmx::power6(v);
88 *c12 = 4 * w * std::pow(v, reppow);
92 /* Interpret negative sigma as c6=0 and c12 with -sigma */
94 *c12 = 4 * w * std::pow(-v, reppow);
104 /* A return value of 0 means parameters were assigned successfully,
105 * returning -1 means this is an all-zero interaction that should not be added.
107 static int assign_param(t_functype ftype, t_iparams* newparam, gmx::ArrayRef<const real> old, int comb, double reppow)
109 bool all_param_zero = true;
112 for (int j = 0; (j < MAXFORCEPARAM); j++)
114 newparam->generic.buf[j] = 0.0;
115 /* If all parameters are zero we might not add some interaction types (selected below).
116 * We cannot apply this to ALL interactions, since many have valid reasons for having
117 * zero parameters (e.g. an index to a Cmap interaction, or LJ parameters), but
118 * we use it for angles and torsions that are typically generated automatically.
120 all_param_zero = all_param_zero && fabs(old[j]) < GMX_REAL_MIN;
125 if (IS_ANGLE(ftype) || IS_RESTRAINT_TYPE(ftype) || ftype == F_IDIHS || ftype == F_PDIHS
126 || ftype == F_PIDIHS || ftype == F_RBDIHS || ftype == F_FOURDIHS)
135 /* Post processing of input data: store cosine iso angle itself */
136 newparam->harmonic.rA = cos(old[0] * DEG2RAD);
137 newparam->harmonic.krA = old[1];
138 newparam->harmonic.rB = cos(old[2] * DEG2RAD);
139 newparam->harmonic.krB = old[3];
142 /* Post processing of input data: store square of length itself */
143 newparam->harmonic.rA = gmx::square(old[0]);
144 newparam->harmonic.krA = old[1];
145 newparam->harmonic.rB = gmx::square(old[2]);
146 newparam->harmonic.krB = old[3];
149 newparam->fene.bm = old[0];
150 newparam->fene.kb = old[1];
153 newparam->restraint.lowA = old[0];
154 newparam->restraint.up1A = old[1];
155 newparam->restraint.up2A = old[2];
156 newparam->restraint.kA = old[3];
157 newparam->restraint.lowB = old[4];
158 newparam->restraint.up1B = old[5];
159 newparam->restraint.up2B = old[6];
160 newparam->restraint.kB = old[7];
166 newparam->tab.table = round_check(old[0], 0, ftype, "table index");
167 newparam->tab.kA = old[1];
168 newparam->tab.kB = old[3];
170 case F_CROSS_BOND_BONDS:
171 newparam->cross_bb.r1e = old[0];
172 newparam->cross_bb.r2e = old[1];
173 newparam->cross_bb.krr = old[2];
175 case F_CROSS_BOND_ANGLES:
176 newparam->cross_ba.r1e = old[0];
177 newparam->cross_ba.r2e = old[1];
178 newparam->cross_ba.r3e = old[2];
179 newparam->cross_ba.krt = old[3];
182 newparam->u_b.thetaA = old[0];
183 newparam->u_b.kthetaA = old[1];
184 newparam->u_b.r13A = old[2];
185 newparam->u_b.kUBA = old[3];
186 newparam->u_b.thetaB = old[4];
187 newparam->u_b.kthetaB = old[5];
188 newparam->u_b.r13B = old[6];
189 newparam->u_b.kUBB = old[7];
191 case F_QUARTIC_ANGLES:
192 newparam->qangle.theta = old[0];
193 for (int i = 0; i < 5; i++)
195 newparam->qangle.c[i] = old[i + 1];
198 case F_LINEAR_ANGLES:
199 newparam->linangle.aA = old[0];
200 newparam->linangle.klinA = old[1];
201 newparam->linangle.aB = old[2];
202 newparam->linangle.klinB = old[3];
208 newparam->harmonic.rA = old[0];
209 newparam->harmonic.krA = old[1];
210 newparam->harmonic.rB = old[2];
211 newparam->harmonic.krB = old[3];
214 newparam->harmonic.rA = old[0];
215 newparam->harmonic.krA = old[1];
218 newparam->morse.b0A = old[0];
219 newparam->morse.cbA = old[1];
220 newparam->morse.betaA = old[2];
221 newparam->morse.b0B = old[3];
222 newparam->morse.cbB = old[4];
223 newparam->morse.betaB = old[5];
226 newparam->cubic.b0 = old[0];
227 newparam->cubic.kb = old[1];
228 newparam->cubic.kcub = old[2];
230 case F_CONNBONDS: break;
231 case F_POLARIZATION: newparam->polarize.alpha = old[0]; break;
233 newparam->anharm_polarize.alpha = old[0];
234 newparam->anharm_polarize.drcut = old[1];
235 newparam->anharm_polarize.khyp = old[2];
238 newparam->wpol.al_x = old[0];
239 newparam->wpol.al_y = old[1];
240 newparam->wpol.al_z = old[2];
241 newparam->wpol.rOH = old[3];
242 newparam->wpol.rHH = old[4];
243 newparam->wpol.rOD = old[5];
246 newparam->thole.a = old[0];
247 newparam->thole.alpha1 = old[1];
248 newparam->thole.alpha2 = old[2];
249 if ((old[1] > 0) && (old[2] > 0))
251 newparam->thole.rfac = old[0] * gmx::invsixthroot(old[1] * old[2]);
255 newparam->thole.rfac = 1;
259 newparam->bham.a = old[0];
260 newparam->bham.b = old[1];
261 newparam->bham.c = old[2];
264 set_ljparams(comb, reppow, old[0], old[1], &newparam->lj14.c6A, &newparam->lj14.c12A);
265 set_ljparams(comb, reppow, old[2], old[3], &newparam->lj14.c6B, &newparam->lj14.c12B);
268 newparam->ljc14.fqq = old[0];
269 newparam->ljc14.qi = old[1];
270 newparam->ljc14.qj = old[2];
271 set_ljparams(comb, reppow, old[3], old[4], &newparam->ljc14.c6, &newparam->ljc14.c12);
274 newparam->ljcnb.qi = old[0];
275 newparam->ljcnb.qj = old[1];
276 set_ljparams(comb, reppow, old[2], old[3], &newparam->ljcnb.c6, &newparam->ljcnb.c12);
279 set_ljparams(comb, reppow, old[0], old[1], &newparam->lj.c6, &newparam->lj.c12);
285 newparam->pdihs.phiA = old[0];
286 newparam->pdihs.cpA = old[1];
288 /* Change 20100720: Amber occasionally uses negative multiplicities (mathematically OK),
289 * so I have changed the lower limit to -99 /EL
291 newparam->pdihs.phiB = old[3];
292 newparam->pdihs.cpB = old[4];
293 /* If both force constants are zero there is no interaction. Return -1 to signal
294 * this entry should NOT be added.
296 if (fabs(newparam->pdihs.cpA) < GMX_REAL_MIN && fabs(newparam->pdihs.cpB) < GMX_REAL_MIN)
301 newparam->pdihs.mult = round_check(old[2], -99, ftype, "multiplicity");
305 newparam->pdihs.phiA = old[0];
306 newparam->pdihs.cpA = old[1];
309 newparam->posres.fcA[XX] = old[0];
310 newparam->posres.fcA[YY] = old[1];
311 newparam->posres.fcA[ZZ] = old[2];
312 newparam->posres.fcB[XX] = old[3];
313 newparam->posres.fcB[YY] = old[4];
314 newparam->posres.fcB[ZZ] = old[5];
315 newparam->posres.pos0A[XX] = old[6];
316 newparam->posres.pos0A[YY] = old[7];
317 newparam->posres.pos0A[ZZ] = old[8];
318 newparam->posres.pos0B[XX] = old[9];
319 newparam->posres.pos0B[YY] = old[10];
320 newparam->posres.pos0B[ZZ] = old[11];
323 newparam->fbposres.geom = round_check(old[0], 0, ftype, "geometry");
324 if (!(newparam->fbposres.geom > efbposresZERO && newparam->fbposres.geom < efbposresNR))
327 "Invalid geometry for flat-bottomed position restraint.\n"
328 "Expected number between 1 and %d. Found %d\n",
329 efbposresNR - 1, newparam->fbposres.geom);
331 newparam->fbposres.r = old[1];
332 newparam->fbposres.k = old[2];
333 newparam->fbposres.pos0[XX] = old[3];
334 newparam->fbposres.pos0[YY] = old[4];
335 newparam->fbposres.pos0[ZZ] = old[5];
338 newparam->disres.label = round_check(old[0], 0, ftype, "label");
339 newparam->disres.type = round_check(old[1], 1, ftype, "type'");
340 newparam->disres.low = old[2];
341 newparam->disres.up1 = old[3];
342 newparam->disres.up2 = old[4];
343 newparam->disres.kfac = old[5];
346 newparam->orires.ex = round_check(old[0], 1, ftype, "experiment") - 1;
347 newparam->orires.label = round_check(old[1], 1, ftype, "label");
348 newparam->orires.power = round_check(old[2], 0, ftype, "power");
349 newparam->orires.c = old[3];
350 newparam->orires.obs = old[4];
351 newparam->orires.kfac = old[5];
354 newparam->dihres.phiA = old[0];
355 newparam->dihres.dphiA = old[1];
356 newparam->dihres.kfacA = old[2];
357 newparam->dihres.phiB = old[3];
358 newparam->dihres.dphiB = old[4];
359 newparam->dihres.kfacB = old[5];
362 for (int i = 0; (i < NR_RBDIHS); i++)
364 newparam->rbdihs.rbcA[i] = old[i];
365 newparam->rbdihs.rbcB[i] = old[NR_RBDIHS + i];
369 for (int i = 0; (i < NR_CBTDIHS); i++)
371 newparam->cbtdihs.cbtcA[i] = old[i];
375 /* Read the dihedral parameters to temporary arrays,
376 * and convert them to the computationally faster
377 * Ryckaert-Bellemans form.
379 /* Use conversion formula for OPLS to Ryckaert-Bellemans: */
380 newparam->rbdihs.rbcA[0] = old[1] + 0.5 * (old[0] + old[2]);
381 newparam->rbdihs.rbcA[1] = 0.5 * (3.0 * old[2] - old[0]);
382 newparam->rbdihs.rbcA[2] = 4.0 * old[3] - old[1];
383 newparam->rbdihs.rbcA[3] = -2.0 * old[2];
384 newparam->rbdihs.rbcA[4] = -4.0 * old[3];
385 newparam->rbdihs.rbcA[5] = 0.0;
387 newparam->rbdihs.rbcB[0] =
388 old[NR_FOURDIHS + 1] + 0.5 * (old[NR_FOURDIHS + 0] + old[NR_FOURDIHS + 2]);
389 newparam->rbdihs.rbcB[1] = 0.5 * (3.0 * old[NR_FOURDIHS + 2] - old[NR_FOURDIHS + 0]);
390 newparam->rbdihs.rbcB[2] = 4.0 * old[NR_FOURDIHS + 3] - old[NR_FOURDIHS + 1];
391 newparam->rbdihs.rbcB[3] = -2.0 * old[NR_FOURDIHS + 2];
392 newparam->rbdihs.rbcB[4] = -4.0 * old[NR_FOURDIHS + 3];
393 newparam->rbdihs.rbcB[5] = 0.0;
397 newparam->constr.dA = old[0];
398 newparam->constr.dB = old[1];
401 newparam->settle.doh = old[0];
402 newparam->settle.dhh = old[1];
411 newparam->vsite.a = old[0];
412 newparam->vsite.b = old[1];
413 newparam->vsite.c = old[2];
414 newparam->vsite.d = old[3];
415 newparam->vsite.e = old[4];
416 newparam->vsite.f = old[5];
419 newparam->vsite.a = old[1] * cos(DEG2RAD * old[0]);
420 newparam->vsite.b = old[1] * sin(DEG2RAD * old[0]);
421 newparam->vsite.c = old[2];
422 newparam->vsite.d = old[3];
423 newparam->vsite.e = old[4];
424 newparam->vsite.f = old[5];
427 newparam->vsiten.n = round_check(old[0], 1, ftype, "number of atoms");
428 newparam->vsiten.a = old[1];
431 newparam->cmap.cmapA = static_cast<int>(old[0]);
432 newparam->cmap.cmapB = static_cast<int>(old[1]);
434 case F_GB12_NOLONGERUSED:
435 case F_GB13_NOLONGERUSED:
436 case F_GB14_NOLONGERUSED: break;
438 gmx_fatal(FARGS, "unknown function type %d in %s line %d", ftype, __FILE__, __LINE__);
443 static int enter_params(gmx_ffparams_t* ffparams,
445 gmx::ArrayRef<const real> forceparams,
455 if ((rc = assign_param(ftype, &newparam, forceparams, comb, reppow)) < 0)
457 /* -1 means this interaction is all-zero and should not be added */
463 for (type = start; (type < ffparams->numTypes()); type++)
465 if (ffparams->functype[type] == ftype)
467 if (memcmp(&newparam, &ffparams->iparams[type], static_cast<size_t>(sizeof(newparam))) == 0)
476 type = ffparams->numTypes();
479 ffparams->iparams.push_back(newparam);
480 ffparams->functype.push_back(ftype);
482 GMX_ASSERT(ffparams->iparams.size() == ffparams->functype.size(), "sizes should match");
487 static void append_interaction(InteractionList* ilist, int type, gmx::ArrayRef<const int> a)
489 ilist->iatoms.push_back(type);
490 for (const auto& atom : a)
492 ilist->iatoms.push_back(atom);
496 static void enter_function(const InteractionsOfType* p,
500 gmx_ffparams_t* ffparams,
505 int start = ffparams->numTypes();
507 for (auto& parm : p->interactionTypes)
509 int type = enter_params(ffparams, ftype, parm.forceParam(), comb, reppow, start, bAppend);
510 /* Type==-1 is used as a signal that this interaction is all-zero and should not be added. */
511 if (!bNB && type >= 0)
513 GMX_RELEASE_ASSERT(il, "Need valid interaction list");
514 GMX_RELEASE_ASSERT(parm.atoms().ssize() == NRAL(ftype),
515 "Need to have correct number of atoms for the parameter");
516 append_interaction(il, type, parm.atoms());
521 void convertInteractionsOfType(int atnr,
522 gmx::ArrayRef<const InteractionsOfType> nbtypes,
523 gmx::ArrayRef<const MoleculeInformation> mi,
524 const MoleculeInformation* intermolecular_interactions,
535 ffp = &mtop->ffparams;
537 ffp->functype.clear();
538 ffp->iparams.clear();
539 ffp->reppow = reppow;
541 enter_function(&(nbtypes[F_LJ]), static_cast<t_functype>(F_LJ), comb, reppow, ffp, nullptr, TRUE, TRUE);
542 enter_function(&(nbtypes[F_BHAM]), static_cast<t_functype>(F_BHAM), comb, reppow, ffp, nullptr,
545 for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
547 molt = &mtop->moltype[mt];
548 for (i = 0; (i < F_NRE); i++)
550 molt->ilist[i].iatoms.clear();
552 gmx::ArrayRef<const InteractionsOfType> interactions = mi[mt].interactions;
554 flags = interaction_function[i].flags;
555 if ((i != F_LJ) && (i != F_BHAM)
556 && ((flags & IF_BOND) || (flags & IF_VSITE) || (flags & IF_CONSTRAINT)))
558 enter_function(&(interactions[i]), static_cast<t_functype>(i), comb, reppow, ffp,
559 &molt->ilist[i], FALSE, (i == F_POSRES || i == F_FBPOSRES));
564 mtop->bIntermolecularInteractions = FALSE;
565 if (intermolecular_interactions != nullptr)
567 /* Process the intermolecular interaction list */
568 mtop->intermolecular_ilist = std::make_unique<InteractionLists>();
570 for (i = 0; (i < F_NRE); i++)
572 (*mtop->intermolecular_ilist)[i].iatoms.clear();
574 gmx::ArrayRef<const InteractionsOfType> interactions = intermolecular_interactions->interactions;
576 if (!interactions[i].interactionTypes.empty())
578 flags = interaction_function[i].flags;
579 /* For intermolecular interactions we (currently)
580 * only support potentials.
581 * Constraints and virtual sites would be possible,
582 * but require a lot of extra (bug-prone) code.
584 if (!(flags & IF_BOND))
587 "The intermolecular_interaction section may only contain bonded "
590 else if (NRAL(i) == 1) /* e.g. position restraints */
593 "Single atom interactions don't make sense in the "
594 "intermolecular_interaction section, you can put them in the "
595 "moleculetype section");
597 else if (flags & IF_CHEMBOND)
600 "The intermolecular_interaction can not contain chemically bonding "
605 enter_function(&(interactions[i]), static_cast<t_functype>(i), comb, reppow,
606 ffp, &(*mtop->intermolecular_ilist)[i], FALSE, FALSE);
608 mtop->bIntermolecularInteractions = TRUE;
613 if (!mtop->bIntermolecularInteractions)
615 mtop->intermolecular_ilist.reset(nullptr);
619 ffp->fudgeQQ = fudgeQQ;