0ed646e2e4263dcd3baaa1ac7f7c77d98bf9d620
[alexxy/gromacs.git] / src / gromacs / topology / forcefieldparameters.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,2019, 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 #include "gmxpre.h"
38
39 #include "forcefieldparameters.h"
40
41 #include "gromacs/utility/basedefinitions.h"
42 #include "gromacs/utility/txtdump.h"
43
44 static void pr_cmap(FILE* fp, int indent, const char* title, const gmx_cmap_t* cmap_grid, gmx_bool bShowNumbers)
45 {
46     int  j, nelem;
47     real dx, idx;
48
49     if (cmap_grid->grid_spacing != 0)
50     {
51         dx = 360.0 / cmap_grid->grid_spacing;
52     }
53     else
54     {
55         dx = 0;
56     }
57     nelem = cmap_grid->grid_spacing * cmap_grid->grid_spacing;
58
59     if (available(fp, cmap_grid, indent, title))
60     {
61         fprintf(fp, "%s\n", title);
62
63         for (gmx::index i = 0; i < gmx::ssize(cmap_grid->cmapdata); i++)
64         {
65             idx = -180.0;
66             fprintf(fp, "%8s %8s %8s %8s\n", "V", "dVdx", "dVdy", "d2dV");
67
68             fprintf(fp, "grid[%3zd]={\n", bShowNumbers ? i : -1);
69
70             for (j = 0; j < nelem; j++)
71             {
72                 if ((j % cmap_grid->grid_spacing) == 0)
73                 {
74                     fprintf(fp, "%8.1f\n", idx);
75                     idx += dx;
76                 }
77
78                 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j * 4]);
79                 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j * 4 + 1]);
80                 fprintf(fp, "%8.3f ", cmap_grid->cmapdata[i].cmap[j * 4 + 2]);
81                 fprintf(fp, "%8.3f\n", cmap_grid->cmapdata[i].cmap[j * 4 + 3]);
82             }
83             fprintf(fp, "\n");
84         }
85     }
86 }
87
88 void pr_ffparams(FILE* fp, int indent, const char* title, const gmx_ffparams_t* ffparams, gmx_bool bShowNumbers)
89 {
90     int i;
91
92     indent = pr_title(fp, indent, title);
93     pr_indent(fp, indent);
94     fprintf(fp, "atnr=%d\n", ffparams->atnr);
95     pr_indent(fp, indent);
96     fprintf(fp, "ntypes=%d\n", ffparams->numTypes());
97     for (i = 0; i < ffparams->numTypes(); i++)
98     {
99         pr_indent(fp, indent + INDENT);
100         fprintf(fp, "functype[%d]=%s, ", bShowNumbers ? i : -1,
101                 interaction_function[ffparams->functype[i]].name);
102         pr_iparams(fp, ffparams->functype[i], &ffparams->iparams[i]);
103     }
104     pr_double(fp, indent, "reppow", ffparams->reppow);
105     pr_real(fp, indent, "fudgeQQ", ffparams->fudgeQQ);
106     pr_cmap(fp, indent, "cmap", &ffparams->cmap_grid, bShowNumbers);
107 }