Merge branch 'release-4-6'
[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  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * 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 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <stdio.h>
43 #include <math.h>
44 #include "typedefs.h"
45 #include "statutil.h"
46 #include "gmx_fatal.h"
47 #include "xvgr.h"
48 #include "pdbio.h"
49 #include "macros.h"
50 #include "smalloc.h"
51 #include "vec.h"
52 #include "pbc.h"
53 #include "physics.h"
54 #include "names.h"
55 #include "txtdump.h"
56 #include "trnio.h"
57 #include "symtab.h"
58 #include "confio.h"
59
60 real pot(real x, real qq, real c6, real cn, int npow)
61 {
62     return cn*pow(x, -npow)-c6*pow(x, -6)+qq*ONE_4PI_EPS0/x;
63 }
64
65 real bhpot(real x, real A, real B, real C)
66 {
67     return A*exp(-B*x) - C*pow(x, -6.0);
68 }
69
70 real dpot(real x, real qq, real c6, real cn, int npow)
71 {
72     return -(npow*cn*pow(x, -npow-1)-6*c6*pow(x, -7)+qq*ONE_4PI_EPS0/sqr(x));
73 }
74
75 int gmx_sigeps(int argc, char *argv[])
76 {
77     const char   *desc[] = {
78         "[TT]g_sigeps[tt] is a simple utility that converts C6/C12 or C6/Cn combinations",
79         "to [GRK]sigma[grk] and [GRK]epsilon[grk], or vice versa. It can also plot the potential",
80         "in  file. In addition, it makes an approximation of a Buckingham potential",
81         "to a Lennard-Jones potential."
82     };
83     static real   c6   = 1.0e-3, cn = 1.0e-6, qi = 0, qj = 0, sig = 0.3, eps = 1, sigfac = 0.7;
84     static real   Abh  = 1e5, Bbh = 32, Cbh = 1e-3;
85     static int    npow = 12;
86     t_pargs       pa[] = {
87         { "-c6",   FALSE,  etREAL,  {&c6},  "C6"   },
88         { "-cn",   FALSE,  etREAL,  {&cn},  "Constant for repulsion"   },
89         { "-pow",  FALSE,  etINT,   {&npow}, "Power of the repulsion term" },
90         { "-sig",  FALSE,  etREAL,  {&sig}, "[GRK]sigma[grk]"  },
91         { "-eps",  FALSE,  etREAL,  {&eps}, "[GRK]epsilon[grk]"  },
92         { "-A",    FALSE,  etREAL,  {&Abh}, "Buckingham A" },
93         { "-B",    FALSE,  etREAL,  {&Bbh}, "Buckingham B" },
94         { "-C",    FALSE,  etREAL,  {&Cbh}, "Buckingham C" },
95         { "-qi",   FALSE,  etREAL,  {&qi},  "qi"   },
96         { "-qj",   FALSE,  etREAL,  {&qj},  "qj"   },
97         { "-sigfac", FALSE, etREAL, {&sigfac}, "Factor in front of [GRK]sigma[grk] for starting the plot" }
98     };
99     t_filenm      fnm[] = {
100         { efXVG, "-o", "potje", ffWRITE }
101     };
102     output_env_t  oenv;
103 #define NFILE asize(fnm)
104     const char   *legend[] = { "Lennard-Jones", "Buckingham" };
105     FILE         *fp;
106     int           i;
107     gmx_bool      bBham;
108     real          qq, x, oldx, minimum, mval, dp[2], pp[2];
109     int           cur = 0;
110 #define next (1-cur)
111
112     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW,
113                            NFILE, fnm, asize(pa), pa, asize(desc),
114                            desc, 0, NULL, &oenv))
115     {
116         return 0;
117     }
118
119     bBham = (opt2parg_bSet("-A", asize(pa), pa) ||
120              opt2parg_bSet("-B", asize(pa), pa) ||
121              opt2parg_bSet("-C", asize(pa), pa));
122
123     if (bBham)
124     {
125         c6  = Cbh;
126         sig = pow((6.0/npow)*pow(npow/Bbh, npow-6.0), 1.0/(npow-6.0));
127         eps = c6/(4*pow(sig, 6.0));
128         cn  = 4*eps*pow(sig, npow);
129     }
130     else
131     {
132         if (opt2parg_bSet("-sig", asize(pa), pa) ||
133             opt2parg_bSet("-eps", asize(pa), pa))
134         {
135             c6  = 4*eps*pow(sig, 6);
136             cn  = 4*eps*pow(sig, npow);
137         }
138         else if (opt2parg_bSet("-c6", asize(pa), pa) ||
139                  opt2parg_bSet("-cn", asize(pa), pa) ||
140                  opt2parg_bSet("-pow", asize(pa), pa))
141         {
142             sig = pow(cn/c6, 1.0/(npow-6.0));
143             eps = 0.25*c6*pow(sig, -6.0);
144         }
145         else
146         {
147             sig = eps = 0;
148         }
149         printf("c6    = %12.5e, c%d    = %12.5e\n", c6, npow, cn);
150         printf("sigma = %12.5f, epsilon = %12.5f\n", sig, eps);
151
152         minimum = pow(npow/6.0*pow(sig, npow-6.0), 1.0/(npow-6));
153         printf("Van der Waals minimum at %g, V = %g\n\n",
154                minimum, pot(minimum, 0, c6, cn, npow));
155         printf("Fit of Lennard Jones (%d-6) to Buckingham:\n", npow);
156         Bbh = npow/minimum;
157         Cbh = c6;
158         Abh = 4*eps*pow(sig/minimum, npow)*exp(npow);
159         printf("A = %g, B = %g, C = %g\n", Abh, Bbh, Cbh);
160     }
161     qq = qi*qj;
162
163     fp = xvgropen(ftp2fn(efXVG, NFILE, fnm), "Potential", "r (nm)", "E (kJ/mol)",
164                   oenv);
165     xvgr_legend(fp, asize(legend), legend,
166                 oenv);
167     if (sig == 0)
168     {
169         sig = 0.25;
170     }
171     minimum = -1;
172     mval    = 0;
173     oldx    = 0;
174     for (i = 0; (i < 100); i++)
175     {
176         x        = sigfac*sig+sig*i*0.02;
177         dp[next] = dpot(x, qq, c6, cn, npow);
178         fprintf(fp, "%10g  %10g  %10g\n", x, pot(x, qq, c6, cn, npow),
179                 bhpot(x, Abh, Bbh, Cbh));
180         if (qq != 0)
181         {
182             if ((i > 0) && (dp[cur]*dp[next] < 0))
183             {
184                 minimum = oldx + dp[cur]*(x-oldx)/(dp[cur]-dp[next]);
185                 mval    = pot(minimum, qq, c6, cn, npow);
186                 printf("Van der Waals + Coulomb minimum at r = %g (nm). Value = %g (kJ/mol)\n",
187                        minimum, mval);
188             }
189         }
190         cur  = next;
191         oldx = x;
192
193     }
194     ffclose(fp);
195
196     do_view(oenv, ftp2fn(efXVG, NFILE, fnm), NULL);
197
198     return 0;
199 }