fdef7ea8155732c56a97cc607e1453081ac8feb1
[alexxy/gromacs.git] / src / kernel / convparm.c
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.
14
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.
19  * 
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.
26  * 
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.
29  * 
30  * For more info, check our website at http://www.gromacs.org
31  * 
32  * And Hey:
33  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
34  */
35 /* This file is completely threadsafe - keep it that way! */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
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"
52 #include "maths.h"
53
54 static int round_check(real r,int limit,int ftype,const char *name)
55 {
56   int i;
57
58   if (r >= 0)
59     i = (int)(r + 0.5);
60   else
61     i = (int)(r - 0.5);
62
63   if (r-i > 0.01 || r-i < -0.01)
64     gmx_fatal(FARGS,"A non-integer value (%f) was supplied for '%s' in %s",
65               r,name,interaction_function[ftype].longname);
66
67   if (i < limit)
68     gmx_fatal(FARGS,"Value of '%s' in %s is %d, which is smaller than the minimum of %d",
69               name,interaction_function[ftype].longname,i,limit);
70
71   return i;
72 }
73
74 static void set_ljparams(int comb,double reppow,real v,real w,
75                          real *c6,real *c12)
76 {
77   if (comb == eCOMB_ARITHMETIC || comb == eCOMB_GEOM_SIG_EPS) {
78     if (v >= 0) {
79       *c6  = 4*w*pow(v,6.0);
80       *c12 = 4*w*pow(v,reppow);
81     } else {
82       /* Interpret negative sigma as c6=0 and c12 with -sigma */
83       *c6  = 0;
84       *c12 = 4*w*pow(-v,reppow);
85     }
86   } else {
87     *c6  = v;
88     *c12 = w;
89   }
90 }
91
92 static void assign_param(t_functype ftype,t_iparams *newparam,
93                          real old[MAXFORCEPARAM],int comb,double reppow)
94 {
95   int  i,j;
96   real tmp;
97
98   /* Set to zero */
99   for(j=0; (j<MAXFORCEPARAM); j++) 
100     {
101       newparam->generic.buf[j]=0.0;
102     }
103   switch (ftype) {
104   case F_G96ANGLES:
105     /* Post processing of input data: store cosine iso angle itself */
106     newparam->harmonic.rA =cos(old[0]*DEG2RAD);
107     newparam->harmonic.krA=old[1];
108     newparam->harmonic.rB =cos(old[2]*DEG2RAD);
109     newparam->harmonic.krB=old[3];
110     break;
111   case F_G96BONDS:
112     /* Post processing of input data: store square of length itself */
113     newparam->harmonic.rA =sqr(old[0]);
114     newparam->harmonic.krA=old[1];
115     newparam->harmonic.rB =sqr(old[2]);
116     newparam->harmonic.krB=old[3];
117     break;
118   case F_FENEBONDS:
119     newparam->fene.bm=old[0];
120     newparam->fene.kb=old[1];
121     break;
122   case F_RESTRBONDS:
123     newparam->restraint.lowA = old[0];
124     newparam->restraint.up1A = old[1];
125     newparam->restraint.up2A = old[2];
126     newparam->restraint.kA   = old[3];
127     newparam->restraint.lowB = old[4];
128     newparam->restraint.up1B = old[5];
129     newparam->restraint.up2B = old[6];
130     newparam->restraint.kB   = old[7];
131     break;
132   case F_TABBONDS:
133   case F_TABBONDSNC:
134   case F_TABANGLES:
135   case F_TABDIHS:
136     newparam->tab.table = round_check(old[0],0,ftype,"table index");
137     newparam->tab.kA    = old[1];
138     newparam->tab.kB    = old[3];
139     break;
140   case F_CROSS_BOND_BONDS:
141     newparam->cross_bb.r1e=old[0];
142     newparam->cross_bb.r2e=old[1];
143     newparam->cross_bb.krr=old[2];
144     break;
145   case F_CROSS_BOND_ANGLES:
146     newparam->cross_ba.r1e=old[0];
147     newparam->cross_ba.r2e=old[1];
148     newparam->cross_ba.r3e=old[2];
149     newparam->cross_ba.krt=old[3];
150     break;
151   case F_UREY_BRADLEY:
152     newparam->u_b.theta=old[0];
153     newparam->u_b.ktheta=old[1];
154     newparam->u_b.r13=old[2];
155     newparam->u_b.kUB=old[3];
156     break;
157   case F_QUARTIC_ANGLES:
158     newparam->qangle.theta=old[0];
159     for(i=0; i<5; i++)
160       newparam->qangle.c[i]=old[i+1];
161     break;
162   case F_ANGLES:
163   case F_BONDS:
164   case F_HARMONIC:
165   case F_IDIHS:
166     newparam->harmonic.rA =old[0];
167     newparam->harmonic.krA=old[1];
168     newparam->harmonic.rB =old[2];
169     newparam->harmonic.krB=old[3];
170     break;
171   case F_MORSE:
172     newparam->morse.b0    =old[0];
173     newparam->morse.cb    =old[1];
174     newparam->morse.beta  =old[2];
175     break;
176   case F_CUBICBONDS:
177     newparam->cubic.b0    =old[0];
178     newparam->cubic.kb    =old[1];
179     newparam->cubic.kcub  =old[2];
180     break;
181   case F_CONNBONDS:
182     break;
183   case F_POLARIZATION:
184     newparam->polarize.alpha = old[0];
185     break;
186   case F_WATER_POL:
187     newparam->wpol.al_x   =old[0];
188     newparam->wpol.al_y   =old[1];
189     newparam->wpol.al_z   =old[2];
190     newparam->wpol.rOH    =old[3];
191     newparam->wpol.rHH    =old[4];
192     newparam->wpol.rOD    =old[5];
193     break;
194   case F_THOLE_POL:
195     newparam->thole.a      = old[0];
196     newparam->thole.alpha1 = old[1];
197     newparam->thole.alpha2 = old[2];
198     if ((old[1] > 0) && (old[2] > 0))
199       newparam->thole.rfac = old[0]*pow(old[1]*old[2],-1.0/6.0);
200     else
201       newparam->thole.rfac = 1;
202     break;
203   case F_BHAM:
204     newparam->bham.a = old[0];
205     newparam->bham.b = old[1];
206     newparam->bham.c = old[2];
207     break;
208   case F_LJ14:
209     set_ljparams(comb,reppow,old[0],old[1],&newparam->lj14.c6A,&newparam->lj14.c12A);
210     set_ljparams(comb,reppow,old[2],old[3],&newparam->lj14.c6B,&newparam->lj14.c12B);
211     break;
212   case F_LJC14_Q:
213     newparam->ljc14.fqq = old[0];
214     newparam->ljc14.qi  = old[1];
215     newparam->ljc14.qj  = old[2];
216     set_ljparams(comb,reppow,old[3],old[4],&newparam->ljc14.c6,&newparam->ljc14.c12);
217     break;
218   case F_LJC_PAIRS_NB:
219     newparam->ljcnb.qi = old[0];
220     newparam->ljcnb.qj = old[1];
221     set_ljparams(comb,reppow,old[2],old[3],&newparam->ljcnb.c6,&newparam->ljcnb.c12);
222     break;
223   case F_LJ:
224     set_ljparams(comb,reppow,old[0],old[1],&newparam->lj.c6,&newparam->lj.c12);
225     break;
226   case F_PDIHS:
227   case F_PIDIHS:
228   case F_ANGRES:
229   case F_ANGRESZ:
230     newparam->pdihs.phiA = old[0];
231     newparam->pdihs.cpA  = old[1];
232                   
233     /* Dont do any checks if all parameters are zero (such interactions will be removed).
234      * Change 20100720: Amber occasionally uses negative multiplicities (mathematically OK),
235      * so I have changed the lower limit to -99 /EL
236      *
237      * Second, if the force constant is zero in both A and B states, we set the phase
238      * and multiplicity to zero too so the interaction gets removed during clean-up.
239      */ 
240     newparam->pdihs.phiB = old[3];
241     newparam->pdihs.cpB  = old[4];
242           
243     if( fabs(newparam->pdihs.cpA) < GMX_REAL_MIN && fabs(newparam->pdihs.cpB) < GMX_REAL_MIN )
244     {
245         newparam->pdihs.phiA = 0.0; 
246         newparam->pdihs.phiB = 0.0; 
247         newparam->pdihs.mult = 0; 
248     } 
249     else
250     {
251         newparam->pdihs.mult = round_check(old[2],-99,ftype,"multiplicity");
252     }
253           
254     break;
255   case F_POSRES:
256     newparam->posres.fcA[XX]   = old[0];
257     newparam->posres.fcA[YY]   = old[1];
258     newparam->posres.fcA[ZZ]   = old[2];
259     newparam->posres.fcB[XX]   = old[3];
260     newparam->posres.fcB[YY]   = old[4];
261     newparam->posres.fcB[ZZ]   = old[5];
262     newparam->posres.pos0A[XX] = old[6];
263     newparam->posres.pos0A[YY] = old[7];
264     newparam->posres.pos0A[ZZ] = old[8];
265     newparam->posres.pos0B[XX] = old[9];
266     newparam->posres.pos0B[YY] = old[10];
267     newparam->posres.pos0B[ZZ] = old[11];
268     break;
269   case F_DISRES:
270     newparam->disres.label = round_check(old[0],0,ftype,"label");
271     newparam->disres.type  = round_check(old[1],1,ftype,"type'");
272     newparam->disres.low   = old[2];
273     newparam->disres.up1   = old[3];
274     newparam->disres.up2   = old[4];
275     newparam->disres.kfac  = old[5];
276     break;
277   case F_ORIRES:
278     newparam->orires.ex    = round_check(old[0],1,ftype,"experiment") - 1;
279     newparam->orires.label = round_check(old[1],1,ftype,"label");
280     newparam->orires.power = round_check(old[2],0,ftype,"power");
281     newparam->orires.c     = old[3];
282     newparam->orires.obs   = old[4];
283     newparam->orires.kfac  = old[5];
284     break;
285   case F_DIHRES:
286     newparam->dihres.label = round_check(old[0],0,ftype,"label");
287     newparam->dihres.phi   = old[1];
288     newparam->dihres.dphi  = old[2];
289     newparam->dihres.kfac  = old[3];
290     newparam->dihres.power = round_check(old[4],0,ftype,"power");
291     break;
292   case F_RBDIHS:
293     for (i=0; (i<NR_RBDIHS); i++) {
294       newparam->rbdihs.rbcA[i]=old[i]; 
295       newparam->rbdihs.rbcB[i]=old[NR_RBDIHS+i]; 
296     }
297     break;
298   case F_FOURDIHS:
299     /* Read the dihedral parameters to temporary arrays,
300      * and convert them to the computationally faster
301      * Ryckaert-Bellemans form.
302      */   
303     /* Use conversion formula for OPLS to Ryckaert-Bellemans: */
304     newparam->rbdihs.rbcA[0]=old[1]+0.5*(old[0]+old[2]);
305     newparam->rbdihs.rbcA[1]=0.5*(3.0*old[2]-old[0]);
306     newparam->rbdihs.rbcA[2]=4.0*old[3]-old[1];
307     newparam->rbdihs.rbcA[3]=-2.0*old[2];
308     newparam->rbdihs.rbcA[4]=-4.0*old[3];
309     newparam->rbdihs.rbcA[5]=0.0;
310
311     newparam->rbdihs.rbcB[0]=old[NR_FOURDIHS+1]+0.5*(old[NR_FOURDIHS+0]+old[NR_FOURDIHS+2]);
312     newparam->rbdihs.rbcB[1]=0.5*(3.0*old[NR_FOURDIHS+2]-old[NR_FOURDIHS+0]);
313     newparam->rbdihs.rbcB[2]=4.0*old[NR_FOURDIHS+3]-old[NR_FOURDIHS+1];
314     newparam->rbdihs.rbcB[3]=-2.0*old[NR_FOURDIHS+2];
315     newparam->rbdihs.rbcB[4]=-4.0*old[NR_FOURDIHS+3];
316     newparam->rbdihs.rbcB[5]=0.0;
317     break;    
318   case F_CONSTR:
319   case F_CONSTRNC:
320     newparam->constr.dA = old[0];
321     newparam->constr.dB = old[1];
322     break;
323   case F_SETTLE:
324     newparam->settle.doh=old[0];
325     newparam->settle.dhh=old[1];
326     break;
327   case F_VSITE2:
328   case F_VSITE3:
329   case F_VSITE3FD:
330   case F_VSITE3OUT:
331   case F_VSITE4FD:
332   case F_VSITE4FDN:
333     newparam->vsite.a=old[0];
334     newparam->vsite.b=old[1];
335     newparam->vsite.c=old[2];
336     newparam->vsite.d=old[3];
337     newparam->vsite.e=old[4];
338     newparam->vsite.f=old[5];
339     break;
340   case F_VSITE3FAD:
341     newparam->vsite.a=old[1] * cos(DEG2RAD * old[0]);
342     newparam->vsite.b=old[1] * sin(DEG2RAD * old[0]);
343     newparam->vsite.c=old[2];
344     newparam->vsite.d=old[3];
345     newparam->vsite.e=old[4];
346     newparam->vsite.f=old[5];
347     break;
348   case F_VSITEN:
349     newparam->vsiten.n = round_check(old[0],1,ftype,"number of atoms");
350     newparam->vsiten.a = old[1];
351     break;
352   case F_CMAP:
353     newparam->cmap.cmapA=old[0];
354     newparam->cmap.cmapB=old[1];
355     break;
356   case F_GB12:
357   case F_GB13:
358   case F_GB14:
359     newparam->gb.sar  = old[0];
360     newparam->gb.st   = old[1];
361     newparam->gb.pi   = old[2];
362     newparam->gb.gbr  = old[3];
363     newparam->gb.bmlt = old[4];
364     break;
365   default:
366     gmx_fatal(FARGS,"unknown function type %d in %s line %d",
367               ftype,__FILE__,__LINE__);
368   }
369 }
370
371 static int enter_params(gmx_ffparams_t *ffparams, t_functype ftype,
372                         real forceparams[MAXFORCEPARAM],int comb,real reppow,
373                         int start,gmx_bool bAppend)
374 {
375   t_iparams newparam;
376   int       type;
377   
378   assign_param(ftype,&newparam,forceparams,comb,reppow);
379   if (!bAppend) {
380     for (type=start; (type<ffparams->ntypes); type++) {
381       if (ffparams->functype[type]==ftype) {
382         if (F_GB13 == ftype) {
383           /* Occasionally, the way the 1-3 reference distance is
384            * computed can lead to non-binary-identical results, but I
385            * don't know why. */
386           if ((gmx_within_tol(newparam.gb.sar,  ffparams->iparams[type].gb.sar,  1e-6)) &&
387               (gmx_within_tol(newparam.gb.st,   ffparams->iparams[type].gb.st,   1e-6)) &&
388               (gmx_within_tol(newparam.gb.pi,   ffparams->iparams[type].gb.pi,   1e-6)) &&
389               (gmx_within_tol(newparam.gb.gbr,  ffparams->iparams[type].gb.gbr,  1e-6)) &&
390               (gmx_within_tol(newparam.gb.bmlt, ffparams->iparams[type].gb.bmlt, 1e-6))) {
391             return type;
392           }
393         }
394         else {
395         if (memcmp(&newparam,&ffparams->iparams[type],(size_t)sizeof(newparam)) == 0)
396           return type;
397         }
398       }
399     }
400   }
401   else {
402     type = ffparams->ntypes;
403   }
404   if (debug)
405     fprintf(debug,"copying newparam to ffparams->iparams[%d] (ntypes=%d)\n",
406             type,ffparams->ntypes);
407   memcpy(&ffparams->iparams[type],&newparam,(size_t)sizeof(newparam));
408   
409   ffparams->ntypes++;
410   ffparams->functype[type]=ftype;
411
412   return type;
413 }
414
415 static void append_interaction(t_ilist *ilist,
416                                int type,int nral,atom_id a[MAXATOMLIST])
417 {
418   int i,where1;
419   
420   where1     = ilist->nr;
421   ilist->nr += nral+1;
422
423   ilist->iatoms[where1++]=type;
424   for (i=0; (i<nral); i++) 
425     ilist->iatoms[where1++]=a[i];
426 }
427
428 static void enter_function(t_params *p,t_functype ftype,int comb,real reppow,
429                            gmx_ffparams_t *ffparams,t_ilist *il,
430                            int *maxtypes,
431                            gmx_bool bNB,gmx_bool bAppend)
432 {
433   int     k,type,nr,nral,delta,start;
434   
435   start = ffparams->ntypes;
436   nr    = p->nr;
437   
438   for (k=0; k<nr; k++) {
439     if (*maxtypes <= ffparams->ntypes) {
440       *maxtypes += 1000;
441       srenew(ffparams->functype,*maxtypes);
442       srenew(ffparams->iparams, *maxtypes);
443       if (debug) 
444         fprintf(debug,"%s, line %d: srenewed idef->functype and idef->iparams to %d\n",
445                 __FILE__,__LINE__,*maxtypes);
446     }
447     type = enter_params(ffparams,ftype,p->param[k].c,comb,reppow,start,bAppend);
448     if (!bNB) {
449       nral  = NRAL(ftype);
450       delta = nr*(nral+1);
451       srenew(il->iatoms,il->nr+delta);
452       append_interaction(il,type,nral,p->param[k].a);
453     }
454   }
455 }
456
457 static void new_interaction_list(t_ilist *ilist)
458 {
459   int i;
460   
461   ilist->nr=0;
462   ilist->iatoms=NULL;
463 }
464
465 void convert_params(int atnr,t_params nbtypes[],
466                     t_molinfo *mi,int comb,double reppow,real fudgeQQ,
467                     gmx_mtop_t *mtop)
468 {
469   int    i,j,maxtypes,mt;
470   unsigned long  flags;
471   gmx_ffparams_t *ffp;
472   gmx_moltype_t *molt;
473   t_params *plist;
474
475   maxtypes=0;
476   
477   ffp = &mtop->ffparams;
478   ffp->ntypes   = 0;
479   ffp->atnr     = atnr;
480   ffp->functype = NULL;
481   ffp->iparams  = NULL;
482   ffp->reppow   = reppow;
483
484   enter_function(&(nbtypes[F_LJ]),  (t_functype)F_LJ,    comb,reppow,ffp,NULL,
485                  &maxtypes,TRUE,TRUE);
486   enter_function(&(nbtypes[F_BHAM]),(t_functype)F_BHAM,  comb,reppow,ffp,NULL,
487                  &maxtypes,TRUE,TRUE);
488
489   for(mt=0; mt<mtop->nmoltype; mt++) {
490     molt = &mtop->moltype[mt];
491     for(i=0; (i<F_NRE); i++) {
492       molt->ilist[i].nr     = 0;
493       molt->ilist[i].iatoms = NULL;
494       
495       plist = mi[mt].plist;
496
497       flags = interaction_function[i].flags;
498       if ((i != F_LJ) && (i != F_BHAM) && ((flags & IF_BOND) ||
499                                            (flags & IF_VSITE) ||
500                                            (flags & IF_CONSTRAINT))) {
501         enter_function(&(plist[i]),(t_functype)i,comb,reppow,
502                        ffp,&molt->ilist[i],
503                        &maxtypes,FALSE,(i == F_POSRES));
504       }
505     }
506   }
507   if (debug) {
508     fprintf(debug,"%s, line %d: There are %d functypes in idef\n",
509             __FILE__,__LINE__,ffp->ntypes);
510   }
511
512   ffp->fudgeQQ = fudgeQQ;
513 }