3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35 /* This file is completely threadsafe - keep it that way! */
46 #include "gmx_fatal.h"
51 #include "gpp_atomtype.h"
53 static int round_check(real r,int limit,int ftype,const char *name)
62 if (r-i > 0.01 || r-i < -0.01)
63 gmx_fatal(FARGS,"A non-integer value (%f) was supplied for '%s' in %s",
64 r,name,interaction_function[ftype].longname);
67 gmx_fatal(FARGS,"Value of '%s' in %s is %d, which is smaller than the minimum of %d",
68 name,interaction_function[ftype].longname,i,limit);
73 static void set_ljparams(int comb,double reppow,real v,real w,
76 if (comb == eCOMB_ARITHMETIC || comb == eCOMB_GEOM_SIG_EPS) {
79 *c12 = 4*w*pow(v,reppow);
81 /* Interpret negative sigma as c6=0 and c12 with -sigma */
83 *c12 = 4*w*pow(-v,reppow);
91 static void assign_param(t_functype ftype,t_iparams *newparam,
92 real old[MAXFORCEPARAM],int comb,double reppow)
98 for(j=0; (j<MAXFORCEPARAM); j++)
100 newparam->generic.buf[j]=0.0;
104 /* Post processing of input data: store cosine iso angle itself */
105 newparam->harmonic.rA =cos(old[0]*DEG2RAD);
106 newparam->harmonic.krA=old[1];
107 newparam->harmonic.rB =cos(old[2]*DEG2RAD);
108 newparam->harmonic.krB=old[3];
111 /* Post processing of input data: store square of length itself */
112 newparam->harmonic.rA =sqr(old[0]);
113 newparam->harmonic.krA=old[1];
114 newparam->harmonic.rB =sqr(old[2]);
115 newparam->harmonic.krB=old[3];
118 newparam->fene.bm=old[0];
119 newparam->fene.kb=old[1];
122 newparam->restraint.lowA = old[0];
123 newparam->restraint.up1A = old[1];
124 newparam->restraint.up2A = old[2];
125 newparam->restraint.kA = old[3];
126 newparam->restraint.lowB = old[4];
127 newparam->restraint.up1B = old[5];
128 newparam->restraint.up2B = old[6];
129 newparam->restraint.kB = old[7];
135 newparam->tab.table = round_check(old[0],0,ftype,"table index");
136 newparam->tab.kA = old[1];
137 newparam->tab.kB = old[3];
139 case F_CROSS_BOND_BONDS:
140 newparam->cross_bb.r1e=old[0];
141 newparam->cross_bb.r2e=old[1];
142 newparam->cross_bb.krr=old[2];
144 case F_CROSS_BOND_ANGLES:
145 newparam->cross_ba.r1e=old[0];
146 newparam->cross_ba.r2e=old[1];
147 newparam->cross_ba.r3e=old[2];
148 newparam->cross_ba.krt=old[3];
151 newparam->u_b.theta=old[0];
152 newparam->u_b.ktheta=old[1];
153 newparam->u_b.r13=old[2];
154 newparam->u_b.kUB=old[3];
156 case F_QUARTIC_ANGLES:
157 newparam->qangle.theta=old[0];
159 newparam->qangle.c[i]=old[i+1];
165 newparam->harmonic.rA =old[0];
166 newparam->harmonic.krA=old[1];
167 newparam->harmonic.rB =old[2];
168 newparam->harmonic.krB=old[3];
171 newparam->morse.b0 =old[0];
172 newparam->morse.cb =old[1];
173 newparam->morse.beta =old[2];
176 newparam->cubic.b0 =old[0];
177 newparam->cubic.kb =old[1];
178 newparam->cubic.kcub =old[2];
183 newparam->polarize.alpha = old[0];
186 newparam->wpol.al_x =old[0];
187 newparam->wpol.al_y =old[1];
188 newparam->wpol.al_z =old[2];
189 newparam->wpol.rOH =old[3];
190 newparam->wpol.rHH =old[4];
191 newparam->wpol.rOD =old[5];
194 newparam->thole.a = old[0];
195 newparam->thole.alpha1 = old[1];
196 newparam->thole.alpha2 = old[2];
197 if ((old[1] > 0) && (old[2] > 0))
198 newparam->thole.rfac = old[0]*pow(old[1]*old[2],-1.0/6.0);
200 newparam->thole.rfac = 1;
203 newparam->bham.a = old[0];
204 newparam->bham.b = old[1];
205 newparam->bham.c = old[2];
208 set_ljparams(comb,reppow,old[0],old[1],&newparam->lj14.c6A,&newparam->lj14.c12A);
209 set_ljparams(comb,reppow,old[2],old[3],&newparam->lj14.c6B,&newparam->lj14.c12B);
212 newparam->ljc14.fqq = old[0];
213 newparam->ljc14.qi = old[1];
214 newparam->ljc14.qj = old[2];
215 set_ljparams(comb,reppow,old[3],old[4],&newparam->ljc14.c6,&newparam->ljc14.c12);
218 newparam->ljcnb.qi = old[0];
219 newparam->ljcnb.qj = old[1];
220 set_ljparams(comb,reppow,old[2],old[3],&newparam->ljcnb.c6,&newparam->ljcnb.c12);
223 set_ljparams(comb,reppow,old[0],old[1],&newparam->lj.c6,&newparam->lj.c12);
229 newparam->pdihs.phiA = old[0];
230 newparam->pdihs.cpA = old[1];
232 /* Dont do any checks if all parameters are zero (such interactions will be removed).
233 * Change 20100720: Amber occasionally uses negative multiplicities (mathematically OK),
234 * so I have changed the lower limit to -99 /EL
236 * Second, if the force constant is zero in both A and B states, we set the phase
237 * and multiplicity to zero too so the interaction gets removed during clean-up.
239 newparam->pdihs.phiB = old[3];
240 newparam->pdihs.cpB = old[4];
242 if( fabs(newparam->pdihs.cpA) < GMX_REAL_MIN && fabs(newparam->pdihs.cpB) < GMX_REAL_MIN )
244 newparam->pdihs.phiA = 0.0;
245 newparam->pdihs.phiB = 0.0;
246 newparam->pdihs.mult = 0;
250 newparam->pdihs.mult = round_check(old[2],-99,ftype,"multiplicity");
255 newparam->posres.fcA[XX] = old[0];
256 newparam->posres.fcA[YY] = old[1];
257 newparam->posres.fcA[ZZ] = old[2];
258 newparam->posres.fcB[XX] = old[3];
259 newparam->posres.fcB[YY] = old[4];
260 newparam->posres.fcB[ZZ] = old[5];
261 newparam->posres.pos0A[XX] = old[6];
262 newparam->posres.pos0A[YY] = old[7];
263 newparam->posres.pos0A[ZZ] = old[8];
264 newparam->posres.pos0B[XX] = old[9];
265 newparam->posres.pos0B[YY] = old[10];
266 newparam->posres.pos0B[ZZ] = old[11];
269 newparam->disres.label = round_check(old[0],0,ftype,"label");
270 newparam->disres.type = round_check(old[1],1,ftype,"type'");
271 newparam->disres.low = old[2];
272 newparam->disres.up1 = old[3];
273 newparam->disres.up2 = old[4];
274 newparam->disres.kfac = old[5];
277 newparam->orires.ex = round_check(old[0],1,ftype,"experiment") - 1;
278 newparam->orires.label = round_check(old[1],1,ftype,"label");
279 newparam->orires.power = round_check(old[2],0,ftype,"power");
280 newparam->orires.c = old[3];
281 newparam->orires.obs = old[4];
282 newparam->orires.kfac = old[5];
285 newparam->dihres.label = round_check(old[0],0,ftype,"label");
286 newparam->dihres.phi = old[1];
287 newparam->dihres.dphi = old[2];
288 newparam->dihres.kfac = old[3];
289 newparam->dihres.power = round_check(old[4],0,ftype,"power");
292 for (i=0; (i<NR_RBDIHS); i++) {
293 newparam->rbdihs.rbcA[i]=old[i];
294 newparam->rbdihs.rbcB[i]=old[NR_RBDIHS+i];
298 /* Read the dihedral parameters to temporary arrays,
299 * and convert them to the computationally faster
300 * Ryckaert-Bellemans form.
302 /* Use conversion formula for OPLS to Ryckaert-Bellemans: */
303 newparam->rbdihs.rbcA[0]=old[1]+0.5*(old[0]+old[2]);
304 newparam->rbdihs.rbcA[1]=0.5*(3.0*old[2]-old[0]);
305 newparam->rbdihs.rbcA[2]=4.0*old[3]-old[1];
306 newparam->rbdihs.rbcA[3]=-2.0*old[2];
307 newparam->rbdihs.rbcA[4]=-4.0*old[3];
308 newparam->rbdihs.rbcA[5]=0.0;
310 newparam->rbdihs.rbcB[0]=old[NR_FOURDIHS+1]+0.5*(old[NR_FOURDIHS+0]+old[NR_FOURDIHS+2]);
311 newparam->rbdihs.rbcB[1]=0.5*(3.0*old[NR_FOURDIHS+2]-old[NR_FOURDIHS+0]);
312 newparam->rbdihs.rbcB[2]=4.0*old[NR_FOURDIHS+3]-old[NR_FOURDIHS+1];
313 newparam->rbdihs.rbcB[3]=-2.0*old[NR_FOURDIHS+2];
314 newparam->rbdihs.rbcB[4]=-4.0*old[NR_FOURDIHS+3];
315 newparam->rbdihs.rbcB[5]=0.0;
319 newparam->constr.dA = old[0];
320 newparam->constr.dB = old[1];
323 newparam->settle.doh=old[0];
324 newparam->settle.dhh=old[1];
332 newparam->vsite.a=old[0];
333 newparam->vsite.b=old[1];
334 newparam->vsite.c=old[2];
335 newparam->vsite.d=old[3];
336 newparam->vsite.e=old[4];
337 newparam->vsite.f=old[5];
340 newparam->vsite.a=old[1] * cos(DEG2RAD * old[0]);
341 newparam->vsite.b=old[1] * sin(DEG2RAD * old[0]);
342 newparam->vsite.c=old[2];
343 newparam->vsite.d=old[3];
344 newparam->vsite.e=old[4];
345 newparam->vsite.f=old[5];
348 newparam->vsiten.n = round_check(old[0],1,ftype,"number of atoms");
349 newparam->vsiten.a = old[1];
352 newparam->cmap.cmapA=old[0];
353 newparam->cmap.cmapB=old[1];
358 newparam->gb.sar = old[0];
359 newparam->gb.st = old[1];
360 newparam->gb.pi = old[2];
361 newparam->gb.gbr = old[3];
362 newparam->gb.bmlt = old[4];
365 gmx_fatal(FARGS,"unknown function type %d in %s line %d",
366 ftype,__FILE__,__LINE__);
370 static int enter_params(gmx_ffparams_t *ffparams, t_functype ftype,
371 real forceparams[MAXFORCEPARAM],int comb,real reppow,
372 int start,gmx_bool bAppend)
377 assign_param(ftype,&newparam,forceparams,comb,reppow);
379 for (type=start; (type<ffparams->ntypes); type++) {
380 if (ffparams->functype[type]==ftype) {
381 if (memcmp(&newparam,&ffparams->iparams[type],(size_t)sizeof(newparam)) == 0)
387 type = ffparams->ntypes;
390 fprintf(debug,"copying newparam to ffparams->iparams[%d] (ntypes=%d)\n",
391 type,ffparams->ntypes);
392 memcpy(&ffparams->iparams[type],&newparam,(size_t)sizeof(newparam));
395 ffparams->functype[type]=ftype;
400 static void append_interaction(t_ilist *ilist,
401 int type,int nral,atom_id a[MAXATOMLIST])
408 ilist->iatoms[where1++]=type;
409 for (i=0; (i<nral); i++)
410 ilist->iatoms[where1++]=a[i];
413 static void enter_function(t_params *p,t_functype ftype,int comb,real reppow,
414 gmx_ffparams_t *ffparams,t_ilist *il,
416 gmx_bool bNB,gmx_bool bAppend)
418 int k,type,nr,nral,delta,start;
420 start = ffparams->ntypes;
423 for (k=0; k<nr; k++) {
424 if (*maxtypes <= ffparams->ntypes) {
426 srenew(ffparams->functype,*maxtypes);
427 srenew(ffparams->iparams, *maxtypes);
429 fprintf(debug,"%s, line %d: srenewed idef->functype and idef->iparams to %d\n",
430 __FILE__,__LINE__,*maxtypes);
432 type = enter_params(ffparams,ftype,p->param[k].c,comb,reppow,start,bAppend);
436 srenew(il->iatoms,il->nr+delta);
437 append_interaction(il,type,nral,p->param[k].a);
442 static void new_interaction_list(t_ilist *ilist)
450 void convert_params(int atnr,t_params nbtypes[],
451 t_molinfo *mi,int comb,double reppow,real fudgeQQ,
462 ffp = &mtop->ffparams;
465 ffp->functype = NULL;
467 ffp->reppow = reppow;
469 enter_function(&(nbtypes[F_LJ]), (t_functype)F_LJ, comb,reppow,ffp,NULL,
470 &maxtypes,TRUE,TRUE);
471 enter_function(&(nbtypes[F_BHAM]),(t_functype)F_BHAM, comb,reppow,ffp,NULL,
472 &maxtypes,TRUE,TRUE);
474 for(mt=0; mt<mtop->nmoltype; mt++) {
475 molt = &mtop->moltype[mt];
476 for(i=0; (i<F_NRE); i++) {
477 molt->ilist[i].nr = 0;
478 molt->ilist[i].iatoms = NULL;
480 plist = mi[mt].plist;
482 flags = interaction_function[i].flags;
483 if ((i != F_LJ) && (i != F_BHAM) && ((flags & IF_BOND) ||
484 (flags & IF_VSITE) ||
485 (flags & IF_CONSTRAINT))) {
486 enter_function(&(plist[i]),(t_functype)i,comb,reppow,
488 &maxtypes,FALSE,(i == F_POSRES));
493 fprintf(debug,"%s, line %d: There are %d functypes in idef\n",
494 __FILE__,__LINE__,ffp->ntypes);
497 ffp->fudgeQQ = fudgeQQ;