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