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 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 /* This file is completely threadsafe - keep it that way! */
49 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
50 #include "gromacs/gmxpreprocess/grompp_impl.h"
51 #include "gromacs/gmxpreprocess/topio.h"
52 #include "gromacs/gmxpreprocess/toputil.h"
53 #include "gromacs/math/functions.h"
54 #include "gromacs/math/units.h"
55 #include "gromacs/math/utilities.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/mdtypes/md_enums.h"
58 #include "gromacs/topology/ifunc.h"
59 #include "gromacs/topology/topology.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/smalloc.h"
63 static int round_check(real r, int limit, int ftype, const char* name)
65 const int i = gmx::roundToInt(r);
67 if (r - i > 0.01 || r - i < -0.01)
70 "A non-integer value (%f) was supplied for '%s' in %s",
73 interaction_function[ftype].longname);
79 "Value of '%s' in %s is %d, which is smaller than the minimum of %d",
81 interaction_function[ftype].longname,
89 static void set_ljparams(CombinationRule comb, double reppow, double v, double w, real* c6, real* c12)
91 if (comb == CombinationRule::Arithmetic || comb == CombinationRule::GeomSigEps)
95 *c6 = 4 * w * gmx::power6(v);
96 *c12 = 4 * w * std::pow(v, reppow);
100 /* Interpret negative sigma as c6=0 and c12 with -sigma */
102 *c12 = 4 * w * std::pow(-v, reppow);
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.
115 static int assign_param(t_functype ftype,
117 gmx::ArrayRef<const real> old,
118 CombinationRule comb,
121 bool all_param_zero = true;
124 for (int j = 0; (j < MAXFORCEPARAM); j++)
126 newparam->generic.buf[j] = 0.0;
127 /* If all parameters are zero we might not add some interaction types (selected below).
128 * We cannot apply this to ALL interactions, since many have valid reasons for having
129 * zero parameters (e.g. an index to a Cmap interaction, or LJ parameters), but
130 * we use it for angles and torsions that are typically generated automatically.
132 all_param_zero = all_param_zero && fabs(old[j]) < GMX_REAL_MIN;
137 if (IS_ANGLE(ftype) || IS_RESTRAINT_TYPE(ftype) || ftype == F_IDIHS || ftype == F_PDIHS
138 || ftype == F_PIDIHS || ftype == F_RBDIHS || ftype == F_FOURDIHS)
147 /* Post processing of input data: store cosine iso angle itself */
148 newparam->harmonic.rA = cos(old[0] * gmx::c_deg2Rad);
149 newparam->harmonic.krA = old[1];
150 newparam->harmonic.rB = cos(old[2] * gmx::c_deg2Rad);
151 newparam->harmonic.krB = old[3];
154 /* Post processing of input data: store square of length itself */
155 newparam->harmonic.rA = gmx::square(old[0]);
156 newparam->harmonic.krA = old[1];
157 newparam->harmonic.rB = gmx::square(old[2]);
158 newparam->harmonic.krB = old[3];
161 newparam->fene.bm = old[0];
162 newparam->fene.kb = old[1];
165 newparam->restraint.lowA = old[0];
166 newparam->restraint.up1A = old[1];
167 newparam->restraint.up2A = old[2];
168 newparam->restraint.kA = old[3];
169 newparam->restraint.lowB = old[4];
170 newparam->restraint.up1B = old[5];
171 newparam->restraint.up2B = old[6];
172 newparam->restraint.kB = old[7];
178 newparam->tab.table = round_check(old[0], 0, ftype, "table index");
179 newparam->tab.kA = old[1];
180 newparam->tab.kB = old[3];
182 case F_CROSS_BOND_BONDS:
183 newparam->cross_bb.r1e = old[0];
184 newparam->cross_bb.r2e = old[1];
185 newparam->cross_bb.krr = old[2];
187 case F_CROSS_BOND_ANGLES:
188 newparam->cross_ba.r1e = old[0];
189 newparam->cross_ba.r2e = old[1];
190 newparam->cross_ba.r3e = old[2];
191 newparam->cross_ba.krt = old[3];
194 newparam->u_b.thetaA = old[0];
195 newparam->u_b.kthetaA = old[1];
196 newparam->u_b.r13A = old[2];
197 newparam->u_b.kUBA = old[3];
198 newparam->u_b.thetaB = old[4];
199 newparam->u_b.kthetaB = old[5];
200 newparam->u_b.r13B = old[6];
201 newparam->u_b.kUBB = old[7];
203 case F_QUARTIC_ANGLES:
204 newparam->qangle.theta = old[0];
205 for (int i = 0; i < 5; i++)
207 newparam->qangle.c[i] = old[i + 1];
210 case F_LINEAR_ANGLES:
211 newparam->linangle.aA = old[0];
212 newparam->linangle.klinA = old[1];
213 newparam->linangle.aB = old[2];
214 newparam->linangle.klinB = old[3];
220 newparam->harmonic.rA = old[0];
221 newparam->harmonic.krA = old[1];
222 newparam->harmonic.rB = old[2];
223 newparam->harmonic.krB = old[3];
226 newparam->harmonic.rA = old[0];
227 newparam->harmonic.krA = old[1];
230 newparam->morse.b0A = old[0];
231 newparam->morse.cbA = old[1];
232 newparam->morse.betaA = old[2];
233 newparam->morse.b0B = old[3];
234 newparam->morse.cbB = old[4];
235 newparam->morse.betaB = old[5];
238 newparam->cubic.b0 = old[0];
239 newparam->cubic.kb = old[1];
240 newparam->cubic.kcub = old[2];
242 case F_CONNBONDS: break;
243 case F_POLARIZATION: newparam->polarize.alpha = old[0]; break;
245 newparam->anharm_polarize.alpha = old[0];
246 newparam->anharm_polarize.drcut = old[1];
247 newparam->anharm_polarize.khyp = old[2];
250 newparam->wpol.al_x = old[0];
251 newparam->wpol.al_y = old[1];
252 newparam->wpol.al_z = old[2];
253 newparam->wpol.rOH = old[3];
254 newparam->wpol.rHH = old[4];
255 newparam->wpol.rOD = old[5];
258 newparam->thole.a = old[0];
259 newparam->thole.alpha1 = old[1];
260 newparam->thole.alpha2 = old[2];
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))
331 "Invalid geometry for flat-bottomed position restraint.\n"
332 "Expected number between 1 and %d. Found %d\n",
334 newparam->fbposres.geom);
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];
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];
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];
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];
367 for (int i = 0; (i < NR_RBDIHS); i++)
369 newparam->rbdihs.rbcA[i] = old[i];
370 newparam->rbdihs.rbcB[i] = old[NR_RBDIHS + i];
374 for (int i = 0; (i < NR_CBTDIHS); i++)
376 newparam->cbtdihs.cbtcA[i] = old[i];
380 /* Read the dihedral parameters to temporary arrays,
381 * and convert them to the computationally faster
382 * Ryckaert-Bellemans form.
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;
392 newparam->rbdihs.rbcB[0] =
393 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];
417 newparam->vsite.a = old[0];
418 newparam->vsite.b = old[1];
419 newparam->vsite.c = old[2];
420 newparam->vsite.d = old[3];
421 newparam->vsite.e = old[4];
422 newparam->vsite.f = old[5];
425 newparam->vsite.a = old[1] * cos(gmx::c_deg2Rad * old[0]);
426 newparam->vsite.b = old[1] * sin(gmx::c_deg2Rad * old[0]);
427 newparam->vsite.c = old[2];
428 newparam->vsite.d = old[3];
429 newparam->vsite.e = old[4];
430 newparam->vsite.f = old[5];
433 newparam->vsiten.n = round_check(old[0], 1, ftype, "number of atoms");
434 newparam->vsiten.a = old[1];
437 newparam->cmap.cmapA = static_cast<int>(old[0]);
438 newparam->cmap.cmapB = static_cast<int>(old[1]);
440 case F_GB12_NOLONGERUSED:
441 case F_GB13_NOLONGERUSED:
442 case F_GB14_NOLONGERUSED: break;
444 gmx_fatal(FARGS, "unknown function type %d in %s line %d", ftype, __FILE__, __LINE__);
449 static int enter_params(gmx_ffparams_t* ffparams,
451 gmx::ArrayRef<const real> forceparams,
452 CombinationRule comb,
460 if ((rc = assign_param(ftype, &newparam, forceparams, comb, reppow)) < 0)
462 /* -1 means this interaction is all-zero and should not be added */
468 if (ftype != F_DISRES)
470 for (int type = start; type < ffparams->numTypes(); type++)
472 // Note that the first condition is always met by starting the loop at start
473 if (ffparams->functype[type] == ftype
474 && memcmp(&newparam, &ffparams->iparams[type], static_cast<size_t>(sizeof(newparam))) == 0)
482 // Distance restraints should have unique labels and pairs with the same label
483 // should be consecutive, so we here we only need to check the last type in the list.
484 // This changes the complexity from quadratic to linear in the number of restraints.
485 const int type = ffparams->numTypes() - 1;
486 if (type >= 0 && ffparams->functype[type] == ftype
487 && memcmp(&newparam, &ffparams->iparams[type], static_cast<size_t>(sizeof(newparam))) == 0)
494 const int type = ffparams->numTypes();
496 ffparams->iparams.push_back(newparam);
497 ffparams->functype.push_back(ftype);
499 GMX_ASSERT(ffparams->iparams.size() == ffparams->functype.size(), "sizes should match");
504 static void append_interaction(InteractionList* ilist, int type, gmx::ArrayRef<const int> a)
506 ilist->iatoms.push_back(type);
507 for (const auto& atom : a)
509 ilist->iatoms.push_back(atom);
513 static void enter_function(const InteractionsOfType* p,
515 CombinationRule comb,
517 gmx_ffparams_t* ffparams,
522 int start = ffparams->numTypes();
524 for (const auto& parm : p->interactionTypes)
526 int type = enter_params(ffparams, ftype, parm.forceParam(), comb, reppow, start, bAppend);
527 /* Type==-1 is used as a signal that this interaction is all-zero and should not be added. */
528 if (!bNB && type >= 0)
530 GMX_RELEASE_ASSERT(il, "Need valid interaction list");
531 GMX_RELEASE_ASSERT(parm.atoms().ssize() == NRAL(ftype),
532 "Need to have correct number of atoms for the parameter");
533 append_interaction(il, type, parm.atoms());
538 void convertInteractionsOfType(int atnr,
539 gmx::ArrayRef<const InteractionsOfType> nbtypes,
540 gmx::ArrayRef<const MoleculeInformation> mi,
541 const MoleculeInformation* intermolecular_interactions,
542 CombinationRule comb,
552 ffp = &mtop->ffparams;
554 ffp->functype.clear();
555 ffp->iparams.clear();
556 ffp->reppow = reppow;
558 enter_function(&(nbtypes[F_LJ]), static_cast<t_functype>(F_LJ), comb, reppow, ffp, nullptr, TRUE, TRUE);
560 &(nbtypes[F_BHAM]), static_cast<t_functype>(F_BHAM), comb, reppow, ffp, nullptr, TRUE, TRUE);
562 for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
564 molt = &mtop->moltype[mt];
565 for (i = 0; (i < F_NRE); i++)
567 molt->ilist[i].iatoms.clear();
569 gmx::ArrayRef<const InteractionsOfType> interactions = mi[mt].interactions;
571 flags = interaction_function[i].flags;
572 if ((i != F_LJ) && (i != F_BHAM)
573 && ((flags & IF_BOND) || (flags & IF_VSITE) || (flags & IF_CONSTRAINT)))
575 enter_function(&(interactions[i]),
576 static_cast<t_functype>(i),
582 (i == F_POSRES || i == F_FBPOSRES));
587 mtop->bIntermolecularInteractions = FALSE;
588 if (intermolecular_interactions != nullptr)
590 /* Process the intermolecular interaction list */
591 mtop->intermolecular_ilist = std::make_unique<InteractionLists>();
593 for (i = 0; (i < F_NRE); i++)
595 (*mtop->intermolecular_ilist)[i].iatoms.clear();
597 gmx::ArrayRef<const InteractionsOfType> interactions = intermolecular_interactions->interactions;
599 if (!interactions[i].interactionTypes.empty())
601 flags = interaction_function[i].flags;
602 /* For intermolecular interactions we (currently)
603 * only support potentials.
604 * Constraints and virtual sites would be possible,
605 * but require a lot of extra (bug-prone) code.
607 if (!(flags & IF_BOND))
610 "The intermolecular_interaction section may only contain bonded "
613 else if (NRAL(i) == 1) /* e.g. position restraints */
616 "Single atom interactions don't make sense in the "
617 "intermolecular_interaction section, you can put them in the "
618 "moleculetype section");
620 else if (flags & IF_CHEMBOND)
623 "The intermolecular_interaction can not contain chemically bonding "
628 enter_function(&(interactions[i]),
629 static_cast<t_functype>(i),
633 &(*mtop->intermolecular_ilist)[i],
637 mtop->bIntermolecularInteractions = TRUE;
642 if (!mtop->bIntermolecularInteractions)
644 mtop->intermolecular_ilist.reset(nullptr);
648 ffp->fudgeQQ = fudgeQQ;