Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_sigeps.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) 2012,2013,2014,2015,2017,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 <cmath>
40 #include <cstdio>
41
42 #include "gromacs/commandline/pargs.h"
43 #include "gromacs/commandline/viewit.h"
44 #include "gromacs/fileio/xvgr.h"
45 #include "gromacs/gmxana/gmx_ana.h"
46 #include "gromacs/math/functions.h"
47 #include "gromacs/math/units.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/mdtypes/md_enums.h"
50 #include "gromacs/utility/arraysize.h"
51 #include "gromacs/utility/smalloc.h"
52
53 static real pot(real x, real qq, real c6, real cn, int npow)
54 {
55     return cn * pow(x, -npow) - c6 / gmx::power6(x) + qq * ONE_4PI_EPS0 / x;
56 }
57
58 static real bhpot(real x, real A, real B, real C)
59 {
60     return A * std::exp(-B * x) - C / gmx::power6(x);
61 }
62
63 static real dpot(real x, real qq, real c6, real cn, int npow)
64 {
65     return -(npow * cn * std::pow(x, -npow - 1) - 6 * c6 / (x * gmx::power6(x))
66              + qq * ONE_4PI_EPS0 / gmx::square(x));
67 }
68
69 int gmx_sigeps(int argc, char* argv[])
70 {
71     const char* desc[] = {
72         "[THISMODULE] is a simple utility that converts C6/C12 or C6/Cn combinations",
73         "to [GRK]sigma[grk] and [GRK]epsilon[grk], or vice versa. It can also plot the potential",
74         "in  file. In addition, it makes an approximation of a Buckingham potential",
75         "to a Lennard-Jones potential."
76     };
77     static real       c6 = 1.0e-3, cn = 1.0e-6, qi = 0, qj = 0, sig = 0.3, eps = 1, sigfac = 0.7;
78     static real       Abh = 1e5, Bbh = 32, Cbh = 1e-3;
79     static int        npow  = 12;
80     t_pargs           pa[]  = { { "-c6", FALSE, etREAL, { &c6 }, "C6" },
81                      { "-cn", FALSE, etREAL, { &cn }, "Constant for repulsion" },
82                      { "-pow", FALSE, etINT, { &npow }, "Power of the repulsion term" },
83                      { "-sig", FALSE, etREAL, { &sig }, "[GRK]sigma[grk]" },
84                      { "-eps", FALSE, etREAL, { &eps }, "[GRK]epsilon[grk]" },
85                      { "-A", FALSE, etREAL, { &Abh }, "Buckingham A" },
86                      { "-B", FALSE, etREAL, { &Bbh }, "Buckingham B" },
87                      { "-C", FALSE, etREAL, { &Cbh }, "Buckingham C" },
88                      { "-qi", FALSE, etREAL, { &qi }, "qi" },
89                      { "-qj", FALSE, etREAL, { &qj }, "qj" },
90                      { "-sigfac",
91                        FALSE,
92                        etREAL,
93                        { &sigfac },
94                        "Factor in front of [GRK]sigma[grk] for starting the plot" } };
95     t_filenm          fnm[] = { { efXVG, "-o", "potje", ffWRITE } };
96     gmx_output_env_t* oenv;
97 #define NFILE asize(fnm)
98     const char* legend[] = { "Lennard-Jones", "Buckingham" };
99     FILE*       fp;
100     int         i;
101     gmx_bool    bBham;
102     real        qq, x, oldx, minimum, mval, dp[2];
103     int         cur = 0;
104 #define next (1 - cur)
105
106     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW, NFILE, fnm, asize(pa), pa, asize(desc), desc,
107                            0, nullptr, &oenv))
108     {
109         return 0;
110     }
111
112     bBham = (opt2parg_bSet("-A", asize(pa), pa) || opt2parg_bSet("-B", asize(pa), pa)
113              || opt2parg_bSet("-C", asize(pa), pa));
114
115     if (bBham)
116     {
117         c6  = Cbh;
118         sig = std::pow((6.0 / npow) * std::pow(npow / Bbh, npow - 6), 1.0 / (npow - 6));
119         eps = c6 / (4 * gmx::power6(sig));
120         cn  = 4 * eps * std::pow(sig, npow);
121     }
122     else
123     {
124         if (opt2parg_bSet("-sig", asize(pa), pa) || opt2parg_bSet("-eps", asize(pa), pa))
125         {
126             c6 = 4 * eps * gmx::power6(sig);
127             cn = 4 * eps * std::pow(sig, npow);
128         }
129         else if (opt2parg_bSet("-c6", asize(pa), pa) || opt2parg_bSet("-cn", asize(pa), pa)
130                  || opt2parg_bSet("-pow", asize(pa), pa))
131         {
132             sig = std::pow(cn / c6, static_cast<real>(1.0 / (npow - 6)));
133             eps = 0.25 * c6 / gmx::power6(sig);
134         }
135         else
136         {
137             sig = eps = 0;
138         }
139         printf("c6    = %12.5e, c%d    = %12.5e\n", c6, npow, cn);
140         printf("sigma = %12.5f, epsilon = %12.5f\n", sig, eps);
141
142         minimum = std::pow(npow / 6.0 * std::pow(sig, npow - 6), 1.0 / (npow - 6));
143         printf("Van der Waals minimum at %g, V = %g\n\n", minimum, pot(minimum, 0, c6, cn, npow));
144         printf("Fit of Lennard Jones (%d-6) to Buckingham:\n", npow);
145         Bbh = npow / minimum;
146         Cbh = c6;
147         Abh = 4 * eps * std::pow(sig / minimum, static_cast<real>(npow))
148               * std::exp(static_cast<real>(npow));
149         printf("A = %g, B = %g, C = %g\n", Abh, Bbh, Cbh);
150     }
151     qq = qi * qj;
152
153     fp = xvgropen(ftp2fn(efXVG, NFILE, fnm), "Potential", "r (nm)", "E (kJ/mol)", oenv);
154     xvgr_legend(fp, asize(legend), legend, oenv);
155     if (sig == 0)
156     {
157         sig = 0.25;
158     }
159     oldx = 0;
160     for (i = 0; (i < 100); i++)
161     {
162         x        = sigfac * sig + sig * i * 0.02;
163         dp[next] = dpot(x, qq, c6, cn, npow);
164         fprintf(fp, "%10g  %10g  %10g\n", x, pot(x, qq, c6, cn, npow), bhpot(x, Abh, Bbh, Cbh));
165         if (qq != 0)
166         {
167             if ((i > 0) && (dp[cur] * dp[next] < 0))
168             {
169                 minimum = oldx + dp[cur] * (x - oldx) / (dp[cur] - dp[next]);
170                 mval    = pot(minimum, qq, c6, cn, npow);
171                 printf("Van der Waals + Coulomb minimum at r = %g (nm). Value = %g (kJ/mol)\n",
172                        minimum, mval);
173             }
174         }
175         cur  = next;
176         oldx = x;
177     }
178     xvgrclose(fp);
179
180     do_view(oenv, ftp2fn(efXVG, NFILE, fnm), nullptr);
181
182     return 0;
183 }