2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2012,2013,2014,2015, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
42 #include "gromacs/commandline/pargs.h"
43 #include "gromacs/commandline/viewit.h"
44 #include "gromacs/fileio/xvgr.h"
45 #include "gromacs/gmxana/gmx_ana.h"
46 #include "gromacs/legacyheaders/names.h"
47 #include "gromacs/legacyheaders/txtdump.h"
48 #include "gromacs/math/units.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/utility/arraysize.h"
51 #include "gromacs/utility/smalloc.h"
53 real
pot(real x
, real qq
, real c6
, real cn
, int npow
)
55 return cn
*pow(x
, -npow
)-c6
*std::pow(x
, -6)+qq
*ONE_4PI_EPS0
/x
;
58 real
bhpot(real x
, real A
, real B
, real C
)
60 return A
*std::exp(-B
*x
) - C
*std::pow(x
, -6);
63 real
dpot(real x
, real qq
, real c6
, real cn
, int npow
)
65 return -(npow
*cn
*std::pow(x
, -npow
-1)-6*c6
*std::pow(x
, -7)+qq
*ONE_4PI_EPS0
/sqr(x
));
68 int gmx_sigeps(int argc
, char *argv
[])
70 const char *desc
[] = {
71 "[THISMODULE] is a simple utility that converts C6/C12 or C6/Cn combinations",
72 "to [GRK]sigma[grk] and [GRK]epsilon[grk], or vice versa. It can also plot the potential",
73 "in file. In addition, it makes an approximation of a Buckingham potential",
74 "to a Lennard-Jones potential."
76 static real c6
= 1.0e-3, cn
= 1.0e-6, qi
= 0, qj
= 0, sig
= 0.3, eps
= 1, sigfac
= 0.7;
77 static real Abh
= 1e5
, Bbh
= 32, Cbh
= 1e-3;
80 { "-c6", FALSE
, etREAL
, {&c6
}, "C6" },
81 { "-cn", FALSE
, etREAL
, {&cn
}, "Constant for repulsion" },
82 { "-pow", FALSE
, etINT
, {&npow
}, "Power of the repulsion term" },
83 { "-sig", FALSE
, etREAL
, {&sig
}, "[GRK]sigma[grk]" },
84 { "-eps", FALSE
, etREAL
, {&eps
}, "[GRK]epsilon[grk]" },
85 { "-A", FALSE
, etREAL
, {&Abh
}, "Buckingham A" },
86 { "-B", FALSE
, etREAL
, {&Bbh
}, "Buckingham B" },
87 { "-C", FALSE
, etREAL
, {&Cbh
}, "Buckingham C" },
88 { "-qi", FALSE
, etREAL
, {&qi
}, "qi" },
89 { "-qj", FALSE
, etREAL
, {&qj
}, "qj" },
90 { "-sigfac", FALSE
, etREAL
, {&sigfac
}, "Factor in front of [GRK]sigma[grk] for starting the plot" }
93 { efXVG
, "-o", "potje", ffWRITE
}
95 gmx_output_env_t
*oenv
;
96 #define NFILE asize(fnm)
97 const char *legend
[] = { "Lennard-Jones", "Buckingham" };
101 real qq
, x
, oldx
, minimum
, mval
, dp
[2];
105 if (!parse_common_args(&argc
, argv
, PCA_CAN_VIEW
,
106 NFILE
, fnm
, asize(pa
), pa
, asize(desc
),
107 desc
, 0, NULL
, &oenv
))
112 bBham
= (opt2parg_bSet("-A", asize(pa
), pa
) ||
113 opt2parg_bSet("-B", asize(pa
), pa
) ||
114 opt2parg_bSet("-C", asize(pa
), pa
));
119 sig
= std::pow((6.0/npow
)*std::pow(npow
/Bbh
, npow
-6), 1.0/(npow
-6));
120 eps
= c6
/(4*std::pow(sig
, 6));
121 cn
= 4*eps
*std::pow(sig
, npow
);
125 if (opt2parg_bSet("-sig", asize(pa
), pa
) ||
126 opt2parg_bSet("-eps", asize(pa
), pa
))
128 c6
= 4*eps
*std::pow(sig
, 6);
129 cn
= 4*eps
*std::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
))
135 sig
= std::pow(cn
/c6
, static_cast<real
>(1.0/(npow
-6)));
136 eps
= 0.25*c6
*std::pow(sig
, -6);
142 printf("c6 = %12.5e, c%d = %12.5e\n", c6
, npow
, cn
);
143 printf("sigma = %12.5f, epsilon = %12.5f\n", sig
, eps
);
145 minimum
= std::pow(npow
/6.0*std::pow(sig
, npow
-6), 1.0/(npow
-6));
146 printf("Van der Waals minimum at %g, V = %g\n\n",
147 minimum
, pot(minimum
, 0, c6
, cn
, npow
));
148 printf("Fit of Lennard Jones (%d-6) to Buckingham:\n", npow
);
151 Abh
= 4*eps
*std::pow(sig
/minimum
, static_cast<real
>(npow
))*std::exp(static_cast<real
>(npow
));
152 printf("A = %g, B = %g, C = %g\n", Abh
, Bbh
, Cbh
);
156 fp
= xvgropen(ftp2fn(efXVG
, NFILE
, fnm
), "Potential", "r (nm)", "E (kJ/mol)",
158 xvgr_legend(fp
, asize(legend
), legend
,
165 for (i
= 0; (i
< 100); i
++)
167 x
= sigfac
*sig
+sig
*i
*0.02;
168 dp
[next
] = dpot(x
, qq
, c6
, cn
, npow
);
169 fprintf(fp
, "%10g %10g %10g\n", x
, pot(x
, qq
, c6
, cn
, npow
),
170 bhpot(x
, Abh
, Bbh
, Cbh
));
173 if ((i
> 0) && (dp
[cur
]*dp
[next
] < 0))
175 minimum
= oldx
+ dp
[cur
]*(x
-oldx
)/(dp
[cur
]-dp
[next
]);
176 mval
= pot(minimum
, qq
, c6
, cn
, npow
);
177 printf("Van der Waals + Coulomb minimum at r = %g (nm). Value = %g (kJ/mol)\n",
187 do_view(oenv
, ftp2fn(efXVG
, NFILE
, fnm
), NULL
);