3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
30 * For more info, check our website at http://www.gromacs.org
33 * Good gRace! Old Maple Actually Chews Slate
44 #include "gmx_fatal.h"
58 real pot(real x,real qq,real c6,real cn,int npow)
60 return cn*pow(x,-npow)-c6*pow(x,-6)+qq*ONE_4PI_EPS0/x;
63 real bhpot(real x,real qq,real A,real B,real C)
65 return A*exp(-B*x) - C*pow(x,-6.0);
68 real dpot(real x,real qq,real c6,real cn,int npow)
70 return -(npow*cn*pow(x,-npow-1)-6*c6*pow(x,-7)+qq*ONE_4PI_EPS0/sqr(x));
73 int main(int argc,char *argv[])
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."
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;
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}, "sig" },
89 { "-eps", FALSE, etREAL, {&eps}, "eps" },
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 sigma for starting the plot" }
98 { efXVG, "-o", "potje", ffWRITE }
101 #define NFILE asize(fnm)
102 const char *legend[] = { "Lennard-Jones", "Buckingham" };
106 real qq,x,oldx,minimum,mval,dp[2],pp[2];
110 /* CopyRight(stdout,argv[0]);*/
111 parse_common_args(&argc,argv,PCA_CAN_VIEW,
112 NFILE,fnm,asize(pa),pa,asize(desc),
115 bBham = (opt2parg_bSet("-A",asize(pa),pa) ||
116 opt2parg_bSet("-B",asize(pa),pa) ||
117 opt2parg_bSet("-C",asize(pa),pa));
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);
126 if (opt2parg_bSet("-sig",asize(pa),pa) ||
127 opt2parg_bSet("-eps",asize(pa),pa)) {
128 c6 = 4*eps*pow(sig,6);
129 cn = 4*eps*pow(sig,npow);
131 else if (opt2parg_bSet("-c6",asize(pa),pa) ||
132 opt2parg_bSet("-cn",asize(pa),pa) ||
133 opt2parg_bSet("-pow",asize(pa),pa)) {
134 sig = pow(cn/c6,1.0/(npow-6.0));
135 eps = 0.25*c6*pow(sig,-6.0);
140 printf("c6 = %12.5e, c%d = %12.5e\n",c6,npow,cn);
141 printf("sigma = %12.5f, epsilon = %12.5f\n",sig,eps);
143 minimum = pow(npow/6.0*pow(sig,npow-6.0),1.0/(npow-6));
144 printf("Van der Waals minimum at %g, V = %g\n\n",
145 minimum,pot(minimum,0,c6,cn,npow));
146 printf("Fit of Lennard Jones (%d-6) to Buckingham:\n",npow);
149 Abh = 4*eps*pow(sig/minimum,npow)*exp(npow);
150 printf("A = %g, B = %g, C = %g\n",Abh,Bbh,Cbh);
154 fp = xvgropen(ftp2fn(efXVG,NFILE,fnm),"Potential","r (nm)","E (kJ/mol)",
156 xvgr_legend(fp,asize(legend),legend,
163 for(i=0; (i<100); i++) {
164 x = sigfac*sig+sig*i*0.02;
165 dp[next] = dpot(x,qq,c6,cn,npow);
166 fprintf(fp,"%10g %10g %10g\n",x,pot(x,qq,c6,cn,npow),
167 bhpot(x,qq,Abh,Bbh,Cbh));
169 if ((i > 0) && (dp[cur]*dp[next] < 0)) {
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",
182 do_view(oenv,ftp2fn(efXVG,NFILE,fnm),NULL);