Merge release-4-6 into master
[alexxy/gromacs.git] / src / programs / grompp / 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 /* A return value of 0 means parameters were assigned successfully,
93  * returning -1 means this is an all-zero interaction that should not be added.
94  */
95 static int
96 assign_param(t_functype ftype,t_iparams *newparam,
97                          real old[MAXFORCEPARAM],int comb,double reppow)
98 {
99   int  i,j;
100   real tmp;
101   gmx_bool all_param_zero=TRUE;
102
103   /* Set to zero */
104   for(j=0; (j<MAXFORCEPARAM); j++) 
105   {
106       newparam->generic.buf[j]=0.0;
107       /* If all parameters are zero we might not add some interaction types (selected below).
108        * We cannot apply this to ALL interactions, since many have valid reasons for having
109        * zero parameters (e.g. an index to a Cmap interaction, or LJ parameters), but
110        * we use it for angles and torsions that are typically generated automatically.
111        */
112       all_param_zero = (all_param_zero==TRUE) && fabs(old[j])<GMX_REAL_MIN;
113   }
114
115   if(all_param_zero==TRUE)
116   {
117       if(IS_ANGLE(ftype) || IS_RESTRAINT_TYPE(ftype) || ftype==F_IDIHS ||
118          ftype==F_PDIHS || ftype==F_PIDIHS || ftype==F_RBDIHS || ftype==F_FOURDIHS)
119       {
120           return -1;
121       }
122   }
123
124   switch (ftype) {
125   case F_G96ANGLES:
126     /* Post processing of input data: store cosine iso angle itself */
127     newparam->harmonic.rA =cos(old[0]*DEG2RAD);
128     newparam->harmonic.krA=old[1];
129     newparam->harmonic.rB =cos(old[2]*DEG2RAD);
130     newparam->harmonic.krB=old[3];
131     break;
132   case F_G96BONDS:
133     /* Post processing of input data: store square of length itself */
134     newparam->harmonic.rA =sqr(old[0]);
135     newparam->harmonic.krA=old[1];
136     newparam->harmonic.rB =sqr(old[2]);
137     newparam->harmonic.krB=old[3];
138     break;
139   case F_FENEBONDS:
140     newparam->fene.bm=old[0];
141     newparam->fene.kb=old[1];
142     break;
143   case F_RESTRBONDS:
144     newparam->restraint.lowA = old[0];
145     newparam->restraint.up1A = old[1];
146     newparam->restraint.up2A = old[2];
147     newparam->restraint.kA   = old[3];
148     newparam->restraint.lowB = old[4];
149     newparam->restraint.up1B = old[5];
150     newparam->restraint.up2B = old[6];
151     newparam->restraint.kB   = old[7];
152     break;
153   case F_TABBONDS:
154   case F_TABBONDSNC:
155   case F_TABANGLES:
156   case F_TABDIHS:
157     newparam->tab.table = round_check(old[0],0,ftype,"table index");
158     newparam->tab.kA    = old[1];
159     newparam->tab.kB    = old[3];
160     break;
161   case F_CROSS_BOND_BONDS:
162     newparam->cross_bb.r1e=old[0];
163     newparam->cross_bb.r2e=old[1];
164     newparam->cross_bb.krr=old[2];
165     break;
166   case F_CROSS_BOND_ANGLES:
167     newparam->cross_ba.r1e=old[0];
168     newparam->cross_ba.r2e=old[1];
169     newparam->cross_ba.r3e=old[2];
170     newparam->cross_ba.krt=old[3];
171     break;
172   case F_UREY_BRADLEY:
173     newparam->u_b.thetaA=old[0];
174     newparam->u_b.kthetaA=old[1];
175     newparam->u_b.r13A=old[2];
176     newparam->u_b.kUBA=old[3];
177     newparam->u_b.thetaB=old[4];
178     newparam->u_b.kthetaB=old[5];
179     newparam->u_b.r13B=old[6];
180     newparam->u_b.kUBB=old[7];
181     break;
182   case F_QUARTIC_ANGLES:
183     newparam->qangle.theta=old[0];
184     for(i=0; i<5; i++)
185       newparam->qangle.c[i]=old[i+1];
186     break;
187   case F_LINEAR_ANGLES:
188     newparam->linangle.aA    = old[0];
189     newparam->linangle.klinA = old[1];
190     newparam->linangle.aB    = old[2];
191     newparam->linangle.klinB = old[3];
192     break;
193   case F_BONDS:
194   case F_ANGLES:
195   case F_HARMONIC:
196   case F_IDIHS:
197     newparam->harmonic.rA =old[0];
198     newparam->harmonic.krA=old[1];
199     newparam->harmonic.rB =old[2];
200     newparam->harmonic.krB=old[3];
201     break;
202   case F_MORSE:
203     newparam->morse.b0A    =old[0];
204     newparam->morse.cbA    =old[1];
205     newparam->morse.betaA  =old[2];
206     newparam->morse.b0B    =old[3];
207     newparam->morse.cbB    =old[4];
208     newparam->morse.betaB  =old[5];
209     break;
210   case F_CUBICBONDS:
211     newparam->cubic.b0    =old[0];
212     newparam->cubic.kb    =old[1];
213     newparam->cubic.kcub  =old[2];
214     break;
215   case F_CONNBONDS:
216     break;
217   case F_POLARIZATION:
218     newparam->polarize.alpha = old[0];
219     break;
220   case F_ANHARM_POL:
221     newparam->anharm_polarize.alpha = old[0];
222     newparam->anharm_polarize.drcut = old[1];
223     newparam->anharm_polarize.khyp  = old[2];
224     break;
225   case F_WATER_POL:
226     newparam->wpol.al_x   =old[0];
227     newparam->wpol.al_y   =old[1];
228     newparam->wpol.al_z   =old[2];
229     newparam->wpol.rOH    =old[3];
230     newparam->wpol.rHH    =old[4];
231     newparam->wpol.rOD    =old[5];
232     break;
233   case F_THOLE_POL:
234     newparam->thole.a      = old[0];
235     newparam->thole.alpha1 = old[1];
236     newparam->thole.alpha2 = old[2];
237     if ((old[1] > 0) && (old[2] > 0))
238       newparam->thole.rfac = old[0]*pow(old[1]*old[2],-1.0/6.0);
239     else
240       newparam->thole.rfac = 1;
241     break;
242   case F_BHAM:
243     newparam->bham.a = old[0];
244     newparam->bham.b = old[1];
245     newparam->bham.c = old[2];
246     break;
247   case F_LJ14:
248     set_ljparams(comb,reppow,old[0],old[1],&newparam->lj14.c6A,&newparam->lj14.c12A);
249     set_ljparams(comb,reppow,old[2],old[3],&newparam->lj14.c6B,&newparam->lj14.c12B);
250     break;
251   case F_LJC14_Q:
252     newparam->ljc14.fqq = old[0];
253     newparam->ljc14.qi  = old[1];
254     newparam->ljc14.qj  = old[2];
255     set_ljparams(comb,reppow,old[3],old[4],&newparam->ljc14.c6,&newparam->ljc14.c12);
256     break;
257   case F_LJC_PAIRS_NB:
258     newparam->ljcnb.qi = old[0];
259     newparam->ljcnb.qj = old[1];
260     set_ljparams(comb,reppow,old[2],old[3],&newparam->ljcnb.c6,&newparam->ljcnb.c12);
261     break;
262   case F_LJ:
263     set_ljparams(comb,reppow,old[0],old[1],&newparam->lj.c6,&newparam->lj.c12);
264     break;
265   case F_PDIHS:
266   case F_PIDIHS:
267   case F_ANGRES:
268   case F_ANGRESZ:
269           newparam->pdihs.phiA = old[0];
270           newparam->pdihs.cpA  = old[1];
271
272           /* Change 20100720: Amber occasionally uses negative multiplicities (mathematically OK),
273            * so I have changed the lower limit to -99 /EL
274            */
275           newparam->pdihs.phiB = old[3];
276           newparam->pdihs.cpB  = old[4];
277           /* If both force constants are zero there is no interaction. Return -1 to signal
278            * this entry should NOT be added.
279            */
280           if( fabs(newparam->pdihs.cpA) < GMX_REAL_MIN && fabs(newparam->pdihs.cpB) < GMX_REAL_MIN )
281           {
282               return -1;
283           }
284     
285           newparam->pdihs.mult = round_check(old[2],-99,ftype,"multiplicity");
286           
287     break;
288   case F_POSRES:
289     newparam->posres.fcA[XX]   = old[0];
290     newparam->posres.fcA[YY]   = old[1];
291     newparam->posres.fcA[ZZ]   = old[2];
292     newparam->posres.fcB[XX]   = old[3];
293     newparam->posres.fcB[YY]   = old[4];
294     newparam->posres.fcB[ZZ]   = old[5];
295     newparam->posres.pos0A[XX] = old[6];
296     newparam->posres.pos0A[YY] = old[7];
297     newparam->posres.pos0A[ZZ] = old[8];
298     newparam->posres.pos0B[XX] = old[9];
299     newparam->posres.pos0B[YY] = old[10];
300     newparam->posres.pos0B[ZZ] = old[11];
301     break;
302   case F_FBPOSRES:
303     newparam->fbposres.geom     = round_check(old[0],0,ftype,"geometry");
304     if ( ! (newparam->fbposres.geom > efbposresZERO && newparam->fbposres.geom < efbposresNR))
305     {
306       gmx_fatal(FARGS,"Invalid geometry for flat-bottomed position restraint.\n"
307                 "Expected number between 1 and %d. Found %d\n", efbposresNR-1,
308                 newparam->fbposres.geom);
309     }
310     newparam->fbposres.r        = old[1];
311     newparam->fbposres.k        = old[2];
312     newparam->fbposres.pos0[XX] = old[3];
313     newparam->fbposres.pos0[YY] = old[4];
314     newparam->fbposres.pos0[ZZ] = old[5];
315     break;
316   case F_DISRES:
317     newparam->disres.label = round_check(old[0],0,ftype,"label");
318     newparam->disres.type  = round_check(old[1],1,ftype,"type'");
319     newparam->disres.low   = old[2];
320     newparam->disres.up1   = old[3];
321     newparam->disres.up2   = old[4];
322     newparam->disres.kfac  = old[5];
323     break;
324   case F_ORIRES:
325     newparam->orires.ex    = round_check(old[0],1,ftype,"experiment") - 1;
326     newparam->orires.label = round_check(old[1],1,ftype,"label");
327     newparam->orires.power = round_check(old[2],0,ftype,"power");
328     newparam->orires.c     = old[3];
329     newparam->orires.obs   = old[4];
330     newparam->orires.kfac  = old[5];
331     break;
332   case F_DIHRES:
333     newparam->dihres.phiA  = old[0];
334     newparam->dihres.dphiA = old[1];
335     newparam->dihres.kfacA = old[2];
336     newparam->dihres.phiB   = old[3];
337     newparam->dihres.dphiB  = old[4];
338     newparam->dihres.kfacB  = old[5];
339     break;
340   case F_RBDIHS:
341     for (i=0; (i<NR_RBDIHS); i++) {
342       newparam->rbdihs.rbcA[i]=old[i]; 
343       newparam->rbdihs.rbcB[i]=old[NR_RBDIHS+i]; 
344     }
345     break;
346   case F_FOURDIHS:
347     /* Read the dihedral parameters to temporary arrays,
348      * and convert them to the computationally faster
349      * Ryckaert-Bellemans form.
350      */   
351     /* Use conversion formula for OPLS to Ryckaert-Bellemans: */
352     newparam->rbdihs.rbcA[0]=old[1]+0.5*(old[0]+old[2]);
353     newparam->rbdihs.rbcA[1]=0.5*(3.0*old[2]-old[0]);
354     newparam->rbdihs.rbcA[2]=4.0*old[3]-old[1];
355     newparam->rbdihs.rbcA[3]=-2.0*old[2];
356     newparam->rbdihs.rbcA[4]=-4.0*old[3];
357     newparam->rbdihs.rbcA[5]=0.0;
358
359     newparam->rbdihs.rbcB[0]=old[NR_FOURDIHS+1]+0.5*(old[NR_FOURDIHS+0]+old[NR_FOURDIHS+2]);
360     newparam->rbdihs.rbcB[1]=0.5*(3.0*old[NR_FOURDIHS+2]-old[NR_FOURDIHS+0]);
361     newparam->rbdihs.rbcB[2]=4.0*old[NR_FOURDIHS+3]-old[NR_FOURDIHS+1];
362     newparam->rbdihs.rbcB[3]=-2.0*old[NR_FOURDIHS+2];
363     newparam->rbdihs.rbcB[4]=-4.0*old[NR_FOURDIHS+3];
364     newparam->rbdihs.rbcB[5]=0.0;
365     break;
366   case F_CONSTR:
367   case F_CONSTRNC:
368     newparam->constr.dA = old[0];
369     newparam->constr.dB = old[1];
370     break;
371   case F_SETTLE:
372     newparam->settle.doh=old[0];
373     newparam->settle.dhh=old[1];
374     break;
375   case F_VSITE2:
376   case F_VSITE3:
377   case F_VSITE3FD:
378   case F_VSITE3OUT:
379   case F_VSITE4FD:
380   case F_VSITE4FDN:
381     newparam->vsite.a=old[0];
382     newparam->vsite.b=old[1];
383     newparam->vsite.c=old[2];
384     newparam->vsite.d=old[3];
385     newparam->vsite.e=old[4];
386     newparam->vsite.f=old[5];
387     break;
388   case F_VSITE3FAD:
389     newparam->vsite.a=old[1] * cos(DEG2RAD * old[0]);
390     newparam->vsite.b=old[1] * sin(DEG2RAD * old[0]);
391     newparam->vsite.c=old[2];
392     newparam->vsite.d=old[3];
393     newparam->vsite.e=old[4];
394     newparam->vsite.f=old[5];
395     break;
396   case F_VSITEN:
397     newparam->vsiten.n = round_check(old[0],1,ftype,"number of atoms");
398     newparam->vsiten.a = old[1];
399     break;
400   case F_CMAP:
401     newparam->cmap.cmapA=old[0];
402     newparam->cmap.cmapB=old[1];
403     break;
404   case F_GB12:
405   case F_GB13:
406   case F_GB14:
407     newparam->gb.sar  = old[0];
408     newparam->gb.st   = old[1];
409     newparam->gb.pi   = old[2];
410     newparam->gb.gbr  = old[3];
411     newparam->gb.bmlt = old[4];
412     break;
413   default:
414     gmx_fatal(FARGS,"unknown function type %d in %s line %d",
415               ftype,__FILE__,__LINE__);
416   }
417     return 0;
418 }
419
420 static int enter_params(gmx_ffparams_t *ffparams, t_functype ftype,
421                         real forceparams[MAXFORCEPARAM],int comb,real reppow,
422                         int start,gmx_bool bAppend)
423 {
424   t_iparams newparam;
425   int       type;
426   int       rc;
427
428   if( (rc=assign_param(ftype,&newparam,forceparams,comb,reppow))<0 )
429   {
430       /* -1 means this interaction is all-zero and should not be added */
431       return rc;
432   }
433
434   if (!bAppend) {
435     for (type=start; (type<ffparams->ntypes); type++) {
436       if (ffparams->functype[type]==ftype) {
437         if (F_GB13 == ftype) {
438           /* Occasionally, the way the 1-3 reference distance is
439            * computed can lead to non-binary-identical results, but I
440            * don't know why. */
441           if ((gmx_within_tol(newparam.gb.sar,  ffparams->iparams[type].gb.sar,  1e-6)) &&
442               (gmx_within_tol(newparam.gb.st,   ffparams->iparams[type].gb.st,   1e-6)) &&
443               (gmx_within_tol(newparam.gb.pi,   ffparams->iparams[type].gb.pi,   1e-6)) &&
444               (gmx_within_tol(newparam.gb.gbr,  ffparams->iparams[type].gb.gbr,  1e-6)) &&
445               (gmx_within_tol(newparam.gb.bmlt, ffparams->iparams[type].gb.bmlt, 1e-6))) {
446             return type;
447           }
448         }
449         else {
450         if (memcmp(&newparam,&ffparams->iparams[type],(size_t)sizeof(newparam)) == 0)
451           return type;
452         }
453       }
454     }
455   }
456   else {
457     type = ffparams->ntypes;
458   }
459   if (debug)
460     fprintf(debug,"copying newparam to ffparams->iparams[%d] (ntypes=%d)\n",
461             type,ffparams->ntypes);
462   memcpy(&ffparams->iparams[type],&newparam,(size_t)sizeof(newparam));
463   
464   ffparams->ntypes++;
465   ffparams->functype[type]=ftype;
466
467   return type;
468 }
469
470 static void append_interaction(t_ilist *ilist,
471                                int type,int nral,atom_id a[MAXATOMLIST])
472 {
473   int i,where1;
474   
475   where1     = ilist->nr;
476   ilist->nr += nral+1;
477
478   ilist->iatoms[where1++]=type;
479   for (i=0; (i<nral); i++) 
480     ilist->iatoms[where1++]=a[i];
481 }
482
483 static void enter_function(t_params *p,t_functype ftype,int comb,real reppow,
484                            gmx_ffparams_t *ffparams,t_ilist *il,
485                            int *maxtypes,
486                            gmx_bool bNB,gmx_bool bAppend)
487 {
488   int     k,type,nr,nral,delta,start;
489   
490   start = ffparams->ntypes;
491   nr    = p->nr;
492   
493   for (k=0; k<nr; k++) {
494     if (*maxtypes <= ffparams->ntypes) {
495       *maxtypes += 1000;
496       srenew(ffparams->functype,*maxtypes);
497       srenew(ffparams->iparams, *maxtypes);
498       if (debug) 
499         fprintf(debug,"%s, line %d: srenewed idef->functype and idef->iparams to %d\n",
500                 __FILE__,__LINE__,*maxtypes);
501     }
502     type = enter_params(ffparams,ftype,p->param[k].c,comb,reppow,start,bAppend);
503     /* Type==-1 is used as a signal that this interaction is all-zero and should not be added. */
504     if (!bNB && type>=0) {
505       nral  = NRAL(ftype);
506       delta = nr*(nral+1);
507       srenew(il->iatoms,il->nr+delta);
508       append_interaction(il,type,nral,p->param[k].a);
509     }
510   }
511 }
512
513 void convert_params(int atnr,t_params nbtypes[],
514                     t_molinfo *mi,int comb,double reppow,real fudgeQQ,
515                     gmx_mtop_t *mtop)
516 {
517   int    i,j,maxtypes,mt;
518   unsigned long  flags;
519   gmx_ffparams_t *ffp;
520   gmx_moltype_t *molt;
521   t_params *plist;
522
523   maxtypes=0;
524   
525   ffp = &mtop->ffparams;
526   ffp->ntypes   = 0;
527   ffp->atnr     = atnr;
528   ffp->functype = NULL;
529   ffp->iparams  = NULL;
530   ffp->reppow   = reppow;
531
532   enter_function(&(nbtypes[F_LJ]),  (t_functype)F_LJ,    comb,reppow,ffp,NULL,
533                  &maxtypes,TRUE,TRUE);
534   enter_function(&(nbtypes[F_BHAM]),(t_functype)F_BHAM,  comb,reppow,ffp,NULL,
535                  &maxtypes,TRUE,TRUE);
536
537   for(mt=0; mt<mtop->nmoltype; mt++) {
538     molt = &mtop->moltype[mt];
539     for(i=0; (i<F_NRE); i++) {
540       molt->ilist[i].nr     = 0;
541       molt->ilist[i].iatoms = NULL;
542       
543       plist = mi[mt].plist;
544
545       flags = interaction_function[i].flags;
546       if ((i != F_LJ) && (i != F_BHAM) && ((flags & IF_BOND) ||
547                                            (flags & IF_VSITE) ||
548                                            (flags & IF_CONSTRAINT))) {
549         enter_function(&(plist[i]),(t_functype)i,comb,reppow,
550                        ffp,&molt->ilist[i],
551                        &maxtypes,FALSE,(i == F_POSRES  || i == F_FBPOSRES));
552       }
553     }
554   }
555   if (debug) {
556     fprintf(debug,"%s, line %d: There are %d functypes in idef\n",
557             __FILE__,__LINE__,ffp->ntypes);
558   }
559
560   ffp->fudgeQQ = fudgeQQ;
561 }