Merge branch release-2018 into release-2019
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / convparm.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
39
40 #include "convparm.h"
41
42 #include <cassert>
43 #include <cmath>
44 #include <cstring>
45
46 #include "gromacs/compat/make_unique.h"
47 #include "gromacs/gmxpreprocess/gpp_atomtype.h"
48 #include "gromacs/gmxpreprocess/topio.h"
49 #include "gromacs/gmxpreprocess/toputil.h"
50 #include "gromacs/math/functions.h"
51 #include "gromacs/math/units.h"
52 #include "gromacs/math/utilities.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/mdtypes/md_enums.h"
55 #include "gromacs/topology/ifunc.h"
56 #include "gromacs/topology/topology.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/smalloc.h"
59
60 static int round_check(real r, int limit, int ftype, const char *name)
61 {
62     const int i = gmx::roundToInt(r);
63
64     if (r-i > 0.01 || r-i < -0.01)
65     {
66         gmx_fatal(FARGS, "A non-integer value (%f) was supplied for '%s' in %s",
67                   r, name, interaction_function[ftype].longname);
68     }
69
70     if (i < limit)
71     {
72         gmx_fatal(FARGS, "Value of '%s' in %s is %d, which is smaller than the minimum of %d",
73                   name, interaction_function[ftype].longname, i, limit);
74     }
75
76     return i;
77 }
78
79 static void set_ljparams(int comb, double reppow, double v, double w,
80                          real *c6, real *c12)
81 {
82     if (comb == eCOMB_ARITHMETIC || comb == eCOMB_GEOM_SIG_EPS)
83     {
84         if (v >= 0)
85         {
86             *c6  = 4*w*gmx::power6(v);
87             *c12 = 4*w*std::pow(v, reppow);
88         }
89         else
90         {
91             /* Interpret negative sigma as c6=0 and c12 with -sigma */
92             *c6  = 0;
93             *c12 = 4*w*std::pow(-v, reppow);
94         }
95     }
96     else
97     {
98         *c6  = v;
99         *c12 = w;
100     }
101 }
102
103 /* A return value of 0 means parameters were assigned successfully,
104  * returning -1 means this is an all-zero interaction that should not be added.
105  */
106 static int
107 assign_param(t_functype ftype, t_iparams *newparam,
108              real old[MAXFORCEPARAM], int comb, double reppow)
109 {
110     int      i, j;
111     bool     all_param_zero = TRUE;
112
113     /* Set to zero */
114     for (j = 0; (j < MAXFORCEPARAM); j++)
115     {
116         newparam->generic.buf[j] = 0.0;
117         /* If all parameters are zero we might not add some interaction types (selected below).
118          * We cannot apply this to ALL interactions, since many have valid reasons for having
119          * zero parameters (e.g. an index to a Cmap interaction, or LJ parameters), but
120          * we use it for angles and torsions that are typically generated automatically.
121          */
122         all_param_zero = all_param_zero && fabs(old[j]) < GMX_REAL_MIN;
123     }
124
125     if (all_param_zero)
126     {
127         if (IS_ANGLE(ftype) || IS_RESTRAINT_TYPE(ftype) || ftype == F_IDIHS ||
128             ftype == F_PDIHS || ftype == F_PIDIHS || ftype == F_RBDIHS || ftype == F_FOURDIHS)
129         {
130             return -1;
131         }
132     }
133
134     switch (ftype)
135     {
136         case F_G96ANGLES:
137             /* Post processing of input data: store cosine iso angle itself */
138             newparam->harmonic.rA  = cos(old[0]*DEG2RAD);
139             newparam->harmonic.krA = old[1];
140             newparam->harmonic.rB  = cos(old[2]*DEG2RAD);
141             newparam->harmonic.krB = old[3];
142             break;
143         case F_G96BONDS:
144             /* Post processing of input data: store square of length itself */
145             newparam->harmonic.rA  = gmx::square(old[0]);
146             newparam->harmonic.krA = old[1];
147             newparam->harmonic.rB  = gmx::square(old[2]);
148             newparam->harmonic.krB = old[3];
149             break;
150         case F_FENEBONDS:
151             newparam->fene.bm = old[0];
152             newparam->fene.kb = old[1];
153             break;
154         case F_RESTRBONDS:
155             newparam->restraint.lowA = old[0];
156             newparam->restraint.up1A = old[1];
157             newparam->restraint.up2A = old[2];
158             newparam->restraint.kA   = old[3];
159             newparam->restraint.lowB = old[4];
160             newparam->restraint.up1B = old[5];
161             newparam->restraint.up2B = old[6];
162             newparam->restraint.kB   = old[7];
163             break;
164         case F_TABBONDS:
165         case F_TABBONDSNC:
166         case F_TABANGLES:
167         case F_TABDIHS:
168             newparam->tab.table = round_check(old[0], 0, ftype, "table index");
169             newparam->tab.kA    = old[1];
170             newparam->tab.kB    = old[3];
171             break;
172         case F_CROSS_BOND_BONDS:
173             newparam->cross_bb.r1e = old[0];
174             newparam->cross_bb.r2e = old[1];
175             newparam->cross_bb.krr = old[2];
176             break;
177         case F_CROSS_BOND_ANGLES:
178             newparam->cross_ba.r1e = old[0];
179             newparam->cross_ba.r2e = old[1];
180             newparam->cross_ba.r3e = old[2];
181             newparam->cross_ba.krt = old[3];
182             break;
183         case F_UREY_BRADLEY:
184             newparam->u_b.thetaA  = old[0];
185             newparam->u_b.kthetaA = old[1];
186             newparam->u_b.r13A    = old[2];
187             newparam->u_b.kUBA    = old[3];
188             newparam->u_b.thetaB  = old[4];
189             newparam->u_b.kthetaB = old[5];
190             newparam->u_b.r13B    = old[6];
191             newparam->u_b.kUBB    = old[7];
192             break;
193         case F_QUARTIC_ANGLES:
194             newparam->qangle.theta = old[0];
195             for (i = 0; i < 5; i++)
196             {
197                 newparam->qangle.c[i] = old[i+1];
198             }
199             break;
200         case F_LINEAR_ANGLES:
201             newparam->linangle.aA    = old[0];
202             newparam->linangle.klinA = old[1];
203             newparam->linangle.aB    = old[2];
204             newparam->linangle.klinB = old[3];
205             break;
206         case F_BONDS:
207         case F_ANGLES:
208         case F_HARMONIC:
209         case F_IDIHS:
210             newparam->harmonic.rA  = old[0];
211             newparam->harmonic.krA = old[1];
212             newparam->harmonic.rB  = old[2];
213             newparam->harmonic.krB = old[3];
214             break;
215         case F_RESTRANGLES:
216             newparam->harmonic.rA  = old[0];
217             newparam->harmonic.krA = old[1];
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]*gmx::invsixthroot(old[1]*old[2]);
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_RESTRDIHS:
310             newparam->pdihs.phiA = old[0];
311             newparam->pdihs.cpA  = old[1];
312             break;
313         case F_POSRES:
314             newparam->posres.fcA[XX]   = old[0];
315             newparam->posres.fcA[YY]   = old[1];
316             newparam->posres.fcA[ZZ]   = old[2];
317             newparam->posres.fcB[XX]   = old[3];
318             newparam->posres.fcB[YY]   = old[4];
319             newparam->posres.fcB[ZZ]   = old[5];
320             newparam->posres.pos0A[XX] = old[6];
321             newparam->posres.pos0A[YY] = old[7];
322             newparam->posres.pos0A[ZZ] = old[8];
323             newparam->posres.pos0B[XX] = old[9];
324             newparam->posres.pos0B[YY] = old[10];
325             newparam->posres.pos0B[ZZ] = old[11];
326             break;
327         case F_FBPOSRES:
328             newparam->fbposres.geom     = round_check(old[0], 0, ftype, "geometry");
329             if (!(newparam->fbposres.geom > efbposresZERO && newparam->fbposres.geom < efbposresNR))
330             {
331                 gmx_fatal(FARGS, "Invalid geometry for flat-bottomed position restraint.\n"
332                           "Expected number between 1 and %d. Found %d\n", efbposresNR-1,
333                           newparam->fbposres.geom);
334             }
335             newparam->fbposres.r        = old[1];
336             newparam->fbposres.k        = old[2];
337             newparam->fbposres.pos0[XX] = old[3];
338             newparam->fbposres.pos0[YY] = old[4];
339             newparam->fbposres.pos0[ZZ] = old[5];
340             break;
341         case F_DISRES:
342             newparam->disres.label = round_check(old[0], 0, ftype, "label");
343             newparam->disres.type  = round_check(old[1], 1, ftype, "type'");
344             newparam->disres.low   = old[2];
345             newparam->disres.up1   = old[3];
346             newparam->disres.up2   = old[4];
347             newparam->disres.kfac  = old[5];
348             break;
349         case F_ORIRES:
350             newparam->orires.ex    = round_check(old[0], 1, ftype, "experiment") - 1;
351             newparam->orires.label = round_check(old[1], 1, ftype, "label");
352             newparam->orires.power = round_check(old[2], 0, ftype, "power");
353             newparam->orires.c     = old[3];
354             newparam->orires.obs   = old[4];
355             newparam->orires.kfac  = old[5];
356             break;
357         case F_DIHRES:
358             newparam->dihres.phiA   = old[0];
359             newparam->dihres.dphiA  = old[1];
360             newparam->dihres.kfacA  = old[2];
361             newparam->dihres.phiB   = old[3];
362             newparam->dihres.dphiB  = old[4];
363             newparam->dihres.kfacB  = old[5];
364             break;
365         case F_RBDIHS:
366             for (i = 0; (i < NR_RBDIHS); i++)
367             {
368                 newparam->rbdihs.rbcA[i] = old[i];
369                 newparam->rbdihs.rbcB[i] = old[NR_RBDIHS+i];
370             }
371             break;
372         case F_CBTDIHS:
373             for (i = 0; (i < NR_CBTDIHS); i++)
374             {
375                 newparam->cbtdihs.cbtcA[i] = old[i];
376             }
377             break;
378         case F_FOURDIHS:
379             /* Read the dihedral parameters to temporary arrays,
380              * and convert them to the computationally faster
381              * Ryckaert-Bellemans form.
382              */
383             /* Use conversion formula for OPLS to Ryckaert-Bellemans: */
384             newparam->rbdihs.rbcA[0] = old[1]+0.5*(old[0]+old[2]);
385             newparam->rbdihs.rbcA[1] = 0.5*(3.0*old[2]-old[0]);
386             newparam->rbdihs.rbcA[2] = 4.0*old[3]-old[1];
387             newparam->rbdihs.rbcA[3] = -2.0*old[2];
388             newparam->rbdihs.rbcA[4] = -4.0*old[3];
389             newparam->rbdihs.rbcA[5] = 0.0;
390
391             newparam->rbdihs.rbcB[0] = old[NR_FOURDIHS+1]+0.5*(old[NR_FOURDIHS+0]+old[NR_FOURDIHS+2]);
392             newparam->rbdihs.rbcB[1] = 0.5*(3.0*old[NR_FOURDIHS+2]-old[NR_FOURDIHS+0]);
393             newparam->rbdihs.rbcB[2] = 4.0*old[NR_FOURDIHS+3]-old[NR_FOURDIHS+1];
394             newparam->rbdihs.rbcB[3] = -2.0*old[NR_FOURDIHS+2];
395             newparam->rbdihs.rbcB[4] = -4.0*old[NR_FOURDIHS+3];
396             newparam->rbdihs.rbcB[5] = 0.0;
397             break;
398         case F_CONSTR:
399         case F_CONSTRNC:
400             newparam->constr.dA = old[0];
401             newparam->constr.dB = old[1];
402             break;
403         case F_SETTLE:
404             newparam->settle.doh = old[0];
405             newparam->settle.dhh = old[1];
406             break;
407         case F_VSITE2:
408         case F_VSITE3:
409         case F_VSITE3FD:
410         case F_VSITE3OUT:
411         case F_VSITE4FD:
412         case F_VSITE4FDN:
413             newparam->vsite.a = old[0];
414             newparam->vsite.b = old[1];
415             newparam->vsite.c = old[2];
416             newparam->vsite.d = old[3];
417             newparam->vsite.e = old[4];
418             newparam->vsite.f = old[5];
419             break;
420         case F_VSITE3FAD:
421             newparam->vsite.a = old[1] * cos(DEG2RAD * old[0]);
422             newparam->vsite.b = old[1] * sin(DEG2RAD * old[0]);
423             newparam->vsite.c = old[2];
424             newparam->vsite.d = old[3];
425             newparam->vsite.e = old[4];
426             newparam->vsite.f = old[5];
427             break;
428         case F_VSITEN:
429             newparam->vsiten.n = round_check(old[0], 1, ftype, "number of atoms");
430             newparam->vsiten.a = old[1];
431             break;
432         case F_CMAP:
433             newparam->cmap.cmapA = static_cast<int>(old[0]);
434             newparam->cmap.cmapB = static_cast<int>(old[1]);
435             break;
436         case F_GB12_NOLONGERUSED:
437         case F_GB13_NOLONGERUSED:
438         case F_GB14_NOLONGERUSED:
439             break;
440         default:
441             gmx_fatal(FARGS, "unknown function type %d in %s line %d",
442                       ftype, __FILE__, __LINE__);
443     }
444     return 0;
445 }
446
447 static int enter_params(gmx_ffparams_t *ffparams, t_functype ftype,
448                         real forceparams[MAXFORCEPARAM], int comb, real reppow,
449                         int start, bool bAppend)
450 {
451     t_iparams newparam;
452     int       type;
453     int       rc;
454
455     if ( (rc = assign_param(ftype, &newparam, forceparams, comb, reppow)) < 0)
456     {
457         /* -1 means this interaction is all-zero and should not be added */
458         return rc;
459     }
460
461     if (!bAppend)
462     {
463         for (type = start; (type < ffparams->numTypes()); type++)
464         {
465             if (ffparams->functype[type] == ftype)
466             {
467                 if (memcmp(&newparam, &ffparams->iparams[type], static_cast<size_t>(sizeof(newparam))) == 0)
468                 {
469                     return type;
470                 }
471             }
472         }
473     }
474     else
475     {
476         type = ffparams->numTypes();
477     }
478
479     ffparams->iparams.push_back(newparam);
480     ffparams->functype.push_back(ftype);
481
482     GMX_ASSERT(ffparams->iparams.size() == ffparams->functype.size(),
483                "sizes should match");
484
485     return type;
486 }
487
488 static void append_interaction(InteractionList *ilist,
489                                int type, int nral, const int a[MAXATOMLIST])
490 {
491     ilist->iatoms.push_back(type);
492     for (int i = 0; (i < nral); i++)
493     {
494         ilist->iatoms.push_back(a[i]);
495     }
496 }
497
498 static void enter_function(t_params *p, t_functype ftype, int comb, real reppow,
499                            gmx_ffparams_t *ffparams, InteractionList *il,
500                            bool bNB, bool bAppend)
501 {
502     int     k, type, nr, nral, start;
503
504     start = ffparams->numTypes();
505     nr    = p->nr;
506
507     for (k = 0; k < nr; k++)
508     {
509         type = enter_params(ffparams, ftype, p->param[k].c, comb, reppow, start, bAppend);
510         /* Type==-1 is used as a signal that this interaction is all-zero and should not be added. */
511         if (!bNB && type >= 0)
512         {
513             assert(il);
514             nral  = NRAL(ftype);
515             append_interaction(il, type, nral, p->param[k].a);
516         }
517     }
518 }
519
520 void convert_params(int atnr, t_params nbtypes[],
521                     t_molinfo *mi, t_molinfo *intermolecular_interactions,
522                     int comb, double reppow, real fudgeQQ,
523                     gmx_mtop_t *mtop)
524 {
525     int             i;
526     unsigned long   flags;
527     gmx_ffparams_t *ffp;
528     gmx_moltype_t  *molt;
529     t_params       *plist;
530
531     ffp           = &mtop->ffparams;
532     ffp->atnr     = atnr;
533     ffp->functype.clear();
534     ffp->iparams.clear();
535     ffp->reppow   = reppow;
536
537     enter_function(&(nbtypes[F_LJ]),  static_cast<t_functype>(F_LJ),    comb, reppow, ffp, nullptr,
538                    TRUE, TRUE);
539     enter_function(&(nbtypes[F_BHAM]), static_cast<t_functype>(F_BHAM),  comb, reppow, ffp, nullptr,
540                    TRUE, TRUE);
541
542     for (size_t mt = 0; mt < mtop->moltype.size(); mt++)
543     {
544         molt = &mtop->moltype[mt];
545         for (i = 0; (i < F_NRE); i++)
546         {
547             molt->ilist[i].iatoms.clear();
548
549             plist = mi[mt].plist;
550
551             flags = interaction_function[i].flags;
552             if ((i != F_LJ) && (i != F_BHAM) && ((flags & IF_BOND) ||
553                                                  (flags & IF_VSITE) ||
554                                                  (flags & IF_CONSTRAINT)))
555             {
556                 enter_function(&(plist[i]), static_cast<t_functype>(i), comb, reppow,
557                                ffp, &molt->ilist[i],
558                                FALSE, (i == F_POSRES  || i == F_FBPOSRES));
559             }
560         }
561     }
562
563     mtop->bIntermolecularInteractions = FALSE;
564     if (intermolecular_interactions != nullptr)
565     {
566         /* Process the intermolecular interaction list */
567         mtop->intermolecular_ilist = gmx::compat::make_unique<InteractionLists>();
568
569         for (i = 0; (i < F_NRE); i++)
570         {
571             (*mtop->intermolecular_ilist)[i].iatoms.clear();
572
573             plist = intermolecular_interactions->plist;
574
575             if (plist[i].nr > 0)
576             {
577                 flags = interaction_function[i].flags;
578                 /* For intermolecular interactions we (currently)
579                  * only support potentials.
580                  * Constraints and virtual sites would be possible,
581                  * but require a lot of extra (bug-prone) code.
582                  */
583                 if (!(flags & IF_BOND))
584                 {
585                     gmx_fatal(FARGS, "The intermolecular_interaction section may only contain bonded potentials");
586                 }
587                 else if (NRAL(i) == 1) /* e.g. position restraints */
588                 {
589                     gmx_fatal(FARGS, "Single atom interactions don't make sense in the intermolecular_interaction section, you can put them in the moleculetype section");
590                 }
591                 else if (flags & IF_CHEMBOND)
592                 {
593                     gmx_fatal(FARGS, "The intermolecular_interaction can not contain chemically bonding interactions");
594                 }
595                 else
596                 {
597                     enter_function(&(plist[i]), static_cast<t_functype>(i), comb, reppow,
598                                    ffp, &(*mtop->intermolecular_ilist)[i],
599                                    FALSE, FALSE);
600
601                     mtop->bIntermolecularInteractions = TRUE;
602                 }
603             }
604         }
605
606         if (!mtop->bIntermolecularInteractions)
607         {
608             mtop->intermolecular_ilist.reset(nullptr);
609         }
610     }
611
612     ffp->fudgeQQ = fudgeQQ;
613 }