Tons of small fixes to documentation.
[alexxy/gromacs.git] / src / tools / g_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 main(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}, "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" }
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   /* CopyRight(stdout,argv[0]);*/
111   parse_common_args(&argc,argv,PCA_CAN_VIEW,
112                     NFILE,fnm,asize(pa),pa,asize(desc),
113                     desc,0,NULL,&oenv);
114
115   bBham = (opt2parg_bSet("-A",asize(pa),pa) || 
116            opt2parg_bSet("-B",asize(pa),pa) ||
117            opt2parg_bSet("-C",asize(pa),pa));
118            
119   if (bBham) {
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     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);
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       sig = pow(cn/c6,1.0/(npow-6.0));
135       eps = 0.25*c6*pow(sig,-6.0);
136     }
137     else {
138       sig = eps = 0;
139     }
140     printf("c6    = %12.5e, c%d    = %12.5e\n",c6,npow,cn);
141     printf("sigma = %12.5f, epsilon = %12.5f\n",sig,eps);
142
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);
147     Bbh = npow/minimum;
148     Cbh = c6;
149     Abh = 4*eps*pow(sig/minimum,npow)*exp(npow);
150     printf("A = %g, B = %g, C = %g\n",Abh,Bbh,Cbh);
151   }
152   qq = qi*qj;
153       
154   fp = xvgropen(ftp2fn(efXVG,NFILE,fnm),"Potential","r (nm)","E (kJ/mol)",
155                 oenv);
156   xvgr_legend(fp,asize(legend),legend,
157                 oenv);
158   if (sig == 0)
159     sig=0.25;
160   minimum = -1;
161   mval    = 0;
162   oldx    = 0;
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));
168     if (qq != 0) {
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",
173                minimum,mval);
174       }
175     }
176     cur = next;
177     oldx = x;
178       
179   }
180   ffclose(fp);
181   
182   do_view(oenv,ftp2fn(efXVG,NFILE,fnm),NULL);
183
184   thanx(stderr);  
185                
186   return 0;
187 }
188
189