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