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