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 * Green Red Orange Magenta Azure Cyan Skyblue
49 #include "gmx_fatal.h"
66 static t_liedata *analyze_names(int nre, gmx_enxnm_t *names, const char *ligand)
72 /* Skip until we come to pressure */
73 for (i = 0; (i < F_NRE); i++)
75 if (strcmp(names[i].name, interaction_function[F_PRES].longname) == 0)
81 /* Now real analysis: find components of energies */
82 sprintf(self, "%s-%s", ligand, ligand);
84 for (; (i < nre); i++)
86 if ((strstr(names[i].name, ligand) != NULL) &&
87 (strstr(names[i].name, self) == NULL))
89 if (strstr(names[i].name, "LJ") != NULL)
92 srenew(ld->lj, ld->nlj);
93 ld->lj[ld->nlj-1] = i;
95 else if (strstr(names[i].name, "Coul") != NULL)
98 srenew(ld->qq, ld->nqq);
99 ld->qq[ld->nqq-1] = i;
103 printf("Using the following energy terms:\n");
105 for (i = 0; (i < ld->nlj); i++)
107 printf(" %12s", names[ld->lj[i]].name);
110 for (i = 0; (i < ld->nqq); i++)
112 printf(" %12s", names[ld->qq[i]].name);
119 real calc_lie(t_liedata *ld, t_energy ee[], real lie_lj, real lie_qq,
120 real fac_lj, real fac_qq)
126 for (i = 0; (i < ld->nlj); i++)
128 lj_tot += ee[ld->lj[i]].e;
131 for (i = 0; (i < ld->nqq); i++)
133 qq_tot += ee[ld->qq[i]].e;
136 /* And now the great LIE formula: */
137 return fac_lj*(lj_tot-lie_lj)+fac_qq*(qq_tot-lie_qq);
140 int gmx_lie(int argc, char *argv[])
142 const char *desc[] = {
143 "[TT]g_lie[tt] computes a free energy estimate based on an energy analysis",
144 "from. One needs an energy file with the following components:",
145 "Coul (A-B) LJ-SR (A-B) etc."
147 static real lie_lj = 0, lie_qq = 0, fac_lj = 0.181, fac_qq = 0.5;
148 static const char *ligand = "none";
150 { "-Elj", FALSE, etREAL, {&lie_lj},
151 "Lennard-Jones interaction between ligand and solvent" },
152 { "-Eqq", FALSE, etREAL, {&lie_qq},
153 "Coulomb interaction between ligand and solvent" },
154 { "-Clj", FALSE, etREAL, {&fac_lj},
155 "Factor in the LIE equation for Lennard-Jones component of energy" },
156 { "-Cqq", FALSE, etREAL, {&fac_qq},
157 "Factor in the LIE equation for Coulomb component of energy" },
158 { "-ligand", FALSE, etSTR, {&ligand},
159 "Name of the ligand in the energy file" }
161 #define NPA asize(pa)
164 int nre, nframes = 0, ct = 0;
168 gmx_enxnm_t *enm = NULL;
171 double lieaver = 0, lieav2 = 0;
175 { efEDR, "-f", "ener", ffREAD },
176 { efXVG, "-o", "lie", ffWRITE }
178 #define NFILE asize(fnm)
180 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
181 NFILE, fnm, NPA, pa, asize(desc), desc, 0, NULL, &oenv))
186 fp = open_enx(ftp2fn(efEDR, NFILE, fnm), "r");
187 do_enxnms(fp, &nre, &enm);
189 ld = analyze_names(nre, enm, ligand);
191 out = xvgropen(ftp2fn(efXVG, NFILE, fnm), "LIE free energy estimate",
192 "Time (ps)", "DGbind (kJ/mol)", oenv);
193 while (do_enx(fp, fr))
195 ct = check_times(fr->t);
198 lie = calc_lie(ld, fr->ener, lie_lj, lie_qq, fac_lj, fac_qq);
202 fprintf(out, "%10g %10g\n", fr->t, lie);
207 fprintf(stderr, "\n");
211 printf("DGbind = %.3f (%.3f)\n", lieaver/nframes,
212 sqrt(lieav2/nframes-sqr(lieaver/nframes)));
215 do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");