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, 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/utility/smalloc.h"
50 #include "gmx_fatal.h"
55 #include "gpp_atomtype.h"
56 #include "gromacs/math/utilities.h"
58 static int round_check(real r, int limit, int ftype, const char *name)
71 if (r-i > 0.01 || r-i < -0.01)
73 gmx_fatal(FARGS, "A non-integer value (%f) was supplied for '%s' in %s",
74 r, name, interaction_function[ftype].longname);
79 gmx_fatal(FARGS, "Value of '%s' in %s is %d, which is smaller than the minimum of %d",
80 name, interaction_function[ftype].longname, i, limit);
86 static void set_ljparams(int comb, double reppow, real v, real w,
89 if (comb == eCOMB_ARITHMETIC || comb == eCOMB_GEOM_SIG_EPS)
93 *c6 = 4*w*pow(v, 6.0);
94 *c12 = 4*w*pow(v, reppow);
98 /* Interpret negative sigma as c6=0 and c12 with -sigma */
100 *c12 = 4*w*pow(-v, reppow);
110 /* A return value of 0 means parameters were assigned successfully,
111 * returning -1 means this is an all-zero interaction that should not be added.
114 assign_param(t_functype ftype, t_iparams *newparam,
115 real old[MAXFORCEPARAM], int comb, double reppow)
119 gmx_bool all_param_zero = TRUE;
122 for (j = 0; (j < MAXFORCEPARAM); j++)
124 newparam->generic.buf[j] = 0.0;
125 /* If all parameters are zero we might not add some interaction types (selected below).
126 * We cannot apply this to ALL interactions, since many have valid reasons for having
127 * zero parameters (e.g. an index to a Cmap interaction, or LJ parameters), but
128 * we use it for angles and torsions that are typically generated automatically.
130 all_param_zero = (all_param_zero == TRUE) && fabs(old[j]) < GMX_REAL_MIN;
133 if (all_param_zero == TRUE)
135 if (IS_ANGLE(ftype) || IS_RESTRAINT_TYPE(ftype) || ftype == F_IDIHS ||
136 ftype == F_PDIHS || ftype == F_PIDIHS || ftype == F_RBDIHS || ftype == F_FOURDIHS)
145 /* Post processing of input data: store cosine iso angle itself */
146 newparam->harmonic.rA = cos(old[0]*DEG2RAD);
147 newparam->harmonic.krA = old[1];
148 newparam->harmonic.rB = cos(old[2]*DEG2RAD);
149 newparam->harmonic.krB = old[3];
152 /* Post processing of input data: store square of length itself */
153 newparam->harmonic.rA = sqr(old[0]);
154 newparam->harmonic.krA = old[1];
155 newparam->harmonic.rB = sqr(old[2]);
156 newparam->harmonic.krB = old[3];
159 newparam->fene.bm = old[0];
160 newparam->fene.kb = old[1];
163 newparam->restraint.lowA = old[0];
164 newparam->restraint.up1A = old[1];
165 newparam->restraint.up2A = old[2];
166 newparam->restraint.kA = old[3];
167 newparam->restraint.lowB = old[4];
168 newparam->restraint.up1B = old[5];
169 newparam->restraint.up2B = old[6];
170 newparam->restraint.kB = old[7];
176 newparam->tab.table = round_check(old[0], 0, ftype, "table index");
177 newparam->tab.kA = old[1];
178 newparam->tab.kB = old[3];
180 case F_CROSS_BOND_BONDS:
181 newparam->cross_bb.r1e = old[0];
182 newparam->cross_bb.r2e = old[1];
183 newparam->cross_bb.krr = old[2];
185 case F_CROSS_BOND_ANGLES:
186 newparam->cross_ba.r1e = old[0];
187 newparam->cross_ba.r2e = old[1];
188 newparam->cross_ba.r3e = old[2];
189 newparam->cross_ba.krt = old[3];
192 newparam->u_b.thetaA = old[0];
193 newparam->u_b.kthetaA = old[1];
194 newparam->u_b.r13A = old[2];
195 newparam->u_b.kUBA = old[3];
196 newparam->u_b.thetaB = old[4];
197 newparam->u_b.kthetaB = old[5];
198 newparam->u_b.r13B = old[6];
199 newparam->u_b.kUBB = old[7];
201 case F_QUARTIC_ANGLES:
202 newparam->qangle.theta = old[0];
203 for (i = 0; i < 5; i++)
205 newparam->qangle.c[i] = old[i+1];
208 case F_LINEAR_ANGLES:
209 newparam->linangle.aA = old[0];
210 newparam->linangle.klinA = old[1];
211 newparam->linangle.aB = old[2];
212 newparam->linangle.klinB = old[3];
218 newparam->harmonic.rA = old[0];
219 newparam->harmonic.krA = old[1];
220 newparam->harmonic.rB = old[2];
221 newparam->harmonic.krB = old[3];
224 newparam->harmonic.rA = old[0];
225 newparam->harmonic.krA = old[1];
228 newparam->morse.b0A = old[0];
229 newparam->morse.cbA = old[1];
230 newparam->morse.betaA = old[2];
231 newparam->morse.b0B = old[3];
232 newparam->morse.cbB = old[4];
233 newparam->morse.betaB = old[5];
236 newparam->cubic.b0 = old[0];
237 newparam->cubic.kb = old[1];
238 newparam->cubic.kcub = old[2];
243 newparam->polarize.alpha = old[0];
246 newparam->anharm_polarize.alpha = old[0];
247 newparam->anharm_polarize.drcut = old[1];
248 newparam->anharm_polarize.khyp = old[2];
251 newparam->wpol.al_x = old[0];
252 newparam->wpol.al_y = old[1];
253 newparam->wpol.al_z = old[2];
254 newparam->wpol.rOH = old[3];
255 newparam->wpol.rHH = old[4];
256 newparam->wpol.rOD = old[5];
259 newparam->thole.a = old[0];
260 newparam->thole.alpha1 = old[1];
261 newparam->thole.alpha2 = old[2];
262 if ((old[1] > 0) && (old[2] > 0))
264 newparam->thole.rfac = old[0]*pow(old[1]*old[2], -1.0/6.0);
268 newparam->thole.rfac = 1;
272 newparam->bham.a = old[0];
273 newparam->bham.b = old[1];
274 newparam->bham.c = old[2];
277 set_ljparams(comb, reppow, old[0], old[1], &newparam->lj14.c6A, &newparam->lj14.c12A);
278 set_ljparams(comb, reppow, old[2], old[3], &newparam->lj14.c6B, &newparam->lj14.c12B);
281 newparam->ljc14.fqq = old[0];
282 newparam->ljc14.qi = old[1];
283 newparam->ljc14.qj = old[2];
284 set_ljparams(comb, reppow, old[3], old[4], &newparam->ljc14.c6, &newparam->ljc14.c12);
287 newparam->ljcnb.qi = old[0];
288 newparam->ljcnb.qj = old[1];
289 set_ljparams(comb, reppow, old[2], old[3], &newparam->ljcnb.c6, &newparam->ljcnb.c12);
292 set_ljparams(comb, reppow, old[0], old[1], &newparam->lj.c6, &newparam->lj.c12);
298 newparam->pdihs.phiA = old[0];
299 newparam->pdihs.cpA = old[1];
301 /* Change 20100720: Amber occasionally uses negative multiplicities (mathematically OK),
302 * so I have changed the lower limit to -99 /EL
304 newparam->pdihs.phiB = old[3];
305 newparam->pdihs.cpB = old[4];
306 /* If both force constants are zero there is no interaction. Return -1 to signal
307 * this entry should NOT be added.
309 if (fabs(newparam->pdihs.cpA) < GMX_REAL_MIN && fabs(newparam->pdihs.cpB) < GMX_REAL_MIN)
314 newparam->pdihs.mult = round_check(old[2], -99, ftype, "multiplicity");
318 newparam->pdihs.phiA = old[0];
319 newparam->pdihs.cpA = old[1];
322 newparam->posres.fcA[XX] = old[0];
323 newparam->posres.fcA[YY] = old[1];
324 newparam->posres.fcA[ZZ] = old[2];
325 newparam->posres.fcB[XX] = old[3];
326 newparam->posres.fcB[YY] = old[4];
327 newparam->posres.fcB[ZZ] = old[5];
328 newparam->posres.pos0A[XX] = old[6];
329 newparam->posres.pos0A[YY] = old[7];
330 newparam->posres.pos0A[ZZ] = old[8];
331 newparam->posres.pos0B[XX] = old[9];
332 newparam->posres.pos0B[YY] = old[10];
333 newparam->posres.pos0B[ZZ] = old[11];
336 newparam->fbposres.geom = round_check(old[0], 0, ftype, "geometry");
337 if (!(newparam->fbposres.geom > efbposresZERO && newparam->fbposres.geom < efbposresNR))
339 gmx_fatal(FARGS, "Invalid geometry for flat-bottomed position restraint.\n"
340 "Expected number between 1 and %d. Found %d\n", efbposresNR-1,
341 newparam->fbposres.geom);
343 newparam->fbposres.r = old[1];
344 newparam->fbposres.k = old[2];
345 newparam->fbposres.pos0[XX] = old[3];
346 newparam->fbposres.pos0[YY] = old[4];
347 newparam->fbposres.pos0[ZZ] = old[5];
350 newparam->disres.label = round_check(old[0], 0, ftype, "label");
351 newparam->disres.type = round_check(old[1], 1, ftype, "type'");
352 newparam->disres.low = old[2];
353 newparam->disres.up1 = old[3];
354 newparam->disres.up2 = old[4];
355 newparam->disres.kfac = old[5];
358 newparam->orires.ex = round_check(old[0], 1, ftype, "experiment") - 1;
359 newparam->orires.label = round_check(old[1], 1, ftype, "label");
360 newparam->orires.power = round_check(old[2], 0, ftype, "power");
361 newparam->orires.c = old[3];
362 newparam->orires.obs = old[4];
363 newparam->orires.kfac = old[5];
366 newparam->dihres.phiA = old[0];
367 newparam->dihres.dphiA = old[1];
368 newparam->dihres.kfacA = old[2];
369 newparam->dihres.phiB = old[3];
370 newparam->dihres.dphiB = old[4];
371 newparam->dihres.kfacB = old[5];
374 for (i = 0; (i < NR_RBDIHS); i++)
376 newparam->rbdihs.rbcA[i] = old[i];
377 newparam->rbdihs.rbcB[i] = old[NR_RBDIHS+i];
381 for (i = 0; (i < NR_CBTDIHS); i++)
383 newparam->cbtdihs.cbtcA[i] = old[i];
387 /* Read the dihedral parameters to temporary arrays,
388 * and convert them to the computationally faster
389 * Ryckaert-Bellemans form.
391 /* Use conversion formula for OPLS to Ryckaert-Bellemans: */
392 newparam->rbdihs.rbcA[0] = old[1]+0.5*(old[0]+old[2]);
393 newparam->rbdihs.rbcA[1] = 0.5*(3.0*old[2]-old[0]);
394 newparam->rbdihs.rbcA[2] = 4.0*old[3]-old[1];
395 newparam->rbdihs.rbcA[3] = -2.0*old[2];
396 newparam->rbdihs.rbcA[4] = -4.0*old[3];
397 newparam->rbdihs.rbcA[5] = 0.0;
399 newparam->rbdihs.rbcB[0] = old[NR_FOURDIHS+1]+0.5*(old[NR_FOURDIHS+0]+old[NR_FOURDIHS+2]);
400 newparam->rbdihs.rbcB[1] = 0.5*(3.0*old[NR_FOURDIHS+2]-old[NR_FOURDIHS+0]);
401 newparam->rbdihs.rbcB[2] = 4.0*old[NR_FOURDIHS+3]-old[NR_FOURDIHS+1];
402 newparam->rbdihs.rbcB[3] = -2.0*old[NR_FOURDIHS+2];
403 newparam->rbdihs.rbcB[4] = -4.0*old[NR_FOURDIHS+3];
404 newparam->rbdihs.rbcB[5] = 0.0;
408 newparam->constr.dA = old[0];
409 newparam->constr.dB = old[1];
412 newparam->settle.doh = old[0];
413 newparam->settle.dhh = old[1];
421 newparam->vsite.a = old[0];
422 newparam->vsite.b = old[1];
423 newparam->vsite.c = old[2];
424 newparam->vsite.d = old[3];
425 newparam->vsite.e = old[4];
426 newparam->vsite.f = old[5];
429 newparam->vsite.a = old[1] * cos(DEG2RAD * old[0]);
430 newparam->vsite.b = old[1] * sin(DEG2RAD * old[0]);
431 newparam->vsite.c = old[2];
432 newparam->vsite.d = old[3];
433 newparam->vsite.e = old[4];
434 newparam->vsite.f = old[5];
437 newparam->vsiten.n = round_check(old[0], 1, ftype, "number of atoms");
438 newparam->vsiten.a = old[1];
441 newparam->cmap.cmapA = old[0];
442 newparam->cmap.cmapB = old[1];
447 newparam->gb.sar = old[0];
448 newparam->gb.st = old[1];
449 newparam->gb.pi = old[2];
450 newparam->gb.gbr = old[3];
451 newparam->gb.bmlt = old[4];
454 gmx_fatal(FARGS, "unknown function type %d in %s line %d",
455 ftype, __FILE__, __LINE__);
460 static int enter_params(gmx_ffparams_t *ffparams, t_functype ftype,
461 real forceparams[MAXFORCEPARAM], int comb, real reppow,
462 int start, gmx_bool bAppend)
468 if ( (rc = assign_param(ftype, &newparam, forceparams, comb, reppow)) < 0)
470 /* -1 means this interaction is all-zero and should not be added */
476 for (type = start; (type < ffparams->ntypes); type++)
478 if (ffparams->functype[type] == ftype)
482 /* Occasionally, the way the 1-3 reference distance is
483 * computed can lead to non-binary-identical results, but I
485 if ((gmx_within_tol(newparam.gb.sar, ffparams->iparams[type].gb.sar, 1e-6)) &&
486 (gmx_within_tol(newparam.gb.st, ffparams->iparams[type].gb.st, 1e-6)) &&
487 (gmx_within_tol(newparam.gb.pi, ffparams->iparams[type].gb.pi, 1e-6)) &&
488 (gmx_within_tol(newparam.gb.gbr, ffparams->iparams[type].gb.gbr, 1e-6)) &&
489 (gmx_within_tol(newparam.gb.bmlt, ffparams->iparams[type].gb.bmlt, 1e-6)))
496 if (memcmp(&newparam, &ffparams->iparams[type], (size_t)sizeof(newparam)) == 0)
506 type = ffparams->ntypes;
510 fprintf(debug, "copying newparam to ffparams->iparams[%d] (ntypes=%d)\n",
511 type, ffparams->ntypes);
513 memcpy(&ffparams->iparams[type], &newparam, (size_t)sizeof(newparam));
516 ffparams->functype[type] = ftype;
521 static void append_interaction(t_ilist *ilist,
522 int type, int nral, atom_id a[MAXATOMLIST])
529 ilist->iatoms[where1++] = type;
530 for (i = 0; (i < nral); i++)
532 ilist->iatoms[where1++] = a[i];
536 static void enter_function(t_params *p, t_functype ftype, int comb, real reppow,
537 gmx_ffparams_t *ffparams, t_ilist *il,
539 gmx_bool bNB, gmx_bool bAppend)
541 int k, type, nr, nral, delta, start;
543 start = ffparams->ntypes;
546 for (k = 0; k < nr; k++)
548 if (*maxtypes <= ffparams->ntypes)
551 srenew(ffparams->functype, *maxtypes);
552 srenew(ffparams->iparams, *maxtypes);
555 fprintf(debug, "%s, line %d: srenewed idef->functype and idef->iparams to %d\n",
556 __FILE__, __LINE__, *maxtypes);
559 type = enter_params(ffparams, ftype, p->param[k].c, comb, reppow, start, bAppend);
560 /* Type==-1 is used as a signal that this interaction is all-zero and should not be added. */
561 if (!bNB && type >= 0)
565 srenew(il->iatoms, il->nr+delta);
566 append_interaction(il, type, nral, p->param[k].a);
571 void convert_params(int atnr, t_params nbtypes[],
572 t_molinfo *mi, int comb, double reppow, real fudgeQQ,
575 int i, j, maxtypes, mt;
583 ffp = &mtop->ffparams;
586 ffp->functype = NULL;
588 ffp->reppow = reppow;
590 enter_function(&(nbtypes[F_LJ]), (t_functype)F_LJ, comb, reppow, ffp, NULL,
591 &maxtypes, TRUE, TRUE);
592 enter_function(&(nbtypes[F_BHAM]), (t_functype)F_BHAM, comb, reppow, ffp, NULL,
593 &maxtypes, TRUE, TRUE);
595 for (mt = 0; mt < mtop->nmoltype; mt++)
597 molt = &mtop->moltype[mt];
598 for (i = 0; (i < F_NRE); i++)
600 molt->ilist[i].nr = 0;
601 molt->ilist[i].iatoms = NULL;
603 plist = mi[mt].plist;
605 flags = interaction_function[i].flags;
606 if ((i != F_LJ) && (i != F_BHAM) && ((flags & IF_BOND) ||
607 (flags & IF_VSITE) ||
608 (flags & IF_CONSTRAINT)))
610 enter_function(&(plist[i]), (t_functype)i, comb, reppow,
611 ffp, &molt->ilist[i],
612 &maxtypes, FALSE, (i == F_POSRES || i == F_FBPOSRES));
618 fprintf(debug, "%s, line %d: There are %d functypes in idef\n",
619 __FILE__, __LINE__, ffp->ntypes);
622 ffp->fudgeQQ = fudgeQQ;