Merge branch 'origin/release-2020' into merge-release-2020
[alexxy/gromacs.git] / src / gromacs / topology / idef.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,2018 by the GROMACS development team.
7  * Copyright (c) 2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include "idef.h"
41
42 #include <cstdio>
43
44 #include "gromacs/topology/forcefieldparameters.h"
45 #include "gromacs/topology/ifunc.h"
46 #include "gromacs/utility/fatalerror.h"
47 #include "gromacs/utility/smalloc.h"
48 #include "gromacs/utility/stringstream.h"
49 #include "gromacs/utility/textwriter.h"
50 #include "gromacs/utility/txtdump.h"
51
52 static void printHarmonicInteraction(gmx::TextWriter* writer,
53                                      const t_iparams& iparams,
54                                      const char*      r,
55                                      const char*      kr)
56 {
57     writer->writeLineFormatted("%sA=%12.5e, %sA=%12.5e, %sB=%12.5e, %sB=%12.5e", r,
58                                iparams.harmonic.rA, kr, iparams.harmonic.krA, r,
59                                iparams.harmonic.rB, kr, iparams.harmonic.krB);
60 }
61
62 void pr_iparams(FILE* fp, t_functype ftype, const t_iparams& iparams)
63 {
64     gmx::StringOutputStream stream;
65     {
66         gmx::TextWriter writer(&stream);
67         printInteractionParameters(&writer, ftype, iparams);
68     }
69     fputs(stream.toString().c_str(), fp);
70 }
71
72 void printInteractionParameters(gmx::TextWriter* writer, t_functype ftype, const t_iparams& iparams)
73 {
74     switch (ftype)
75     {
76         case F_ANGLES:
77         case F_G96ANGLES: printHarmonicInteraction(writer, iparams, "th", "ct"); break;
78         case F_CROSS_BOND_BONDS:
79             writer->writeLineFormatted("r1e=%15.8e, r2e=%15.8e, krr=%15.8e", iparams.cross_bb.r1e,
80                                        iparams.cross_bb.r2e, iparams.cross_bb.krr);
81             break;
82         case F_CROSS_BOND_ANGLES:
83             writer->writeLineFormatted("r1e=%15.8e, r1e=%15.8e, r3e=%15.8e, krt=%15.8e",
84                                        iparams.cross_ba.r1e, iparams.cross_ba.r2e,
85                                        iparams.cross_ba.r3e, iparams.cross_ba.krt);
86             break;
87         case F_LINEAR_ANGLES:
88             writer->writeLineFormatted("klinA=%15.8e, aA=%15.8e, klinB=%15.8e, aB=%15.8e",
89                                        iparams.linangle.klinA, iparams.linangle.aA,
90                                        iparams.linangle.klinB, iparams.linangle.aB);
91             break;
92         case F_UREY_BRADLEY:
93             writer->writeLineFormatted(
94                     "thetaA=%15.8e, kthetaA=%15.8e, r13A=%15.8e, kUBA=%15.8e, thetaB=%15.8e, "
95                     "kthetaB=%15.8e, r13B=%15.8e, kUBB=%15.8e",
96                     iparams.u_b.thetaA, iparams.u_b.kthetaA, iparams.u_b.r13A, iparams.u_b.kUBA,
97                     iparams.u_b.thetaB, iparams.u_b.kthetaB, iparams.u_b.r13B, iparams.u_b.kUBB);
98             break;
99         case F_QUARTIC_ANGLES:
100             writer->writeStringFormatted("theta=%15.8e", iparams.qangle.theta);
101             for (int i = 0; i < 5; i++)
102             {
103                 writer->writeStringFormatted(", c%c=%15.8e", '0' + i, iparams.qangle.c[i]);
104             }
105             writer->ensureLineBreak();
106             break;
107         case F_BHAM:
108             writer->writeLineFormatted("a=%15.8e, b=%15.8e, c=%15.8e", iparams.bham.a,
109                                        iparams.bham.b, iparams.bham.c);
110             break;
111         case F_BONDS:
112         case F_G96BONDS:
113         case F_HARMONIC: printHarmonicInteraction(writer, iparams, "b0", "cb"); break;
114         case F_IDIHS: printHarmonicInteraction(writer, iparams, "xi", "cx"); break;
115         case F_MORSE:
116             writer->writeLineFormatted(
117                     "b0A=%15.8e, cbA=%15.8e, betaA=%15.8e, b0B=%15.8e, cbB=%15.8e, betaB=%15.8e",
118                     iparams.morse.b0A, iparams.morse.cbA, iparams.morse.betaA, iparams.morse.b0B,
119                     iparams.morse.cbB, iparams.morse.betaB);
120             break;
121         case F_CUBICBONDS:
122             writer->writeLineFormatted("b0=%15.8e, kb=%15.8e, kcub=%15.8e", iparams.cubic.b0,
123                                        iparams.cubic.kb, iparams.cubic.kcub);
124             break;
125         case F_CONNBONDS: writer->ensureEmptyLine(); break;
126         case F_FENEBONDS:
127             writer->writeLineFormatted("bm=%15.8e, kb=%15.8e", iparams.fene.bm, iparams.fene.kb);
128             break;
129         case F_RESTRBONDS:
130             writer->writeLineFormatted(
131                     "lowA=%15.8e, up1A=%15.8e, up2A=%15.8e, kA=%15.8e, lowB=%15.8e, up1B=%15.8e, "
132                     "up2B=%15.8e, kB=%15.8e,",
133                     iparams.restraint.lowA, iparams.restraint.up1A, iparams.restraint.up2A,
134                     iparams.restraint.kA, iparams.restraint.lowB, iparams.restraint.up1B,
135                     iparams.restraint.up2B, iparams.restraint.kB);
136             break;
137         case F_TABBONDS:
138         case F_TABBONDSNC:
139         case F_TABANGLES:
140         case F_TABDIHS:
141             writer->writeLineFormatted("tab=%d, kA=%15.8e, kB=%15.8e", iparams.tab.table,
142                                        iparams.tab.kA, iparams.tab.kB);
143             break;
144         case F_POLARIZATION:
145             writer->writeLineFormatted("alpha=%15.8e", iparams.polarize.alpha);
146             break;
147         case F_ANHARM_POL:
148             writer->writeLineFormatted("alpha=%15.8e drcut=%15.8e khyp=%15.8e",
149                                        iparams.anharm_polarize.alpha, iparams.anharm_polarize.drcut,
150                                        iparams.anharm_polarize.khyp);
151             break;
152         case F_THOLE_POL:
153             writer->writeLineFormatted("a=%15.8e, alpha1=%15.8e, alpha2=%15.8e, rfac=%15.8e",
154                                        iparams.thole.a, iparams.thole.alpha1, iparams.thole.alpha2,
155                                        iparams.thole.rfac);
156             break;
157         case F_WATER_POL:
158             writer->writeLineFormatted(
159                     "al_x=%15.8e, al_y=%15.8e, al_z=%15.8e, rOH=%9.6f, rHH=%9.6f, rOD=%9.6f",
160                     iparams.wpol.al_x, iparams.wpol.al_y, iparams.wpol.al_z, iparams.wpol.rOH,
161                     iparams.wpol.rHH, iparams.wpol.rOD);
162             break;
163         case F_LJ:
164             writer->writeLineFormatted("c6=%15.8e, c12=%15.8e", iparams.lj.c6, iparams.lj.c12);
165             break;
166         case F_LJ14:
167             writer->writeLineFormatted("c6A=%15.8e, c12A=%15.8e, c6B=%15.8e, c12B=%15.8e",
168                                        iparams.lj14.c6A, iparams.lj14.c12A, iparams.lj14.c6B,
169                                        iparams.lj14.c12B);
170             break;
171         case F_LJC14_Q:
172             writer->writeLineFormatted("fqq=%15.8e, qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e",
173                                        iparams.ljc14.fqq, iparams.ljc14.qi, iparams.ljc14.qj,
174                                        iparams.ljc14.c6, iparams.ljc14.c12);
175             break;
176         case F_LJC_PAIRS_NB:
177             writer->writeLineFormatted("qi=%15.8e, qj=%15.8e, c6=%15.8e, c12=%15.8e", iparams.ljcnb.qi,
178                                        iparams.ljcnb.qj, iparams.ljcnb.c6, iparams.ljcnb.c12);
179             break;
180         case F_PDIHS:
181         case F_PIDIHS:
182         case F_ANGRES:
183         case F_ANGRESZ:
184             writer->writeLineFormatted("phiA=%15.8e, cpA=%15.8e, phiB=%15.8e, cpB=%15.8e, mult=%d",
185                                        iparams.pdihs.phiA, iparams.pdihs.cpA, iparams.pdihs.phiB,
186                                        iparams.pdihs.cpB, iparams.pdihs.mult);
187             break;
188         case F_DISRES:
189             writer->writeLineFormatted(
190                     "label=%4d, type=%1d, low=%15.8e, up1=%15.8e, up2=%15.8e, fac=%15.8e)",
191                     iparams.disres.label, iparams.disres.type, iparams.disres.low,
192                     iparams.disres.up1, iparams.disres.up2, iparams.disres.kfac);
193             break;
194         case F_ORIRES:
195             writer->writeLineFormatted(
196                     "ex=%4d, label=%d, power=%4d, c=%15.8e, obs=%15.8e, kfac=%15.8e)",
197                     iparams.orires.ex, iparams.orires.label, iparams.orires.power, iparams.orires.c,
198                     iparams.orires.obs, iparams.orires.kfac);
199             break;
200         case F_DIHRES:
201             writer->writeLineFormatted(
202                     "phiA=%15.8e, dphiA=%15.8e, kfacA=%15.8e, phiB=%15.8e, dphiB=%15.8e, "
203                     "kfacB=%15.8e",
204                     iparams.dihres.phiA, iparams.dihres.dphiA, iparams.dihres.kfacA,
205                     iparams.dihres.phiB, iparams.dihres.dphiB, iparams.dihres.kfacB);
206             break;
207         case F_POSRES:
208             writer->writeLineFormatted(
209                     "pos0A=(%15.8e,%15.8e,%15.8e), fcA=(%15.8e,%15.8e,%15.8e), "
210                     "pos0B=(%15.8e,%15.8e,%15.8e), fcB=(%15.8e,%15.8e,%15.8e)",
211                     iparams.posres.pos0A[XX], iparams.posres.pos0A[YY], iparams.posres.pos0A[ZZ],
212                     iparams.posres.fcA[XX], iparams.posres.fcA[YY], iparams.posres.fcA[ZZ],
213                     iparams.posres.pos0B[XX], iparams.posres.pos0B[YY], iparams.posres.pos0B[ZZ],
214                     iparams.posres.fcB[XX], iparams.posres.fcB[YY], iparams.posres.fcB[ZZ]);
215             break;
216         case F_FBPOSRES:
217             writer->writeLineFormatted(
218                     "pos0=(%15.8e,%15.8e,%15.8e), geometry=%d, r=%15.8e, k=%15.8e",
219                     iparams.fbposres.pos0[XX], iparams.fbposres.pos0[YY], iparams.fbposres.pos0[ZZ],
220                     iparams.fbposres.geom, iparams.fbposres.r, iparams.fbposres.k);
221             break;
222         case F_RBDIHS:
223             for (int i = 0; i < NR_RBDIHS; i++)
224             {
225                 writer->writeStringFormatted("%srbcA[%d]=%15.8e", i == 0 ? "" : ", ", i,
226                                              iparams.rbdihs.rbcA[i]);
227             }
228             writer->ensureLineBreak();
229             for (int i = 0; i < NR_RBDIHS; i++)
230             {
231                 writer->writeStringFormatted("%srbcB[%d]=%15.8e", i == 0 ? "" : ", ", i,
232                                              iparams.rbdihs.rbcB[i]);
233             }
234             writer->ensureLineBreak();
235             break;
236         case F_FOURDIHS:
237         {
238             /* Use the OPLS -> Ryckaert-Bellemans formula backwards to get
239              * the OPLS potential constants back.
240              */
241             const real* rbcA = iparams.rbdihs.rbcA;
242             const real* rbcB = iparams.rbdihs.rbcB;
243             real        VA[4], VB[4];
244
245             VA[3] = -0.25 * rbcA[4];
246             VA[2] = -0.5 * rbcA[3];
247             VA[1] = 4.0 * VA[3] - rbcA[2];
248             VA[0] = 3.0 * VA[2] - 2.0 * rbcA[1];
249
250             VB[3] = -0.25 * rbcB[4];
251             VB[2] = -0.5 * rbcB[3];
252             VB[1] = 4.0 * VB[3] - rbcB[2];
253             VB[0] = 3.0 * VB[2] - 2.0 * rbcB[1];
254
255             for (int i = 0; i < NR_FOURDIHS; i++)
256             {
257                 writer->writeStringFormatted("%sFourA[%d]=%15.8e", i == 0 ? "" : ", ", i, VA[i]);
258             }
259             writer->ensureLineBreak();
260             for (int i = 0; i < NR_FOURDIHS; i++)
261             {
262                 writer->writeStringFormatted("%sFourB[%d]=%15.8e", i == 0 ? "" : ", ", i, VB[i]);
263             }
264             writer->ensureLineBreak();
265             break;
266         }
267
268         case F_CONSTR:
269         case F_CONSTRNC:
270             writer->writeLineFormatted("dA=%15.8e, dB=%15.8e", iparams.constr.dA, iparams.constr.dB);
271             break;
272         case F_SETTLE:
273             writer->writeLineFormatted("doh=%15.8e, dhh=%15.8e", iparams.settle.doh, iparams.settle.dhh);
274             break;
275         case F_VSITE2: writer->writeLineFormatted("a=%15.8e", iparams.vsite.a); break;
276         case F_VSITE3:
277         case F_VSITE3FD:
278         case F_VSITE3FAD:
279             writer->writeLineFormatted("a=%15.8e, b=%15.8e", iparams.vsite.a, iparams.vsite.b);
280             break;
281         case F_VSITE3OUT:
282         case F_VSITE4FD:
283         case F_VSITE4FDN:
284             writer->writeLineFormatted("a=%15.8e, b=%15.8e, c=%15.8e", iparams.vsite.a,
285                                        iparams.vsite.b, iparams.vsite.c);
286             break;
287         case F_VSITEN:
288             writer->writeLineFormatted("n=%2d, a=%15.8e", iparams.vsiten.n, iparams.vsiten.a);
289             break;
290         case F_GB12_NOLONGERUSED:
291         case F_GB13_NOLONGERUSED:
292         case F_GB14_NOLONGERUSED:
293             // These could only be generated by grompp, not written in
294             // a .top file. Now that implicit solvent is not
295             // supported, they can't be generated, and the values are
296             // ignored if read from an old .tpr file. So there is
297             // nothing to print.
298             break;
299         case F_CMAP:
300             writer->writeLineFormatted("cmapA=%1d, cmapB=%1d", iparams.cmap.cmapA, iparams.cmap.cmapB);
301             break;
302         case F_RESTRANGLES: printHarmonicInteraction(writer, iparams, "ktheta", "costheta0"); break;
303         case F_RESTRDIHS:
304             writer->writeLineFormatted("phiA=%15.8e, cpA=%15.8e", iparams.pdihs.phiA, iparams.pdihs.cpA);
305             break;
306         case F_CBTDIHS:
307             writer->writeLineFormatted("kphi=%15.8e", iparams.cbtdihs.cbtcA[0]);
308             for (int i = 1; i < NR_CBTDIHS; i++)
309             {
310                 writer->writeStringFormatted(", cbtcA[%d]=%15.8e", i - 1, iparams.cbtdihs.cbtcA[i]);
311             }
312             writer->ensureLineBreak();
313             break;
314         default:
315             gmx_fatal(FARGS, "unknown function type %d (%s) in %s line %d", ftype,
316                       interaction_function[ftype].name, __FILE__, __LINE__);
317     }
318 }
319
320 template<typename T>
321 static void printIlist(FILE*             fp,
322                        int               indent,
323                        const char*       title,
324                        const t_functype* functype,
325                        const T&          ilist,
326                        gmx_bool          bShowNumbers,
327                        gmx_bool          bShowParameters,
328                        const t_iparams*  iparams)
329 {
330     int i, j, k, type, ftype;
331
332     indent = pr_title(fp, indent, title);
333     pr_indent(fp, indent);
334     fprintf(fp, "nr: %d\n", ilist.size());
335     if (!ilist.empty())
336     {
337         pr_indent(fp, indent);
338         fprintf(fp, "iatoms:\n");
339         for (i = j = 0; i < ilist.size();)
340         {
341             pr_indent(fp, indent + INDENT);
342             type  = ilist.iatoms[i];
343             ftype = functype[type];
344             if (bShowNumbers)
345             {
346                 fprintf(fp, "%d type=%d ", j, type);
347             }
348             j++;
349             printf("(%s)", interaction_function[ftype].name);
350             for (k = 0; k < interaction_function[ftype].nratoms; k++)
351             {
352                 fprintf(fp, " %3d", ilist.iatoms[i + 1 + k]);
353             }
354             if (bShowParameters)
355             {
356                 fprintf(fp, "  ");
357                 pr_iparams(fp, ftype, iparams[type]);
358             }
359             fprintf(fp, "\n");
360             i += 1 + interaction_function[ftype].nratoms;
361         }
362     }
363 }
364
365 void pr_ilist(FILE*                  fp,
366               int                    indent,
367               const char*            title,
368               const t_functype*      functype,
369               const InteractionList& ilist,
370               gmx_bool               bShowNumbers,
371               gmx_bool               bShowParameters,
372               const t_iparams*       iparams)
373 {
374     printIlist(fp, indent, title, functype, ilist, bShowNumbers, bShowParameters, iparams);
375 }
376
377 void pr_idef(FILE* fp, int indent, const char* title, const t_idef* idef, gmx_bool bShowNumbers, gmx_bool bShowParameters)
378 {
379     int i, j;
380
381     if (available(fp, idef, indent, title))
382     {
383         indent = pr_title(fp, indent, title);
384         pr_indent(fp, indent);
385         fprintf(fp, "atnr=%d\n", idef->atnr);
386         pr_indent(fp, indent);
387         fprintf(fp, "ntypes=%d\n", idef->ntypes);
388         for (i = 0; i < idef->ntypes; i++)
389         {
390             pr_indent(fp, indent + INDENT);
391             fprintf(fp, "functype[%d]=%s, ", bShowNumbers ? i : -1,
392                     interaction_function[idef->functype[i]].name);
393             pr_iparams(fp, idef->functype[i], idef->iparams[i]);
394         }
395         pr_real(fp, indent, "fudgeQQ", idef->fudgeQQ);
396
397         for (j = 0; (j < F_NRE); j++)
398         {
399             printIlist(fp, indent, interaction_function[j].longname, idef->functype, idef->il[j],
400                        bShowNumbers, bShowParameters, idef->iparams);
401         }
402     }
403 }
404
405 void init_idef(t_idef* idef)
406 {
407     idef->ntypes           = 0;
408     idef->atnr             = 0;
409     idef->functype         = nullptr;
410     idef->iparams          = nullptr;
411     idef->fudgeQQ          = 0.0;
412     idef->iparams_posres   = nullptr;
413     idef->iparams_fbposres = nullptr;
414     for (int f = 0; f < F_NRE; ++f)
415     {
416         idef->il[f].iatoms = nullptr;
417         idef->il[f].nalloc = 0;
418         idef->il[f].nr     = 0;
419     }
420 }
421
422 InteractionDefinitions::InteractionDefinitions(const gmx_ffparams_t& ffparams) :
423     iparams(ffparams.iparams),
424     functype(ffparams.functype),
425     cmap_grid(ffparams.cmap_grid)
426 {
427 }
428
429 void InteractionDefinitions::clear()
430 {
431     /* Clear the counts */
432     for (auto& ilist : il)
433     {
434         ilist.clear();
435     }
436     iparams_posres.clear();
437     iparams_fbposres.clear();
438 }
439
440 void done_idef(t_idef* idef)
441 {
442     sfree(idef->functype);
443     sfree(idef->iparams);
444     sfree(idef->iparams_posres);
445     sfree(idef->iparams_fbposres);
446     for (int f = 0; f < F_NRE; ++f)
447     {
448         sfree(idef->il[f].iatoms);
449     }
450
451     init_idef(idef);
452 }