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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * 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 "gmx_fatal.h"
54 #include "gpp_atomtype.h"
57 static int round_check(real r,int limit,int ftype,const char *name)
66 if (r-i > 0.01 || r-i < -0.01)
67 gmx_fatal(FARGS,"A non-integer value (%f) was supplied for '%s' in %s",
68 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);
77 static void set_ljparams(int comb,double reppow,real v,real w,
80 if (comb == eCOMB_ARITHMETIC || comb == eCOMB_GEOM_SIG_EPS) {
83 *c12 = 4*w*pow(v,reppow);
85 /* Interpret negative sigma as c6=0 and c12 with -sigma */
87 *c12 = 4*w*pow(-v,reppow);
95 static void assign_param(t_functype ftype,t_iparams *newparam,
96 real old[MAXFORCEPARAM],int comb,double reppow)
102 for(j=0; (j<MAXFORCEPARAM); j++)
104 newparam->generic.buf[j]=0.0;
108 /* Post processing of input data: store cosine iso angle itself */
109 newparam->harmonic.rA =cos(old[0]*DEG2RAD);
110 newparam->harmonic.krA=old[1];
111 newparam->harmonic.rB =cos(old[2]*DEG2RAD);
112 newparam->harmonic.krB=old[3];
115 /* Post processing of input data: store square of length itself */
116 newparam->harmonic.rA =sqr(old[0]);
117 newparam->harmonic.krA=old[1];
118 newparam->harmonic.rB =sqr(old[2]);
119 newparam->harmonic.krB=old[3];
122 newparam->fene.bm=old[0];
123 newparam->fene.kb=old[1];
126 newparam->restraint.lowA = old[0];
127 newparam->restraint.up1A = old[1];
128 newparam->restraint.up2A = old[2];
129 newparam->restraint.kA = old[3];
130 newparam->restraint.lowB = old[4];
131 newparam->restraint.up1B = old[5];
132 newparam->restraint.up2B = old[6];
133 newparam->restraint.kB = old[7];
139 newparam->tab.table = round_check(old[0],0,ftype,"table index");
140 newparam->tab.kA = old[1];
141 newparam->tab.kB = old[3];
143 case F_CROSS_BOND_BONDS:
144 newparam->cross_bb.r1e=old[0];
145 newparam->cross_bb.r2e=old[1];
146 newparam->cross_bb.krr=old[2];
148 case F_CROSS_BOND_ANGLES:
149 newparam->cross_ba.r1e=old[0];
150 newparam->cross_ba.r2e=old[1];
151 newparam->cross_ba.r3e=old[2];
152 newparam->cross_ba.krt=old[3];
155 newparam->u_b.thetaA=old[0];
156 newparam->u_b.kthetaA=old[1];
157 newparam->u_b.r13A=old[2];
158 newparam->u_b.kUBA=old[3];
159 newparam->u_b.thetaB=old[4];
160 newparam->u_b.kthetaB=old[5];
161 newparam->u_b.r13B=old[6];
162 newparam->u_b.kUBB=old[7];
164 case F_QUARTIC_ANGLES:
165 newparam->qangle.theta=old[0];
167 newparam->qangle.c[i]=old[i+1];
169 case F_LINEAR_ANGLES:
170 newparam->linangle.aA = old[0];
171 newparam->linangle.klinA = old[1];
172 newparam->linangle.aB = old[2];
173 newparam->linangle.klinB = old[3];
179 newparam->harmonic.rA =old[0];
180 newparam->harmonic.krA=old[1];
181 newparam->harmonic.rB =old[2];
182 newparam->harmonic.krB=old[3];
185 newparam->morse.b0A =old[0];
186 newparam->morse.cbA =old[1];
187 newparam->morse.betaA =old[2];
188 newparam->morse.b0B =old[3];
189 newparam->morse.cbB =old[4];
190 newparam->morse.betaB =old[5];
193 newparam->cubic.b0 =old[0];
194 newparam->cubic.kb =old[1];
195 newparam->cubic.kcub =old[2];
200 newparam->polarize.alpha = old[0];
203 newparam->anharm_polarize.alpha = old[0];
204 newparam->anharm_polarize.drcut = old[1];
205 newparam->anharm_polarize.khyp = old[2];
208 newparam->wpol.al_x =old[0];
209 newparam->wpol.al_y =old[1];
210 newparam->wpol.al_z =old[2];
211 newparam->wpol.rOH =old[3];
212 newparam->wpol.rHH =old[4];
213 newparam->wpol.rOD =old[5];
216 newparam->thole.a = old[0];
217 newparam->thole.alpha1 = old[1];
218 newparam->thole.alpha2 = old[2];
219 if ((old[1] > 0) && (old[2] > 0))
220 newparam->thole.rfac = old[0]*pow(old[1]*old[2],-1.0/6.0);
222 newparam->thole.rfac = 1;
225 newparam->bham.a = old[0];
226 newparam->bham.b = old[1];
227 newparam->bham.c = old[2];
230 set_ljparams(comb,reppow,old[0],old[1],&newparam->lj14.c6A,&newparam->lj14.c12A);
231 set_ljparams(comb,reppow,old[2],old[3],&newparam->lj14.c6B,&newparam->lj14.c12B);
234 newparam->ljc14.fqq = old[0];
235 newparam->ljc14.qi = old[1];
236 newparam->ljc14.qj = old[2];
237 set_ljparams(comb,reppow,old[3],old[4],&newparam->ljc14.c6,&newparam->ljc14.c12);
240 newparam->ljcnb.qi = old[0];
241 newparam->ljcnb.qj = old[1];
242 set_ljparams(comb,reppow,old[2],old[3],&newparam->ljcnb.c6,&newparam->ljcnb.c12);
245 set_ljparams(comb,reppow,old[0],old[1],&newparam->lj.c6,&newparam->lj.c12);
251 newparam->pdihs.phiA = old[0];
252 newparam->pdihs.cpA = old[1];
254 /* Dont do any checks if all parameters are zero (such interactions will be removed).
255 * Change 20100720: Amber occasionally uses negative multiplicities (mathematically OK),
256 * so I have changed the lower limit to -99 /EL
258 * Second, if the force constant is zero in both A and B states, we set the phase
259 * and multiplicity to zero too so the interaction gets removed during clean-up.
261 newparam->pdihs.phiB = old[3];
262 newparam->pdihs.cpB = old[4];
264 if( fabs(newparam->pdihs.cpA) < GMX_REAL_MIN && fabs(newparam->pdihs.cpB) < GMX_REAL_MIN )
266 newparam->pdihs.phiA = 0.0;
267 newparam->pdihs.phiB = 0.0;
268 newparam->pdihs.mult = 0;
272 newparam->pdihs.mult = round_check(old[2],-99,ftype,"multiplicity");
277 newparam->posres.fcA[XX] = old[0];
278 newparam->posres.fcA[YY] = old[1];
279 newparam->posres.fcA[ZZ] = old[2];
280 newparam->posres.fcB[XX] = old[3];
281 newparam->posres.fcB[YY] = old[4];
282 newparam->posres.fcB[ZZ] = old[5];
283 newparam->posres.pos0A[XX] = old[6];
284 newparam->posres.pos0A[YY] = old[7];
285 newparam->posres.pos0A[ZZ] = old[8];
286 newparam->posres.pos0B[XX] = old[9];
287 newparam->posres.pos0B[YY] = old[10];
288 newparam->posres.pos0B[ZZ] = old[11];
291 newparam->disres.label = round_check(old[0],0,ftype,"label");
292 newparam->disres.type = round_check(old[1],1,ftype,"type'");
293 newparam->disres.low = old[2];
294 newparam->disres.up1 = old[3];
295 newparam->disres.up2 = old[4];
296 newparam->disres.kfac = old[5];
299 newparam->orires.ex = round_check(old[0],1,ftype,"experiment") - 1;
300 newparam->orires.label = round_check(old[1],1,ftype,"label");
301 newparam->orires.power = round_check(old[2],0,ftype,"power");
302 newparam->orires.c = old[3];
303 newparam->orires.obs = old[4];
304 newparam->orires.kfac = old[5];
307 newparam->dihres.phiA = old[0];
308 newparam->dihres.dphiA = old[1];
309 newparam->dihres.kfacA = old[2];
310 newparam->dihres.phiB = old[3];
311 newparam->dihres.dphiB = old[4];
312 newparam->dihres.kfacB = old[5];
315 for (i=0; (i<NR_RBDIHS); i++) {
316 newparam->rbdihs.rbcA[i]=old[i];
317 newparam->rbdihs.rbcB[i]=old[NR_RBDIHS+i];
321 /* Read the dihedral parameters to temporary arrays,
322 * and convert them to the computationally faster
323 * Ryckaert-Bellemans form.
325 /* Use conversion formula for OPLS to Ryckaert-Bellemans: */
326 newparam->rbdihs.rbcA[0]=old[1]+0.5*(old[0]+old[2]);
327 newparam->rbdihs.rbcA[1]=0.5*(3.0*old[2]-old[0]);
328 newparam->rbdihs.rbcA[2]=4.0*old[3]-old[1];
329 newparam->rbdihs.rbcA[3]=-2.0*old[2];
330 newparam->rbdihs.rbcA[4]=-4.0*old[3];
331 newparam->rbdihs.rbcA[5]=0.0;
333 newparam->rbdihs.rbcB[0]=old[NR_FOURDIHS+1]+0.5*(old[NR_FOURDIHS+0]+old[NR_FOURDIHS+2]);
334 newparam->rbdihs.rbcB[1]=0.5*(3.0*old[NR_FOURDIHS+2]-old[NR_FOURDIHS+0]);
335 newparam->rbdihs.rbcB[2]=4.0*old[NR_FOURDIHS+3]-old[NR_FOURDIHS+1];
336 newparam->rbdihs.rbcB[3]=-2.0*old[NR_FOURDIHS+2];
337 newparam->rbdihs.rbcB[4]=-4.0*old[NR_FOURDIHS+3];
338 newparam->rbdihs.rbcB[5]=0.0;
342 newparam->constr.dA = old[0];
343 newparam->constr.dB = old[1];
346 newparam->settle.doh=old[0];
347 newparam->settle.dhh=old[1];
355 newparam->vsite.a=old[0];
356 newparam->vsite.b=old[1];
357 newparam->vsite.c=old[2];
358 newparam->vsite.d=old[3];
359 newparam->vsite.e=old[4];
360 newparam->vsite.f=old[5];
363 newparam->vsite.a=old[1] * cos(DEG2RAD * old[0]);
364 newparam->vsite.b=old[1] * sin(DEG2RAD * old[0]);
365 newparam->vsite.c=old[2];
366 newparam->vsite.d=old[3];
367 newparam->vsite.e=old[4];
368 newparam->vsite.f=old[5];
371 newparam->vsiten.n = round_check(old[0],1,ftype,"number of atoms");
372 newparam->vsiten.a = old[1];
375 newparam->cmap.cmapA=old[0];
376 newparam->cmap.cmapB=old[1];
381 newparam->gb.sar = old[0];
382 newparam->gb.st = old[1];
383 newparam->gb.pi = old[2];
384 newparam->gb.gbr = old[3];
385 newparam->gb.bmlt = old[4];
388 gmx_fatal(FARGS,"unknown function type %d in %s line %d",
389 ftype,__FILE__,__LINE__);
393 static int enter_params(gmx_ffparams_t *ffparams, t_functype ftype,
394 real forceparams[MAXFORCEPARAM],int comb,real reppow,
395 int start,gmx_bool bAppend)
400 assign_param(ftype,&newparam,forceparams,comb,reppow);
402 for (type=start; (type<ffparams->ntypes); type++) {
403 if (ffparams->functype[type]==ftype) {
404 if (F_GB13 == ftype) {
405 /* Occasionally, the way the 1-3 reference distance is
406 * computed can lead to non-binary-identical results, but I
408 if ((gmx_within_tol(newparam.gb.sar, ffparams->iparams[type].gb.sar, 1e-6)) &&
409 (gmx_within_tol(newparam.gb.st, ffparams->iparams[type].gb.st, 1e-6)) &&
410 (gmx_within_tol(newparam.gb.pi, ffparams->iparams[type].gb.pi, 1e-6)) &&
411 (gmx_within_tol(newparam.gb.gbr, ffparams->iparams[type].gb.gbr, 1e-6)) &&
412 (gmx_within_tol(newparam.gb.bmlt, ffparams->iparams[type].gb.bmlt, 1e-6))) {
417 if (memcmp(&newparam,&ffparams->iparams[type],(size_t)sizeof(newparam)) == 0)
424 type = ffparams->ntypes;
427 fprintf(debug,"copying newparam to ffparams->iparams[%d] (ntypes=%d)\n",
428 type,ffparams->ntypes);
429 memcpy(&ffparams->iparams[type],&newparam,(size_t)sizeof(newparam));
432 ffparams->functype[type]=ftype;
437 static void append_interaction(t_ilist *ilist,
438 int type,int nral,atom_id a[MAXATOMLIST])
445 ilist->iatoms[where1++]=type;
446 for (i=0; (i<nral); i++)
447 ilist->iatoms[where1++]=a[i];
450 static void enter_function(t_params *p,t_functype ftype,int comb,real reppow,
451 gmx_ffparams_t *ffparams,t_ilist *il,
453 gmx_bool bNB,gmx_bool bAppend)
455 int k,type,nr,nral,delta,start;
457 start = ffparams->ntypes;
460 for (k=0; k<nr; k++) {
461 if (*maxtypes <= ffparams->ntypes) {
463 srenew(ffparams->functype,*maxtypes);
464 srenew(ffparams->iparams, *maxtypes);
466 fprintf(debug,"%s, line %d: srenewed idef->functype and idef->iparams to %d\n",
467 __FILE__,__LINE__,*maxtypes);
469 type = enter_params(ffparams,ftype,p->param[k].c,comb,reppow,start,bAppend);
473 srenew(il->iatoms,il->nr+delta);
474 append_interaction(il,type,nral,p->param[k].a);
479 void convert_params(int atnr,t_params nbtypes[],
480 t_molinfo *mi,int comb,double reppow,real fudgeQQ,
491 ffp = &mtop->ffparams;
494 ffp->functype = NULL;
496 ffp->reppow = reppow;
498 enter_function(&(nbtypes[F_LJ]), (t_functype)F_LJ, comb,reppow,ffp,NULL,
499 &maxtypes,TRUE,TRUE);
500 enter_function(&(nbtypes[F_BHAM]),(t_functype)F_BHAM, comb,reppow,ffp,NULL,
501 &maxtypes,TRUE,TRUE);
503 for(mt=0; mt<mtop->nmoltype; mt++) {
504 molt = &mtop->moltype[mt];
505 for(i=0; (i<F_NRE); i++) {
506 molt->ilist[i].nr = 0;
507 molt->ilist[i].iatoms = NULL;
509 plist = mi[mt].plist;
511 flags = interaction_function[i].flags;
512 if ((i != F_LJ) && (i != F_BHAM) && ((flags & IF_BOND) ||
513 (flags & IF_VSITE) ||
514 (flags & IF_CONSTRAINT))) {
515 enter_function(&(plist[i]),(t_functype)i,comb,reppow,
517 &maxtypes,FALSE,(i == F_POSRES));
522 fprintf(debug,"%s, line %d: There are %d functypes in idef\n",
523 __FILE__,__LINE__,ffp->ntypes);
526 ffp->fudgeQQ = fudgeQQ;