Don't use POSIX fnmatch() for pattern matching.
[gromacs/qmmm-gamess-us.git] / src / kernel / convparm.c
blob6472ba49a4324501650c866a97cea0c8376bb5b3
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++)
99 newparam->generic.buf[j]=0.0;
101 switch (ftype) {
102 case F_G96ANGLES:
103 /* Post processing of input data: store cosine iso angle itself */
104 newparam->harmonic.rA =cos(old[0]*DEG2RAD);
105 newparam->harmonic.krA=old[1];
106 newparam->harmonic.rB =cos(old[2]*DEG2RAD);
107 newparam->harmonic.krB=old[3];
108 break;
109 case F_G96BONDS:
110 /* Post processing of input data: store square of length itself */
111 newparam->harmonic.rA =sqr(old[0]);
112 newparam->harmonic.krA=old[1];
113 newparam->harmonic.rB =sqr(old[2]);
114 newparam->harmonic.krB=old[3];
115 break;
116 case F_FENEBONDS:
117 newparam->fene.bm=old[0];
118 newparam->fene.kb=old[1];
119 break;
120 case F_TABBONDS:
121 case F_TABBONDSNC:
122 case F_TABANGLES:
123 case F_TABDIHS:
124 newparam->tab.table = round_check(old[0],0,ftype,"table index");
125 newparam->tab.kA = old[1];
126 newparam->tab.kB = old[3];
127 break;
128 case F_CROSS_BOND_BONDS:
129 newparam->cross_bb.r1e=old[0];
130 newparam->cross_bb.r2e=old[1];
131 newparam->cross_bb.krr=old[2];
132 break;
133 case F_CROSS_BOND_ANGLES:
134 newparam->cross_ba.r1e=old[0];
135 newparam->cross_ba.r2e=old[1];
136 newparam->cross_ba.r3e=old[2];
137 newparam->cross_ba.krt=old[3];
138 break;
139 case F_UREY_BRADLEY:
140 newparam->u_b.theta=old[0];
141 newparam->u_b.ktheta=old[1];
142 newparam->u_b.r13=old[2];
143 newparam->u_b.kUB=old[3];
144 break;
145 case F_QUARTIC_ANGLES:
146 newparam->qangle.theta=old[0];
147 for(i=0; i<5; i++)
148 newparam->qangle.c[i]=old[i+1];
149 break;
150 case F_ANGLES:
151 case F_BONDS:
152 case F_HARMONIC:
153 case F_IDIHS:
154 newparam->harmonic.rA =old[0];
155 newparam->harmonic.krA=old[1];
156 newparam->harmonic.rB =old[2];
157 newparam->harmonic.krB=old[3];
158 break;
159 case F_MORSE:
160 newparam->morse.b0 =old[0];
161 newparam->morse.cb =old[1];
162 newparam->morse.beta =old[2];
163 break;
164 case F_CUBICBONDS:
165 newparam->cubic.b0 =old[0];
166 newparam->cubic.kb =old[1];
167 newparam->cubic.kcub =old[2];
168 break;
169 case F_CONNBONDS:
170 break;
171 case F_POLARIZATION:
172 newparam->polarize.alpha = old[0];
173 break;
174 case F_WATER_POL:
175 newparam->wpol.al_x =old[0];
176 newparam->wpol.al_y =old[1];
177 newparam->wpol.al_z =old[2];
178 newparam->wpol.rOH =old[3];
179 newparam->wpol.rHH =old[4];
180 newparam->wpol.rOD =old[5];
181 break;
182 case F_THOLE_POL:
183 newparam->thole.a = old[0];
184 newparam->thole.alpha1 = old[1];
185 newparam->thole.alpha2 = old[2];
186 if ((old[1] > 0) && (old[2] > 0))
187 newparam->thole.rfac = old[0]*pow(old[1]*old[2],-1.0/6.0);
188 else
189 newparam->thole.rfac = 1;
190 break;
191 case F_BHAM:
192 newparam->bham.a = old[0];
193 newparam->bham.b = old[1];
194 newparam->bham.c = old[2];
195 break;
196 case F_LJ14:
197 set_ljparams(comb,reppow,old[0],old[1],&newparam->lj14.c6A,&newparam->lj14.c12A);
198 set_ljparams(comb,reppow,old[2],old[3],&newparam->lj14.c6B,&newparam->lj14.c12B);
199 break;
200 case F_LJC14_Q:
201 newparam->ljc14.fqq = old[0];
202 newparam->ljc14.qi = old[1];
203 newparam->ljc14.qj = old[2];
204 set_ljparams(comb,reppow,old[3],old[4],&newparam->ljc14.c6,&newparam->ljc14.c12);
205 break;
206 case F_LJC_PAIRS_NB:
207 newparam->ljcnb.qi = old[0];
208 newparam->ljcnb.qj = old[1];
209 set_ljparams(comb,reppow,old[2],old[3],&newparam->ljcnb.c6,&newparam->ljcnb.c12);
210 break;
211 case F_LJ:
212 set_ljparams(comb,reppow,old[0],old[1],&newparam->lj.c6,&newparam->lj.c12);
213 break;
214 case F_PDIHS:
215 case F_ANGRES:
216 case F_ANGRESZ:
217 newparam->pdihs.phiA = old[0];
218 newparam->pdihs.cpA = old[1];
220 /* Dont do any checks if all parameters are zero (such interactions will be removed) */
221 tmp=fabs(old[0])+fabs(old[1])+fabs(old[2])+fabs(old[3])+fabs(old[4]);
222 newparam->pdihs.mult = (tmp < GMX_REAL_MIN) ? 0 : round_check(old[2],1,ftype,"multiplicity");
224 newparam->pdihs.phiB = old[3];
225 newparam->pdihs.cpB = old[4];
226 break;
227 case F_POSRES:
228 newparam->posres.fcA[XX] = old[0];
229 newparam->posres.fcA[YY] = old[1];
230 newparam->posres.fcA[ZZ] = old[2];
231 newparam->posres.fcB[XX] = old[3];
232 newparam->posres.fcB[YY] = old[4];
233 newparam->posres.fcB[ZZ] = old[5];
234 newparam->posres.pos0A[XX] = old[6];
235 newparam->posres.pos0A[YY] = old[7];
236 newparam->posres.pos0A[ZZ] = old[8];
237 newparam->posres.pos0B[XX] = old[9];
238 newparam->posres.pos0B[YY] = old[10];
239 newparam->posres.pos0B[ZZ] = old[11];
240 break;
241 case F_DISRES:
242 newparam->disres.label = round_check(old[0],0,ftype,"label");
243 newparam->disres.type = round_check(old[1],1,ftype,"type'");
244 newparam->disres.low = old[2];
245 newparam->disres.up1 = old[3];
246 newparam->disres.up2 = old[4];
247 newparam->disres.kfac = old[5];
248 break;
249 case F_ORIRES:
250 newparam->orires.ex = round_check(old[0],1,ftype,"experiment") - 1;
251 newparam->orires.label = round_check(old[1],1,ftype,"label");
252 newparam->orires.power = round_check(old[2],0,ftype,"power");
253 newparam->orires.c = old[3];
254 newparam->orires.obs = old[4];
255 newparam->orires.kfac = old[5];
256 break;
257 case F_DIHRES:
258 newparam->dihres.label = round_check(old[0],0,ftype,"label");
259 newparam->dihres.phi = old[1];
260 newparam->dihres.dphi = old[2];
261 newparam->dihres.kfac = old[3];
262 newparam->dihres.power = round_check(old[4],0,ftype,"power");
263 break;
264 case F_RBDIHS:
265 for (i=0; (i<NR_RBDIHS); i++) {
266 newparam->rbdihs.rbcA[i]=old[i];
267 newparam->rbdihs.rbcB[i]=old[NR_RBDIHS+i];
269 break;
270 case F_FOURDIHS:
271 /* Read the dihedral parameters to temporary arrays,
272 * and convert them to the computationally faster
273 * Ryckaert-Bellemans form.
275 /* Use conversion formula for OPLS to Ryckaert-Bellemans: */
276 newparam->rbdihs.rbcA[0]=old[1]+0.5*(old[0]+old[2]);
277 newparam->rbdihs.rbcA[1]=0.5*(3.0*old[2]-old[0]);
278 newparam->rbdihs.rbcA[2]=4.0*old[3]-old[1];
279 newparam->rbdihs.rbcA[3]=-2.0*old[2];
280 newparam->rbdihs.rbcA[4]=-4.0*old[3];
281 newparam->rbdihs.rbcA[5]=0.0;
283 newparam->rbdihs.rbcB[0]=old[NR_FOURDIHS+1]+0.5*(old[NR_FOURDIHS+0]+old[NR_FOURDIHS+2]);
284 newparam->rbdihs.rbcB[1]=0.5*(3.0*old[NR_FOURDIHS+2]-old[NR_FOURDIHS+0]);
285 newparam->rbdihs.rbcB[2]=4.0*old[NR_FOURDIHS+3]-old[NR_FOURDIHS+1];
286 newparam->rbdihs.rbcB[3]=-2.0*old[NR_FOURDIHS+2];
287 newparam->rbdihs.rbcB[4]=-4.0*old[NR_FOURDIHS+3];
288 newparam->rbdihs.rbcB[5]=0.0;
289 break;
290 case F_CONSTR:
291 case F_CONSTRNC:
292 newparam->constr.dA = old[0];
293 newparam->constr.dB = old[1];
294 break;
295 case F_SETTLE:
296 newparam->settle.doh=old[0];
297 newparam->settle.dhh=old[1];
298 break;
299 case F_VSITE2:
300 case F_VSITE3:
301 case F_VSITE3FD:
302 case F_VSITE3OUT:
303 case F_VSITE4FD:
304 case F_VSITE4FDN:
305 newparam->vsite.a=old[0];
306 newparam->vsite.b=old[1];
307 newparam->vsite.c=old[2];
308 newparam->vsite.d=old[3];
309 newparam->vsite.e=old[4];
310 newparam->vsite.f=old[5];
311 break;
312 case F_VSITE3FAD:
313 newparam->vsite.a=old[1] * cos(DEG2RAD * old[0]);
314 newparam->vsite.b=old[1] * sin(DEG2RAD * old[0]);
315 newparam->vsite.c=old[2];
316 newparam->vsite.d=old[3];
317 newparam->vsite.e=old[4];
318 newparam->vsite.f=old[5];
319 break;
320 case F_VSITEN:
321 newparam->vsiten.n = round_check(old[0],1,ftype,"number of atoms");
322 newparam->vsiten.a = old[1];
323 break;
324 case F_CMAP:
325 newparam->cmap.cmapA=old[0];
326 newparam->cmap.cmapB=old[1];
327 break;
328 case F_GB12:
329 case F_GB13:
330 case F_GB14:
331 newparam->gb.sar = old[0];
332 newparam->gb.st = old[1];
333 newparam->gb.pi = old[2];
334 newparam->gb.gbr = old[3];
335 newparam->gb.bmlt = old[4];
336 break;
337 default:
338 gmx_fatal(FARGS,"unknown function type %d in %s line %d",
339 ftype,__FILE__,__LINE__);
343 static int enter_params(gmx_ffparams_t *ffparams, t_functype ftype,
344 real forceparams[MAXFORCEPARAM],int comb,real reppow,
345 int start,bool bAppend)
347 t_iparams newparam;
348 int type;
350 assign_param(ftype,&newparam,forceparams,comb,reppow);
351 if (!bAppend) {
352 for (type=start; (type<ffparams->ntypes); type++) {
353 if (ffparams->functype[type]==ftype) {
354 if (memcmp(&newparam,&ffparams->iparams[type],(size_t)sizeof(newparam)) == 0)
355 return type;
359 else {
360 type = ffparams->ntypes;
362 if (debug)
363 fprintf(debug,"copying newparam to ffparams->iparams[%d] (ntypes=%d)\n",
364 type,ffparams->ntypes);
365 memcpy(&ffparams->iparams[type],&newparam,(size_t)sizeof(newparam));
367 ffparams->ntypes++;
368 ffparams->functype[type]=ftype;
370 return type;
373 static void append_interaction(t_ilist *ilist,
374 int type,int nral,atom_id a[MAXATOMLIST])
376 int i,where1;
378 where1 = ilist->nr;
379 ilist->nr += nral+1;
381 ilist->iatoms[where1++]=type;
382 for (i=0; (i<nral); i++)
383 ilist->iatoms[where1++]=a[i];
386 static void enter_function(t_params *p,t_functype ftype,int comb,real reppow,
387 gmx_ffparams_t *ffparams,t_ilist *il,
388 int *maxtypes,
389 bool bNB,bool bAppend)
391 int k,type,nr,nral,delta,start;
393 start = ffparams->ntypes;
394 nr = p->nr;
396 for (k=0; k<nr; k++) {
397 if (*maxtypes <= ffparams->ntypes) {
398 *maxtypes += 1000;
399 srenew(ffparams->functype,*maxtypes);
400 srenew(ffparams->iparams, *maxtypes);
401 if (debug)
402 fprintf(debug,"%s, line %d: srenewed idef->functype and idef->iparams to %d\n",
403 __FILE__,__LINE__,*maxtypes);
405 type = enter_params(ffparams,ftype,p->param[k].c,comb,reppow,start,bAppend);
406 if (!bNB) {
407 nral = NRAL(ftype);
408 delta = nr*(nral+1);
409 srenew(il->iatoms,il->nr+delta);
410 append_interaction(il,type,nral,p->param[k].a);
415 static void new_interaction_list(t_ilist *ilist)
417 int i;
419 ilist->nr=0;
420 ilist->iatoms=NULL;
423 void convert_params(int atnr,t_params nbtypes[],
424 t_molinfo *mi,int comb,double reppow,real fudgeQQ,
425 gmx_mtop_t *mtop)
427 int i,j,maxtypes,mt;
428 unsigned long flags;
429 gmx_ffparams_t *ffp;
430 gmx_moltype_t *molt;
431 t_params *plist;
433 maxtypes=0;
435 ffp = &mtop->ffparams;
436 ffp->ntypes = 0;
437 ffp->atnr = atnr;
438 ffp->functype = NULL;
439 ffp->iparams = NULL;
440 ffp->reppow = reppow;
442 enter_function(&(nbtypes[F_LJ]), (t_functype)F_LJ, comb,reppow,ffp,NULL,
443 &maxtypes,TRUE,TRUE);
444 enter_function(&(nbtypes[F_BHAM]),(t_functype)F_BHAM, comb,reppow,ffp,NULL,
445 &maxtypes,TRUE,TRUE);
447 for(mt=0; mt<mtop->nmoltype; mt++) {
448 molt = &mtop->moltype[mt];
449 for(i=0; (i<F_NRE); i++) {
450 molt->ilist[i].nr = 0;
451 molt->ilist[i].iatoms = NULL;
453 plist = mi[mt].plist;
455 flags = interaction_function[i].flags;
456 if ((i != F_LJ) && (i != F_BHAM) && ((flags & IF_BOND) ||
457 (flags & IF_VSITE) ||
458 (flags & IF_CONSTRAINT))) {
459 enter_function(&(plist[i]),(t_functype)i,comb,reppow,
460 ffp,&molt->ilist[i],
461 &maxtypes,FALSE,(i == F_POSRES));
465 if (debug) {
466 fprintf(debug,"%s, line %d: There are %d functypes in idef\n",
467 __FILE__,__LINE__,ffp->ntypes);
470 ffp->fudgeQQ = fudgeQQ;