9a871fb6c1bdd17d7dc99c2b6a7504a64c0c28f7
[alexxy/gromacs.git] / src / tools / gmx_sigeps.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Good gRace! Old Maple Actually Chews Slate
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <stdio.h>
40 #include <math.h>
41 #include "typedefs.h"
42 #include "statutil.h"
43 #include "copyrite.h"
44 #include "gmx_fatal.h"
45 #include "xvgr.h"
46 #include "pdbio.h"
47 #include "macros.h"
48 #include "smalloc.h"
49 #include "vec.h"
50 #include "pbc.h"
51 #include "physics.h"
52 #include "names.h"
53 #include "txtdump.h"
54 #include "trnio.h"
55 #include "symtab.h"
56 #include "confio.h"
57 #include "gmx_ana.h"
58
59 real pot(real x, real qq, real c6, real cn, int npow)
60 {
61     return cn*pow(x, -npow)-c6*pow(x, -6)+qq*ONE_4PI_EPS0/x;
62 }
63
64 real bhpot(real x, real qq, real A, real B, real C)
65 {
66     return A*exp(-B*x) - C*pow(x, -6.0);
67 }
68
69 real dpot(real x, real qq, real c6, real cn, int npow)
70 {
71     return -(npow*cn*pow(x, -npow-1)-6*c6*pow(x, -7)+qq*ONE_4PI_EPS0/sqr(x));
72 }
73
74 int gmx_sigeps(int argc, char *argv[])
75 {
76     const char   *desc[] = {
77         "[TT]g_sigeps[tt] is a simple utility that converts C6/C12 or C6/Cn combinations",
78         "to [GRK]sigma[grk] and [GRK]epsilon[grk], or vice versa. It can also plot the potential",
79         "in  file. In addition, it makes an approximation of a Buckingham potential",
80         "to a Lennard-Jones potential."
81     };
82     static real   c6   = 1.0e-3, cn = 1.0e-6, qi = 0, qj = 0, sig = 0.3, eps = 1, sigfac = 0.7;
83     static real   Abh  = 1e5, Bbh = 32, Cbh = 1e-3;
84     static int    npow = 12;
85     t_pargs       pa[] = {
86         { "-c6",   FALSE,  etREAL,  {&c6},  "C6"   },
87         { "-cn",   FALSE,  etREAL,  {&cn},  "Constant for repulsion"   },
88         { "-pow",  FALSE,  etINT,   {&npow}, "Power of the repulsion term" },
89         { "-sig",  FALSE,  etREAL,  {&sig}, "[GRK]sigma[grk]"  },
90         { "-eps",  FALSE,  etREAL,  {&eps}, "[GRK]epsilon[grk]"  },
91         { "-A",    FALSE,  etREAL,  {&Abh}, "Buckingham A" },
92         { "-B",    FALSE,  etREAL,  {&Bbh}, "Buckingham B" },
93         { "-C",    FALSE,  etREAL,  {&Cbh}, "Buckingham C" },
94         { "-qi",   FALSE,  etREAL,  {&qi},  "qi"   },
95         { "-qj",   FALSE,  etREAL,  {&qj},  "qj"   },
96         { "-sigfac", FALSE, etREAL, {&sigfac}, "Factor in front of [GRK]sigma[grk] for starting the plot" }
97     };
98     t_filenm      fnm[] = {
99         { efXVG, "-o", "potje", ffWRITE }
100     };
101     output_env_t  oenv;
102 #define NFILE asize(fnm)
103     const char   *legend[] = { "Lennard-Jones", "Buckingham" };
104     FILE         *fp;
105     int           i;
106     gmx_bool      bBham;
107     real          qq, x, oldx, minimum, mval, dp[2], pp[2];
108     int           cur = 0;
109 #define next (1-cur)
110
111     /* CopyRight(stdout,argv[0]);*/
112     parse_common_args(&argc, argv, PCA_CAN_VIEW,
113                       NFILE, fnm, asize(pa), pa, asize(desc),
114                       desc, 0, NULL, &oenv);
115
116     bBham = (opt2parg_bSet("-A", asize(pa), pa) ||
117              opt2parg_bSet("-B", asize(pa), pa) ||
118              opt2parg_bSet("-C", asize(pa), pa));
119
120     if (bBham)
121     {
122         c6  = Cbh;
123         sig = pow((6.0/npow)*pow(npow/Bbh, npow-6.0), 1.0/(npow-6.0));
124         eps = c6/(4*pow(sig, 6.0));
125         cn  = 4*eps*pow(sig, npow);
126     }
127     else
128     {
129         if (opt2parg_bSet("-sig", asize(pa), pa) ||
130             opt2parg_bSet("-eps", asize(pa), pa))
131         {
132             c6  = 4*eps*pow(sig, 6);
133             cn  = 4*eps*pow(sig, npow);
134         }
135         else if (opt2parg_bSet("-c6", asize(pa), pa) ||
136                  opt2parg_bSet("-cn", asize(pa), pa) ||
137                  opt2parg_bSet("-pow", asize(pa), pa))
138         {
139             sig = pow(cn/c6, 1.0/(npow-6.0));
140             eps = 0.25*c6*pow(sig, -6.0);
141         }
142         else
143         {
144             sig = eps = 0;
145         }
146         printf("c6    = %12.5e, c%d    = %12.5e\n", c6, npow, cn);
147         printf("sigma = %12.5f, epsilon = %12.5f\n", sig, eps);
148
149         minimum = pow(npow/6.0*pow(sig, npow-6.0), 1.0/(npow-6));
150         printf("Van der Waals minimum at %g, V = %g\n\n",
151                minimum, pot(minimum, 0, c6, cn, npow));
152         printf("Fit of Lennard Jones (%d-6) to Buckingham:\n", npow);
153         Bbh = npow/minimum;
154         Cbh = c6;
155         Abh = 4*eps*pow(sig/minimum, npow)*exp(npow);
156         printf("A = %g, B = %g, C = %g\n", Abh, Bbh, Cbh);
157     }
158     qq = qi*qj;
159
160     fp = xvgropen(ftp2fn(efXVG, NFILE, fnm), "Potential", "r (nm)", "E (kJ/mol)",
161                   oenv);
162     xvgr_legend(fp, asize(legend), legend,
163                 oenv);
164     if (sig == 0)
165     {
166         sig = 0.25;
167     }
168     minimum = -1;
169     mval    = 0;
170     oldx    = 0;
171     for (i = 0; (i < 100); i++)
172     {
173         x        = sigfac*sig+sig*i*0.02;
174         dp[next] = dpot(x, qq, c6, cn, npow);
175         fprintf(fp, "%10g  %10g  %10g\n", x, pot(x, qq, c6, cn, npow),
176                 bhpot(x, qq, Abh, Bbh, Cbh));
177         if (qq != 0)
178         {
179             if ((i > 0) && (dp[cur]*dp[next] < 0))
180             {
181                 minimum = oldx + dp[cur]*(x-oldx)/(dp[cur]-dp[next]);
182                 mval    = pot(minimum, qq, c6, cn, npow);
183                 printf("Van der Waals + Coulomb minimum at r = %g (nm). Value = %g (kJ/mol)\n",
184                        minimum, mval);
185             }
186         }
187         cur  = next;
188         oldx = x;
189
190     }
191     ffclose(fp);
192
193     do_view(oenv, ftp2fn(efXVG, NFILE, fnm), NULL);
194
195     thanx(stderr);
196
197     return 0;
198 }