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) 2013,2014,2015,2017,2018 by the GROMACS development team.
7 * Copyright (c) 2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "gromacs/commandline/pargs.h"
46 #include "gromacs/commandline/viewit.h"
47 #include "gromacs/fileio/enxio.h"
48 #include "gromacs/fileio/trxio.h"
49 #include "gromacs/fileio/xvgr.h"
50 #include "gromacs/gmxana/gmx_ana.h"
51 #include "gromacs/gmxana/gstat.h"
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/topology/ifunc.h"
55 #include "gromacs/trajectory/energyframe.h"
56 #include "gromacs/utility/arraysize.h"
57 #include "gromacs/utility/futil.h"
58 #include "gromacs/utility/smalloc.h"
67 static t_liedata* analyze_names(int nre, gmx_enxnm_t* names, const char* ligand)
73 /* Skip until we come to pressure */
74 for (i = 0; (i < nre); i++)
76 if (std::strcmp(names[i].name, interaction_function[F_PRES].longname) == 0)
82 /* Now real analysis: find components of energies */
83 sprintf(self, "%s-%s", ligand, ligand);
85 for (; (i < nre); i++)
87 if ((std::strstr(names[i].name, ligand) != nullptr) && (std::strstr(names[i].name, self) == nullptr))
89 if (std::strstr(names[i].name, "LJ") != nullptr)
92 srenew(ld->lj, ld->nlj);
93 ld->lj[ld->nlj - 1] = i;
95 else if (std::strstr(names[i].name, "Coul") != nullptr)
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 static real calc_lie(t_liedata* ld, t_energy ee[], real lie_lj, real lie_qq, real fac_lj, real fac_qq)
125 for (i = 0; (i < ld->nlj); i++)
127 lj_tot += ee[ld->lj[i]].e;
130 for (i = 0; (i < ld->nqq); i++)
132 qq_tot += ee[ld->qq[i]].e;
135 /* And now the great LIE formula: */
136 return fac_lj * (lj_tot - lie_lj) + fac_qq * (qq_tot - lie_qq);
139 int gmx_lie(int argc, char* argv[])
141 const char* desc[] = {
142 "[THISMODULE] computes a free energy estimate based on an energy analysis",
143 "from nonbonded energies. One needs an energy file with the following components:",
144 "Coul-(A-B) LJ-SR (A-B) etc.[PAR]",
145 "To utilize [TT]g_lie[tt] correctly, two simulations are required: one with the",
146 "molecule of interest bound to its receptor and one with the molecule in water.",
147 "Both need to utilize [TT]energygrps[tt] such that Coul-SR(A-B), LJ-SR(A-B), etc. terms",
148 "are written to the [REF].edr[ref] file. Values from the molecule-in-water simulation",
149 "are necessary for supplying suitable values for -Elj and -Eqq."
151 static real lie_lj = 0, lie_qq = 0, fac_lj = 0.181, fac_qq = 0.5;
152 static const char* ligand = "none";
158 "Lennard-Jones interaction between ligand and solvent" },
159 { "-Eqq", FALSE, etREAL, { &lie_qq }, "Coulomb interaction between ligand and solvent" },
164 "Factor in the LIE equation for Lennard-Jones component of energy" },
169 "Factor in the LIE equation for Coulomb component of energy" },
170 { "-ligand", FALSE, etSTR, { &ligand }, "Name of the ligand in the energy file" }
172 #define NPA asize(pa)
175 int nre, nframes = 0, ct = 0;
178 gmx_enxnm_t* enm = nullptr;
181 double lieaver = 0, lieav2 = 0;
182 gmx_output_env_t* oenv;
184 t_filenm fnm[] = { { efEDR, "-f", "ener", ffREAD }, { efXVG, "-o", "lie", ffWRITE } };
185 #define NFILE asize(fnm)
187 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME, NFILE, fnm, NPA, pa,
188 asize(desc), desc, 0, nullptr, &oenv))
193 fp = open_enx(ftp2fn(efEDR, NFILE, fnm), "r");
194 do_enxnms(fp, &nre, &enm);
196 ld = analyze_names(nre, enm, ligand);
199 out = xvgropen(ftp2fn(efXVG, NFILE, fnm), "LIE free energy estimate", "Time (ps)",
200 "DGbind (kJ/mol)", oenv);
201 while (do_enx(fp, fr))
203 ct = check_times(fr->t);
206 lie = calc_lie(ld, fr->ener, lie_lj, lie_qq, fac_lj, fac_qq);
210 fprintf(out, "%10g %10g\n", fr->t, lie);
215 fprintf(stderr, "\n");
219 printf("DGbind = %.3f (%.3f)\n", lieaver / nframes,
220 std::sqrt(lieav2 / nframes - gmx::square(lieaver / nframes)));
223 do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");