d19ee95a309108b1a5f1de74e799329ea1ad571e
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_lie.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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,2019, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include <cmath>
40 #include <cstdio>
41 #include <cstdlib>
42 #include <cstring>
43
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"
58
59 typedef struct
60 {
61     int  nlj, nqq;
62     int* lj;
63     int* qq;
64 } t_liedata;
65
66 static t_liedata* analyze_names(int nre, gmx_enxnm_t* names, const char* ligand)
67 {
68     int        i;
69     t_liedata* ld;
70     char       self[256];
71
72     /* Skip until we come to pressure */
73     for (i = 0; (i < F_NRE); i++)
74     {
75         if (std::strcmp(names[i].name, interaction_function[F_PRES].longname) == 0)
76         {
77             break;
78         }
79     }
80
81     /* Now real analysis: find components of energies */
82     sprintf(self, "%s-%s", ligand, ligand);
83     snew(ld, 1);
84     for (; (i < nre); i++)
85     {
86         if ((std::strstr(names[i].name, ligand) != nullptr) && (std::strstr(names[i].name, self) == nullptr))
87         {
88             if (std::strstr(names[i].name, "LJ") != nullptr)
89             {
90                 ld->nlj++;
91                 srenew(ld->lj, ld->nlj);
92                 ld->lj[ld->nlj - 1] = i;
93             }
94             else if (std::strstr(names[i].name, "Coul") != nullptr)
95             {
96                 ld->nqq++;
97                 srenew(ld->qq, ld->nqq);
98                 ld->qq[ld->nqq - 1] = i;
99             }
100         }
101     }
102     printf("Using the following energy terms:\n");
103     printf("LJ:  ");
104     for (i = 0; (i < ld->nlj); i++)
105     {
106         printf("  %12s", names[ld->lj[i]].name);
107     }
108     printf("\nCoul:");
109     for (i = 0; (i < ld->nqq); i++)
110     {
111         printf("  %12s", names[ld->qq[i]].name);
112     }
113     printf("\n");
114
115     return ld;
116 }
117
118 static real calc_lie(t_liedata* ld, t_energy ee[], real lie_lj, real lie_qq, real fac_lj, real fac_qq)
119 {
120     int  i;
121     real lj_tot, qq_tot;
122
123     lj_tot = 0;
124     for (i = 0; (i < ld->nlj); i++)
125     {
126         lj_tot += ee[ld->lj[i]].e;
127     }
128     qq_tot = 0;
129     for (i = 0; (i < ld->nqq); i++)
130     {
131         qq_tot += ee[ld->qq[i]].e;
132     }
133
134     /* And now the great LIE formula: */
135     return fac_lj * (lj_tot - lie_lj) + fac_qq * (qq_tot - lie_qq);
136 }
137
138 int gmx_lie(int argc, char* argv[])
139 {
140     const char* desc[] = {
141         "[THISMODULE] computes a free energy estimate based on an energy analysis",
142         "from nonbonded energies. One needs an energy file with the following components:",
143         "Coul-(A-B) LJ-SR (A-B) etc.[PAR]",
144         "To utilize [TT]g_lie[tt] correctly, two simulations are required: one with the",
145         "molecule of interest bound to its receptor and one with the molecule in water.",
146         "Both need to utilize [TT]energygrps[tt] such that Coul-SR(A-B), LJ-SR(A-B), etc. terms",
147         "are written to the [REF].edr[ref] file. Values from the molecule-in-water simulation",
148         "are necessary for supplying suitable values for -Elj and -Eqq."
149     };
150     static real        lie_lj = 0, lie_qq = 0, fac_lj = 0.181, fac_qq = 0.5;
151     static const char* ligand = "none";
152     t_pargs            pa[]   = {
153         { "-Elj",
154           FALSE,
155           etREAL,
156           { &lie_lj },
157           "Lennard-Jones interaction between ligand and solvent" },
158         { "-Eqq", FALSE, etREAL, { &lie_qq }, "Coulomb interaction between ligand and solvent" },
159         { "-Clj",
160           FALSE,
161           etREAL,
162           { &fac_lj },
163           "Factor in the LIE equation for Lennard-Jones component of energy" },
164         { "-Cqq",
165           FALSE,
166           etREAL,
167           { &fac_qq },
168           "Factor in the LIE equation for Coulomb component of energy" },
169         { "-ligand", FALSE, etSTR, { &ligand }, "Name of the ligand in the energy file" }
170     };
171 #define NPA asize(pa)
172
173     FILE*             out;
174     int               nre, nframes = 0, ct = 0;
175     ener_file_t       fp;
176     t_liedata*        ld;
177     gmx_enxnm_t*      enm = nullptr;
178     t_enxframe*       fr;
179     real              lie;
180     double            lieaver = 0, lieav2 = 0;
181     gmx_output_env_t* oenv;
182
183     t_filenm fnm[] = { { efEDR, "-f", "ener", ffREAD }, { efXVG, "-o", "lie", ffWRITE } };
184 #define NFILE asize(fnm)
185
186     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME, NFILE, fnm, NPA, pa,
187                            asize(desc), desc, 0, nullptr, &oenv))
188     {
189         return 0;
190     }
191
192     fp = open_enx(ftp2fn(efEDR, NFILE, fnm), "r");
193     do_enxnms(fp, &nre, &enm);
194
195     ld = analyze_names(nre, enm, ligand);
196     snew(fr, 1);
197
198     out = xvgropen(ftp2fn(efXVG, NFILE, fnm), "LIE free energy estimate", "Time (ps)",
199                    "DGbind (kJ/mol)", oenv);
200     while (do_enx(fp, fr))
201     {
202         ct = check_times(fr->t);
203         if (ct == 0)
204         {
205             lie = calc_lie(ld, fr->ener, lie_lj, lie_qq, fac_lj, fac_qq);
206             lieaver += lie;
207             lieav2 += lie * lie;
208             nframes++;
209             fprintf(out, "%10g  %10g\n", fr->t, lie);
210         }
211     }
212     close_enx(fp);
213     xvgrclose(out);
214     fprintf(stderr, "\n");
215
216     if (nframes > 0)
217     {
218         printf("DGbind = %.3f (%.3f)\n", lieaver / nframes,
219                std::sqrt(lieav2 / nframes - gmx::square(lieaver / nframes)));
220     }
221
222     do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");
223
224     return 0;
225 }