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