Merge branch 'master' of git@git.gromacs.org:gromacs
[gromacs/rigid-bodies.git] / src / kernel / convparm.c
blob0af4ceac029c7f3cda1516acfe1fe41506352190
1 /*
2 *
3 * This source code is part of
4 *
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
9 * VERSION 3.2.0
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
32 * And Hey:
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35 /* This file is completely threadsafe - keep it that way! */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 #include <math.h>
41 #include "sysstuff.h"
42 #include "physics.h"
43 #include "vec.h"
44 #include "smalloc.h"
45 #include "typedefs.h"
46 #include "gmx_fatal.h"
47 #include "topio.h"
48 #include "toputil.h"
49 #include "convparm.h"
50 #include "names.h"
51 #include "gpp_atomtype.h"
53 static int round_check(real r,int limit,int ftype,const char *name)
55 int i;
57 if (r >= 0)
58 i = (int)(r + 0.5);
59 else
60 i = (int)(r - 0.5);
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);
66 if (i < limit)
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);
70 return i;
73 static void set_ljparams(int comb,double reppow,real v,real w,
74 real *c6,real *c12)
76 if (comb == eCOMB_ARITHMETIC || comb == eCOMB_GEOM_SIG_EPS) {
77 if (v >= 0) {
78 *c6 = 4*w*pow(v,6.0);
79 *c12 = 4*w*pow(v,reppow);
80 } else {
81 /* Interpret negative sigma as c6=0 and c12 with -sigma */
82 *c6 = 0;
83 *c12 = 4*w*pow(-v,reppow);
85 } else {
86 *c6 = v;
87 *c12 = w;
91 static void assign_param(t_functype ftype,t_iparams *newparam,
92 real old[MAXFORCEPARAM],int comb,double reppow)
94 int i,j;
95 real tmp;
97 /* Set to zero */
98 for(j=0; (j<MAXFORCEPARAM); j++)
100 newparam->generic.buf[j]=0.0;
102 switch (ftype) {
103 case F_G96ANGLES:
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];
109 break;
110 case F_G96BONDS:
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];
116 break;
117 case F_FENEBONDS:
118 newparam->fene.bm=old[0];
119 newparam->fene.kb=old[1];
120 break;
121 case F_RESTRBONDS:
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];
130 break;
131 case F_TABBONDS:
132 case F_TABBONDSNC:
133 case F_TABANGLES:
134 case F_TABDIHS:
135 newparam->tab.table = round_check(old[0],0,ftype,"table index");
136 newparam->tab.kA = old[1];
137 newparam->tab.kB = old[3];
138 break;
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];
143 break;
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];
149 break;
150 case F_UREY_BRADLEY:
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];
155 break;
156 case F_QUARTIC_ANGLES:
157 newparam->qangle.theta=old[0];
158 for(i=0; i<5; i++)
159 newparam->qangle.c[i]=old[i+1];
160 break;
161 case F_ANGLES:
162 case F_BONDS:
163 case F_HARMONIC:
164 case F_IDIHS:
165 newparam->harmonic.rA =old[0];
166 newparam->harmonic.krA=old[1];
167 newparam->harmonic.rB =old[2];
168 newparam->harmonic.krB=old[3];
169 break;
170 case F_MORSE:
171 newparam->morse.b0 =old[0];
172 newparam->morse.cb =old[1];
173 newparam->morse.beta =old[2];
174 break;
175 case F_CUBICBONDS:
176 newparam->cubic.b0 =old[0];
177 newparam->cubic.kb =old[1];
178 newparam->cubic.kcub =old[2];
179 break;
180 case F_CONNBONDS:
181 break;
182 case F_POLARIZATION:
183 newparam->polarize.alpha = old[0];
184 break;
185 case F_WATER_POL:
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];
192 break;
193 case F_THOLE_POL:
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);
199 else
200 newparam->thole.rfac = 1;
201 break;
202 case F_BHAM:
203 newparam->bham.a = old[0];
204 newparam->bham.b = old[1];
205 newparam->bham.c = old[2];
206 break;
207 case F_LJ14:
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);
210 break;
211 case F_LJC14_Q:
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);
216 break;
217 case F_LJC_PAIRS_NB:
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);
221 break;
222 case F_LJ:
223 set_ljparams(comb,reppow,old[0],old[1],&newparam->lj.c6,&newparam->lj.c12);
224 break;
225 case F_PDIHS:
226 case F_ANGRES:
227 case F_ANGRESZ:
228 newparam->pdihs.phiA = old[0];
229 newparam->pdihs.cpA = old[1];
231 /* Dont do any checks if all parameters are zero (such interactions will be removed) */
232 tmp=fabs(old[0])+fabs(old[1])+fabs(old[2])+fabs(old[3])+fabs(old[4]);
233 newparam->pdihs.mult = (tmp < GMX_REAL_MIN) ? 0 : round_check(old[2],1,ftype,"multiplicity");
235 newparam->pdihs.phiB = old[3];
236 newparam->pdihs.cpB = old[4];
237 break;
238 case F_POSRES:
239 newparam->posres.fcA[XX] = old[0];
240 newparam->posres.fcA[YY] = old[1];
241 newparam->posres.fcA[ZZ] = old[2];
242 newparam->posres.fcB[XX] = old[3];
243 newparam->posres.fcB[YY] = old[4];
244 newparam->posres.fcB[ZZ] = old[5];
245 newparam->posres.pos0A[XX] = old[6];
246 newparam->posres.pos0A[YY] = old[7];
247 newparam->posres.pos0A[ZZ] = old[8];
248 newparam->posres.pos0B[XX] = old[9];
249 newparam->posres.pos0B[YY] = old[10];
250 newparam->posres.pos0B[ZZ] = old[11];
251 break;
252 case F_DISRES:
253 newparam->disres.label = round_check(old[0],0,ftype,"label");
254 newparam->disres.type = round_check(old[1],1,ftype,"type'");
255 newparam->disres.low = old[2];
256 newparam->disres.up1 = old[3];
257 newparam->disres.up2 = old[4];
258 newparam->disres.kfac = old[5];
259 break;
260 case F_ORIRES:
261 newparam->orires.ex = round_check(old[0],1,ftype,"experiment") - 1;
262 newparam->orires.label = round_check(old[1],1,ftype,"label");
263 newparam->orires.power = round_check(old[2],0,ftype,"power");
264 newparam->orires.c = old[3];
265 newparam->orires.obs = old[4];
266 newparam->orires.kfac = old[5];
267 break;
268 case F_DIHRES:
269 newparam->dihres.label = round_check(old[0],0,ftype,"label");
270 newparam->dihres.phi = old[1];
271 newparam->dihres.dphi = old[2];
272 newparam->dihres.kfac = old[3];
273 newparam->dihres.power = round_check(old[4],0,ftype,"power");
274 break;
275 case F_RBDIHS:
276 for (i=0; (i<NR_RBDIHS); i++) {
277 newparam->rbdihs.rbcA[i]=old[i];
278 newparam->rbdihs.rbcB[i]=old[NR_RBDIHS+i];
280 break;
281 case F_FOURDIHS:
282 /* Read the dihedral parameters to temporary arrays,
283 * and convert them to the computationally faster
284 * Ryckaert-Bellemans form.
286 /* Use conversion formula for OPLS to Ryckaert-Bellemans: */
287 newparam->rbdihs.rbcA[0]=old[1]+0.5*(old[0]+old[2]);
288 newparam->rbdihs.rbcA[1]=0.5*(3.0*old[2]-old[0]);
289 newparam->rbdihs.rbcA[2]=4.0*old[3]-old[1];
290 newparam->rbdihs.rbcA[3]=-2.0*old[2];
291 newparam->rbdihs.rbcA[4]=-4.0*old[3];
292 newparam->rbdihs.rbcA[5]=0.0;
294 newparam->rbdihs.rbcB[0]=old[NR_FOURDIHS+1]+0.5*(old[NR_FOURDIHS+0]+old[NR_FOURDIHS+2]);
295 newparam->rbdihs.rbcB[1]=0.5*(3.0*old[NR_FOURDIHS+2]-old[NR_FOURDIHS+0]);
296 newparam->rbdihs.rbcB[2]=4.0*old[NR_FOURDIHS+3]-old[NR_FOURDIHS+1];
297 newparam->rbdihs.rbcB[3]=-2.0*old[NR_FOURDIHS+2];
298 newparam->rbdihs.rbcB[4]=-4.0*old[NR_FOURDIHS+3];
299 newparam->rbdihs.rbcB[5]=0.0;
300 break;
301 case F_CONSTR:
302 case F_CONSTRNC:
303 newparam->constr.dA = old[0];
304 newparam->constr.dB = old[1];
305 break;
306 case F_SETTLE:
307 newparam->settle.doh=old[0];
308 newparam->settle.dhh=old[1];
309 break;
310 case F_VSITE2:
311 case F_VSITE3:
312 case F_VSITE3FD:
313 case F_VSITE3OUT:
314 case F_VSITE4FD:
315 case F_VSITE4FDN:
316 newparam->vsite.a=old[0];
317 newparam->vsite.b=old[1];
318 newparam->vsite.c=old[2];
319 newparam->vsite.d=old[3];
320 newparam->vsite.e=old[4];
321 newparam->vsite.f=old[5];
322 break;
323 case F_VSITE3FAD:
324 newparam->vsite.a=old[1] * cos(DEG2RAD * old[0]);
325 newparam->vsite.b=old[1] * sin(DEG2RAD * old[0]);
326 newparam->vsite.c=old[2];
327 newparam->vsite.d=old[3];
328 newparam->vsite.e=old[4];
329 newparam->vsite.f=old[5];
330 break;
331 case F_VSITEN:
332 newparam->vsiten.n = round_check(old[0],1,ftype,"number of atoms");
333 newparam->vsiten.a = old[1];
334 break;
335 case F_CMAP:
336 newparam->cmap.cmapA=old[0];
337 newparam->cmap.cmapB=old[1];
338 break;
339 case F_GB12:
340 case F_GB13:
341 case F_GB14:
342 newparam->gb.sar = old[0];
343 newparam->gb.st = old[1];
344 newparam->gb.pi = old[2];
345 newparam->gb.gbr = old[3];
346 newparam->gb.bmlt = old[4];
347 break;
348 default:
349 gmx_fatal(FARGS,"unknown function type %d in %s line %d",
350 ftype,__FILE__,__LINE__);
354 static int enter_params(gmx_ffparams_t *ffparams, t_functype ftype,
355 real forceparams[MAXFORCEPARAM],int comb,real reppow,
356 int start,bool bAppend)
358 t_iparams newparam;
359 int type;
361 assign_param(ftype,&newparam,forceparams,comb,reppow);
362 if (!bAppend) {
363 for (type=start; (type<ffparams->ntypes); type++) {
364 if (ffparams->functype[type]==ftype) {
365 if (memcmp(&newparam,&ffparams->iparams[type],(size_t)sizeof(newparam)) == 0)
366 return type;
370 else {
371 type = ffparams->ntypes;
373 if (debug)
374 fprintf(debug,"copying newparam to ffparams->iparams[%d] (ntypes=%d)\n",
375 type,ffparams->ntypes);
376 memcpy(&ffparams->iparams[type],&newparam,(size_t)sizeof(newparam));
378 ffparams->ntypes++;
379 ffparams->functype[type]=ftype;
381 return type;
384 static void append_interaction(t_ilist *ilist,
385 int type,int nral,atom_id a[MAXATOMLIST])
387 int i,where1;
389 where1 = ilist->nr;
390 ilist->nr += nral+1;
392 ilist->iatoms[where1++]=type;
393 for (i=0; (i<nral); i++)
394 ilist->iatoms[where1++]=a[i];
397 static void enter_function(t_params *p,t_functype ftype,int comb,real reppow,
398 gmx_ffparams_t *ffparams,t_ilist *il,
399 int *maxtypes,
400 bool bNB,bool bAppend)
402 int k,type,nr,nral,delta,start;
404 start = ffparams->ntypes;
405 nr = p->nr;
407 for (k=0; k<nr; k++) {
408 if (*maxtypes <= ffparams->ntypes) {
409 *maxtypes += 1000;
410 srenew(ffparams->functype,*maxtypes);
411 srenew(ffparams->iparams, *maxtypes);
412 if (debug)
413 fprintf(debug,"%s, line %d: srenewed idef->functype and idef->iparams to %d\n",
414 __FILE__,__LINE__,*maxtypes);
416 type = enter_params(ffparams,ftype,p->param[k].c,comb,reppow,start,bAppend);
417 if (!bNB) {
418 nral = NRAL(ftype);
419 delta = nr*(nral+1);
420 srenew(il->iatoms,il->nr+delta);
421 append_interaction(il,type,nral,p->param[k].a);
426 static void new_interaction_list(t_ilist *ilist)
428 int i;
430 ilist->nr=0;
431 ilist->iatoms=NULL;
434 void convert_params(int atnr,t_params nbtypes[],
435 t_molinfo *mi,int comb,double reppow,real fudgeQQ,
436 gmx_mtop_t *mtop)
438 int i,j,maxtypes,mt;
439 unsigned long flags;
440 gmx_ffparams_t *ffp;
441 gmx_moltype_t *molt;
442 t_params *plist;
444 maxtypes=0;
446 ffp = &mtop->ffparams;
447 ffp->ntypes = 0;
448 ffp->atnr = atnr;
449 ffp->functype = NULL;
450 ffp->iparams = NULL;
451 ffp->reppow = reppow;
453 enter_function(&(nbtypes[F_LJ]), (t_functype)F_LJ, comb,reppow,ffp,NULL,
454 &maxtypes,TRUE,TRUE);
455 enter_function(&(nbtypes[F_BHAM]),(t_functype)F_BHAM, comb,reppow,ffp,NULL,
456 &maxtypes,TRUE,TRUE);
458 for(mt=0; mt<mtop->nmoltype; mt++) {
459 molt = &mtop->moltype[mt];
460 for(i=0; (i<F_NRE); i++) {
461 molt->ilist[i].nr = 0;
462 molt->ilist[i].iatoms = NULL;
464 plist = mi[mt].plist;
466 flags = interaction_function[i].flags;
467 if ((i != F_LJ) && (i != F_BHAM) && ((flags & IF_BOND) ||
468 (flags & IF_VSITE) ||
469 (flags & IF_CONSTRAINT))) {
470 enter_function(&(plist[i]),(t_functype)i,comb,reppow,
471 ffp,&molt->ilist[i],
472 &maxtypes,FALSE,(i == F_POSRES));
476 if (debug) {
477 fprintf(debug,"%s, line %d: There are %d functypes in idef\n",
478 __FILE__,__LINE__,ffp->ntypes);
481 ffp->fudgeQQ = fudgeQQ;