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