06328d5f3a50914e6b40fe3046d6c88cbfe75df1
[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.thetaA=old[0];
153     newparam->u_b.kthetaA=old[1];
154     newparam->u_b.r13A=old[2];
155     newparam->u_b.kUBA=old[3];
156     newparam->u_b.thetaB=old[4];
157     newparam->u_b.kthetaB=old[5];
158     newparam->u_b.r13B=old[6];
159     newparam->u_b.kUBB=old[7];
160     break;
161   case F_QUARTIC_ANGLES:
162     newparam->qangle.theta=old[0];
163     for(i=0; i<5; i++)
164       newparam->qangle.c[i]=old[i+1];
165     break;
166   case F_LINEAR_ANGLES:
167     newparam->linangle.aA    = old[0];
168     newparam->linangle.klinA = old[1];
169     newparam->linangle.aB    = old[2];
170     newparam->linangle.klinB = old[3];
171     break;
172   case F_BONDS:
173   case F_ANGLES:
174   case F_HARMONIC:
175   case F_IDIHS:
176     newparam->harmonic.rA =old[0];
177     newparam->harmonic.krA=old[1];
178     newparam->harmonic.rB =old[2];
179     newparam->harmonic.krB=old[3];
180     break;
181   case F_MORSE:
182     newparam->morse.b0A    =old[0];
183     newparam->morse.cbA    =old[1];
184     newparam->morse.betaA  =old[2];
185     newparam->morse.b0B    =old[3];
186     newparam->morse.cbB    =old[4];
187     newparam->morse.betaB  =old[5];
188     break;
189   case F_CUBICBONDS:
190     newparam->cubic.b0    =old[0];
191     newparam->cubic.kb    =old[1];
192     newparam->cubic.kcub  =old[2];
193     break;
194   case F_CONNBONDS:
195     break;
196   case F_POLARIZATION:
197     newparam->polarize.alpha = old[0];
198     break;
199   case F_ANHARM_POL:
200     newparam->anharm_polarize.alpha = old[0];
201     newparam->anharm_polarize.drcut = old[1];
202     newparam->anharm_polarize.khyp  = old[2];
203     break;
204   case F_WATER_POL:
205     newparam->wpol.al_x   =old[0];
206     newparam->wpol.al_y   =old[1];
207     newparam->wpol.al_z   =old[2];
208     newparam->wpol.rOH    =old[3];
209     newparam->wpol.rHH    =old[4];
210     newparam->wpol.rOD    =old[5];
211     break;
212   case F_THOLE_POL:
213     newparam->thole.a      = old[0];
214     newparam->thole.alpha1 = old[1];
215     newparam->thole.alpha2 = old[2];
216     if ((old[1] > 0) && (old[2] > 0))
217       newparam->thole.rfac = old[0]*pow(old[1]*old[2],-1.0/6.0);
218     else
219       newparam->thole.rfac = 1;
220     break;
221   case F_BHAM:
222     newparam->bham.a = old[0];
223     newparam->bham.b = old[1];
224     newparam->bham.c = old[2];
225     break;
226   case F_LJ14:
227     set_ljparams(comb,reppow,old[0],old[1],&newparam->lj14.c6A,&newparam->lj14.c12A);
228     set_ljparams(comb,reppow,old[2],old[3],&newparam->lj14.c6B,&newparam->lj14.c12B);
229     break;
230   case F_LJC14_Q:
231     newparam->ljc14.fqq = old[0];
232     newparam->ljc14.qi  = old[1];
233     newparam->ljc14.qj  = old[2];
234     set_ljparams(comb,reppow,old[3],old[4],&newparam->ljc14.c6,&newparam->ljc14.c12);
235     break;
236   case F_LJC_PAIRS_NB:
237     newparam->ljcnb.qi = old[0];
238     newparam->ljcnb.qj = old[1];
239     set_ljparams(comb,reppow,old[2],old[3],&newparam->ljcnb.c6,&newparam->ljcnb.c12);
240     break;
241   case F_LJ:
242     set_ljparams(comb,reppow,old[0],old[1],&newparam->lj.c6,&newparam->lj.c12);
243     break;
244   case F_PDIHS:
245   case F_PIDIHS:
246   case F_ANGRES:
247   case F_ANGRESZ:
248     newparam->pdihs.phiA = old[0];
249     newparam->pdihs.cpA  = old[1];
250                   
251     /* Dont do any checks if all parameters are zero (such interactions will be removed).
252      * Change 20100720: Amber occasionally uses negative multiplicities (mathematically OK),
253      * so I have changed the lower limit to -99 /EL
254      *
255      * Second, if the force constant is zero in both A and B states, we set the phase
256      * and multiplicity to zero too so the interaction gets removed during clean-up.
257      */ 
258     newparam->pdihs.phiB = old[3];
259     newparam->pdihs.cpB  = old[4];
260           
261     if( fabs(newparam->pdihs.cpA) < GMX_REAL_MIN && fabs(newparam->pdihs.cpB) < GMX_REAL_MIN )
262     {
263         newparam->pdihs.phiA = 0.0; 
264         newparam->pdihs.phiB = 0.0; 
265         newparam->pdihs.mult = 0; 
266     } 
267     else
268     {
269         newparam->pdihs.mult = round_check(old[2],-99,ftype,"multiplicity");
270     }
271           
272     break;
273   case F_POSRES:
274     newparam->posres.fcA[XX]   = old[0];
275     newparam->posres.fcA[YY]   = old[1];
276     newparam->posres.fcA[ZZ]   = old[2];
277     newparam->posres.fcB[XX]   = old[3];
278     newparam->posres.fcB[YY]   = old[4];
279     newparam->posres.fcB[ZZ]   = old[5];
280     newparam->posres.pos0A[XX] = old[6];
281     newparam->posres.pos0A[YY] = old[7];
282     newparam->posres.pos0A[ZZ] = old[8];
283     newparam->posres.pos0B[XX] = old[9];
284     newparam->posres.pos0B[YY] = old[10];
285     newparam->posres.pos0B[ZZ] = old[11];
286     break;
287   case F_DISRES:
288     newparam->disres.label = round_check(old[0],0,ftype,"label");
289     newparam->disres.type  = round_check(old[1],1,ftype,"type'");
290     newparam->disres.low   = old[2];
291     newparam->disres.up1   = old[3];
292     newparam->disres.up2   = old[4];
293     newparam->disres.kfac  = old[5];
294     break;
295   case F_ORIRES:
296     newparam->orires.ex    = round_check(old[0],1,ftype,"experiment") - 1;
297     newparam->orires.label = round_check(old[1],1,ftype,"label");
298     newparam->orires.power = round_check(old[2],0,ftype,"power");
299     newparam->orires.c     = old[3];
300     newparam->orires.obs   = old[4];
301     newparam->orires.kfac  = old[5];
302     break;
303   case F_DIHRES:
304     newparam->dihres.phiA  = old[0];
305     newparam->dihres.dphiA = old[1];
306     newparam->dihres.kfacA = old[2];
307     newparam->dihres.phiB   = old[3];
308     newparam->dihres.dphiB  = old[4];
309     newparam->dihres.kfacB  = old[5];
310     break;
311   case F_RBDIHS:
312     for (i=0; (i<NR_RBDIHS); i++) {
313       newparam->rbdihs.rbcA[i]=old[i]; 
314       newparam->rbdihs.rbcB[i]=old[NR_RBDIHS+i]; 
315     }
316     break;
317   case F_FOURDIHS:
318     /* Read the dihedral parameters to temporary arrays,
319      * and convert them to the computationally faster
320      * Ryckaert-Bellemans form.
321      */   
322     /* Use conversion formula for OPLS to Ryckaert-Bellemans: */
323     newparam->rbdihs.rbcA[0]=old[1]+0.5*(old[0]+old[2]);
324     newparam->rbdihs.rbcA[1]=0.5*(3.0*old[2]-old[0]);
325     newparam->rbdihs.rbcA[2]=4.0*old[3]-old[1];
326     newparam->rbdihs.rbcA[3]=-2.0*old[2];
327     newparam->rbdihs.rbcA[4]=-4.0*old[3];
328     newparam->rbdihs.rbcA[5]=0.0;
329
330     newparam->rbdihs.rbcB[0]=old[NR_FOURDIHS+1]+0.5*(old[NR_FOURDIHS+0]+old[NR_FOURDIHS+2]);
331     newparam->rbdihs.rbcB[1]=0.5*(3.0*old[NR_FOURDIHS+2]-old[NR_FOURDIHS+0]);
332     newparam->rbdihs.rbcB[2]=4.0*old[NR_FOURDIHS+3]-old[NR_FOURDIHS+1];
333     newparam->rbdihs.rbcB[3]=-2.0*old[NR_FOURDIHS+2];
334     newparam->rbdihs.rbcB[4]=-4.0*old[NR_FOURDIHS+3];
335     newparam->rbdihs.rbcB[5]=0.0;
336     break;    
337   case F_CONSTR:
338   case F_CONSTRNC:
339     newparam->constr.dA = old[0];
340     newparam->constr.dB = old[1];
341     break;
342   case F_SETTLE:
343     newparam->settle.doh=old[0];
344     newparam->settle.dhh=old[1];
345     break;
346   case F_VSITE2:
347   case F_VSITE3:
348   case F_VSITE3FD:
349   case F_VSITE3OUT:
350   case F_VSITE4FD:
351   case F_VSITE4FDN:
352     newparam->vsite.a=old[0];
353     newparam->vsite.b=old[1];
354     newparam->vsite.c=old[2];
355     newparam->vsite.d=old[3];
356     newparam->vsite.e=old[4];
357     newparam->vsite.f=old[5];
358     break;
359   case F_VSITE3FAD:
360     newparam->vsite.a=old[1] * cos(DEG2RAD * old[0]);
361     newparam->vsite.b=old[1] * sin(DEG2RAD * old[0]);
362     newparam->vsite.c=old[2];
363     newparam->vsite.d=old[3];
364     newparam->vsite.e=old[4];
365     newparam->vsite.f=old[5];
366     break;
367   case F_VSITEN:
368     newparam->vsiten.n = round_check(old[0],1,ftype,"number of atoms");
369     newparam->vsiten.a = old[1];
370     break;
371   case F_CMAP:
372     newparam->cmap.cmapA=old[0];
373     newparam->cmap.cmapB=old[1];
374     break;
375   case F_GB12:
376   case F_GB13:
377   case F_GB14:
378     newparam->gb.sar  = old[0];
379     newparam->gb.st   = old[1];
380     newparam->gb.pi   = old[2];
381     newparam->gb.gbr  = old[3];
382     newparam->gb.bmlt = old[4];
383     break;
384   default:
385     gmx_fatal(FARGS,"unknown function type %d in %s line %d",
386               ftype,__FILE__,__LINE__);
387   }
388 }
389
390 static int enter_params(gmx_ffparams_t *ffparams, t_functype ftype,
391                         real forceparams[MAXFORCEPARAM],int comb,real reppow,
392                         int start,gmx_bool bAppend)
393 {
394   t_iparams newparam;
395   int       type;
396   
397   assign_param(ftype,&newparam,forceparams,comb,reppow);
398   if (!bAppend) {
399     for (type=start; (type<ffparams->ntypes); type++) {
400       if (ffparams->functype[type]==ftype) {
401         if (F_GB13 == ftype) {
402           /* Occasionally, the way the 1-3 reference distance is
403            * computed can lead to non-binary-identical results, but I
404            * don't know why. */
405           if ((gmx_within_tol(newparam.gb.sar,  ffparams->iparams[type].gb.sar,  1e-6)) &&
406               (gmx_within_tol(newparam.gb.st,   ffparams->iparams[type].gb.st,   1e-6)) &&
407               (gmx_within_tol(newparam.gb.pi,   ffparams->iparams[type].gb.pi,   1e-6)) &&
408               (gmx_within_tol(newparam.gb.gbr,  ffparams->iparams[type].gb.gbr,  1e-6)) &&
409               (gmx_within_tol(newparam.gb.bmlt, ffparams->iparams[type].gb.bmlt, 1e-6))) {
410             return type;
411           }
412         }
413         else {
414         if (memcmp(&newparam,&ffparams->iparams[type],(size_t)sizeof(newparam)) == 0)
415           return type;
416         }
417       }
418     }
419   }
420   else {
421     type = ffparams->ntypes;
422   }
423   if (debug)
424     fprintf(debug,"copying newparam to ffparams->iparams[%d] (ntypes=%d)\n",
425             type,ffparams->ntypes);
426   memcpy(&ffparams->iparams[type],&newparam,(size_t)sizeof(newparam));
427   
428   ffparams->ntypes++;
429   ffparams->functype[type]=ftype;
430
431   return type;
432 }
433
434 static void append_interaction(t_ilist *ilist,
435                                int type,int nral,atom_id a[MAXATOMLIST])
436 {
437   int i,where1;
438   
439   where1     = ilist->nr;
440   ilist->nr += nral+1;
441
442   ilist->iatoms[where1++]=type;
443   for (i=0; (i<nral); i++) 
444     ilist->iatoms[where1++]=a[i];
445 }
446
447 static void enter_function(t_params *p,t_functype ftype,int comb,real reppow,
448                            gmx_ffparams_t *ffparams,t_ilist *il,
449                            int *maxtypes,
450                            gmx_bool bNB,gmx_bool bAppend)
451 {
452   int     k,type,nr,nral,delta,start;
453   
454   start = ffparams->ntypes;
455   nr    = p->nr;
456   
457   for (k=0; k<nr; k++) {
458     if (*maxtypes <= ffparams->ntypes) {
459       *maxtypes += 1000;
460       srenew(ffparams->functype,*maxtypes);
461       srenew(ffparams->iparams, *maxtypes);
462       if (debug) 
463         fprintf(debug,"%s, line %d: srenewed idef->functype and idef->iparams to %d\n",
464                 __FILE__,__LINE__,*maxtypes);
465     }
466     type = enter_params(ffparams,ftype,p->param[k].c,comb,reppow,start,bAppend);
467     if (!bNB) {
468       nral  = NRAL(ftype);
469       delta = nr*(nral+1);
470       srenew(il->iatoms,il->nr+delta);
471       append_interaction(il,type,nral,p->param[k].a);
472     }
473   }
474 }
475
476 void convert_params(int atnr,t_params nbtypes[],
477                     t_molinfo *mi,int comb,double reppow,real fudgeQQ,
478                     gmx_mtop_t *mtop)
479 {
480   int    i,j,maxtypes,mt;
481   unsigned long  flags;
482   gmx_ffparams_t *ffp;
483   gmx_moltype_t *molt;
484   t_params *plist;
485
486   maxtypes=0;
487   
488   ffp = &mtop->ffparams;
489   ffp->ntypes   = 0;
490   ffp->atnr     = atnr;
491   ffp->functype = NULL;
492   ffp->iparams  = NULL;
493   ffp->reppow   = reppow;
494
495   enter_function(&(nbtypes[F_LJ]),  (t_functype)F_LJ,    comb,reppow,ffp,NULL,
496                  &maxtypes,TRUE,TRUE);
497   enter_function(&(nbtypes[F_BHAM]),(t_functype)F_BHAM,  comb,reppow,ffp,NULL,
498                  &maxtypes,TRUE,TRUE);
499
500   for(mt=0; mt<mtop->nmoltype; mt++) {
501     molt = &mtop->moltype[mt];
502     for(i=0; (i<F_NRE); i++) {
503       molt->ilist[i].nr     = 0;
504       molt->ilist[i].iatoms = NULL;
505       
506       plist = mi[mt].plist;
507
508       flags = interaction_function[i].flags;
509       if ((i != F_LJ) && (i != F_BHAM) && ((flags & IF_BOND) ||
510                                            (flags & IF_VSITE) ||
511                                            (flags & IF_CONSTRAINT))) {
512         enter_function(&(plist[i]),(t_functype)i,comb,reppow,
513                        ffp,&molt->ilist[i],
514                        &maxtypes,FALSE,(i == F_POSRES));
515       }
516     }
517   }
518   if (debug) {
519     fprintf(debug,"%s, line %d: There are %d functypes in idef\n",
520             __FILE__,__LINE__,ffp->ntypes);
521   }
522
523   ffp->fudgeQQ = fudgeQQ;
524 }