Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / gmxana / gmx_lie.cpp
blob17dba3971f2d9a444d63aac4bdf47bb51589b066
1 /*
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, 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.
37 #include "gmxpre.h"
39 #include <cmath>
40 #include <cstdio>
41 #include <cstdlib>
42 #include <cstring>
44 #include "gromacs/commandline/pargs.h"
45 #include "gromacs/commandline/viewit.h"
46 #include "gromacs/fileio/enxio.h"
47 #include "gromacs/fileio/trxio.h"
48 #include "gromacs/fileio/xvgr.h"
49 #include "gromacs/gmxana/gmx_ana.h"
50 #include "gromacs/gmxana/gstat.h"
51 #include "gromacs/math/functions.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/topology/ifunc.h"
54 #include "gromacs/trajectory/energyframe.h"
55 #include "gromacs/utility/arraysize.h"
56 #include "gromacs/utility/futil.h"
57 #include "gromacs/utility/smalloc.h"
59 typedef struct {
60 int nlj, nqq;
61 int *lj;
62 int *qq;
63 } t_liedata;
65 static t_liedata *analyze_names(int nre, gmx_enxnm_t *names, const char *ligand)
67 int i;
68 t_liedata *ld;
69 char self[256];
71 /* Skip until we come to pressure */
72 for (i = 0; (i < F_NRE); i++)
74 if (std::strcmp(names[i].name, interaction_function[F_PRES].longname) == 0)
76 break;
80 /* Now real analysis: find components of energies */
81 sprintf(self, "%s-%s", ligand, ligand);
82 snew(ld, 1);
83 for (; (i < nre); i++)
85 if ((std::strstr(names[i].name, ligand) != nullptr) &&
86 (std::strstr(names[i].name, self) == nullptr))
88 if (std::strstr(names[i].name, "LJ") != nullptr)
90 ld->nlj++;
91 srenew(ld->lj, ld->nlj);
92 ld->lj[ld->nlj-1] = i;
94 else if (std::strstr(names[i].name, "Coul") != nullptr)
96 ld->nqq++;
97 srenew(ld->qq, ld->nqq);
98 ld->qq[ld->nqq-1] = i;
102 printf("Using the following energy terms:\n");
103 printf("LJ: ");
104 for (i = 0; (i < ld->nlj); i++)
106 printf(" %12s", names[ld->lj[i]].name);
108 printf("\nCoul:");
109 for (i = 0; (i < ld->nqq); i++)
111 printf(" %12s", names[ld->qq[i]].name);
113 printf("\n");
115 return ld;
118 static real calc_lie(t_liedata *ld, t_energy ee[], real lie_lj, real lie_qq,
119 real fac_lj, real fac_qq)
121 int i;
122 real lj_tot, qq_tot;
124 lj_tot = 0;
125 for (i = 0; (i < ld->nlj); i++)
127 lj_tot += ee[ld->lj[i]].e;
129 qq_tot = 0;
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";
153 t_pargs pa[] = {
154 { "-Elj", FALSE, etREAL, {&lie_lj},
155 "Lennard-Jones interaction between ligand and solvent" },
156 { "-Eqq", FALSE, etREAL, {&lie_qq},
157 "Coulomb interaction between ligand and solvent" },
158 { "-Clj", FALSE, etREAL, {&fac_lj},
159 "Factor in the LIE equation for Lennard-Jones component of energy" },
160 { "-Cqq", FALSE, etREAL, {&fac_qq},
161 "Factor in the LIE equation for Coulomb component of energy" },
162 { "-ligand", FALSE, etSTR, {&ligand},
163 "Name of the ligand in the energy file" }
165 #define NPA asize(pa)
167 FILE *out;
168 int nre, nframes = 0, ct = 0;
169 ener_file_t fp;
170 t_liedata *ld;
171 gmx_enxnm_t *enm = nullptr;
172 t_enxframe *fr;
173 real lie;
174 double lieaver = 0, lieav2 = 0;
175 gmx_output_env_t *oenv;
177 t_filenm fnm[] = {
178 { efEDR, "-f", "ener", ffREAD },
179 { efXVG, "-o", "lie", ffWRITE }
181 #define NFILE asize(fnm)
183 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME,
184 NFILE, fnm, NPA, pa, asize(desc), desc, 0, nullptr, &oenv))
186 return 0;
189 fp = open_enx(ftp2fn(efEDR, NFILE, fnm), "r");
190 do_enxnms(fp, &nre, &enm);
192 ld = analyze_names(nre, enm, ligand);
193 snew(fr, 1);
195 out = xvgropen(ftp2fn(efXVG, NFILE, fnm), "LIE free energy estimate",
196 "Time (ps)", "DGbind (kJ/mol)", oenv);
197 while (do_enx(fp, fr))
199 ct = check_times(fr->t);
200 if (ct == 0)
202 lie = calc_lie(ld, fr->ener, lie_lj, lie_qq, fac_lj, fac_qq);
203 lieaver += lie;
204 lieav2 += lie*lie;
205 nframes++;
206 fprintf(out, "%10g %10g\n", fr->t, lie);
209 close_enx(fp);
210 xvgrclose(out);
211 fprintf(stderr, "\n");
213 if (nframes > 0)
215 printf("DGbind = %.3f (%.3f)\n", lieaver/nframes,
216 std::sqrt(lieav2/nframes-gmx::square(lieaver/nframes)));
219 do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");
221 return 0;