290f3d2b8497f72e6e6013cc4b3e27d9b5a74cc3
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_lie.c
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, 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 "config.h"
40
41 #include <math.h>
42 #include <stdio.h>
43 #include <stdlib.h>
44 #include <string.h>
45
46 #include "gromacs/commandline/pargs.h"
47 #include "gromacs/legacyheaders/typedefs.h"
48 #include "gromacs/utility/smalloc.h"
49 #include "gromacs/legacyheaders/macros.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/utility/futil.h"
52 #include "gromacs/legacyheaders/txtdump.h"
53 #include "gromacs/fileio/enxio.h"
54 #include "gstat.h"
55 #include "gromacs/fileio/xvgr.h"
56 #include "gromacs/legacyheaders/viewit.h"
57 #include "gmx_ana.h"
58 #include "gromacs/fileio/trxio.h"
59
60 typedef struct {
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 (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 ((strstr(names[i].name, ligand) != NULL) &&
87             (strstr(names[i].name, self) == NULL))
88         {
89             if (strstr(names[i].name, "LJ") != NULL)
90             {
91                 ld->nlj++;
92                 srenew(ld->lj, ld->nlj);
93                 ld->lj[ld->nlj-1] = i;
94             }
95             else if (strstr(names[i].name, "Coul") != NULL)
96             {
97                 ld->nqq++;
98                 srenew(ld->qq, ld->nqq);
99                 ld->qq[ld->nqq-1] = i;
100             }
101         }
102     }
103     printf("Using the following energy terms:\n");
104     printf("LJ:  ");
105     for (i = 0; (i < ld->nlj); i++)
106     {
107         printf("  %12s", names[ld->lj[i]].name);
108     }
109     printf("\nCoul:");
110     for (i = 0; (i < ld->nqq); i++)
111     {
112         printf("  %12s", names[ld->qq[i]].name);
113     }
114     printf("\n");
115
116     return ld;
117 }
118
119 real calc_lie(t_liedata *ld, t_energy ee[], real lie_lj, real lie_qq,
120               real fac_lj, real fac_qq)
121 {
122     int  i;
123     real lj_tot, qq_tot;
124
125     lj_tot = 0;
126     for (i = 0; (i < ld->nlj); i++)
127     {
128         lj_tot += ee[ld->lj[i]].e;
129     }
130     qq_tot = 0;
131     for (i = 0; (i < ld->nqq); i++)
132     {
133         qq_tot += ee[ld->qq[i]].e;
134     }
135
136     /* And now the great LIE formula: */
137     return fac_lj*(lj_tot-lie_lj)+fac_qq*(qq_tot-lie_qq);
138 }
139
140 int gmx_lie(int argc, char *argv[])
141 {
142     const char        *desc[] = {
143         "[THISMODULE] computes a free energy estimate based on an energy analysis",
144         "from nonbonded energies. One needs an energy file with the following components:",
145         "Coul-(A-B) LJ-SR (A-B) etc.[PAR]",
146         "To utilize [TT]g_lie[tt] correctly, two simulations are required: one with the",
147         "molecule of interest bound to its receptor and one with the molecule in water.",
148         "Both need to utilize [TT]energygrps[tt] such that Coul-SR(A-B), LJ-SR(A-B), etc. terms",
149         "are written to the [TT].edr[tt] file. Values from the molecule-in-water simulation",
150         "are necessary for supplying suitable values for -Elj and -Eqq."
151     };
152     static real        lie_lj = 0, lie_qq = 0, fac_lj = 0.181, fac_qq = 0.5;
153     static const char *ligand = "none";
154     t_pargs            pa[]   = {
155         { "-Elj",  FALSE, etREAL, {&lie_lj},
156           "Lennard-Jones interaction between ligand and solvent" },
157         { "-Eqq",  FALSE, etREAL, {&lie_qq},
158           "Coulomb interaction between ligand and solvent" },
159         { "-Clj",  FALSE, etREAL, {&fac_lj},
160           "Factor in the LIE equation for Lennard-Jones component of energy" },
161         { "-Cqq",  FALSE, etREAL, {&fac_qq},
162           "Factor in the LIE equation for Coulomb component of energy" },
163         { "-ligand",  FALSE, etSTR, {&ligand},
164           "Name of the ligand in the energy file" }
165     };
166 #define NPA asize(pa)
167
168     FILE         *out;
169     int           nre, nframes = 0, ct = 0;
170     ener_file_t   fp;
171     gmx_bool      bCont;
172     t_liedata    *ld;
173     gmx_enxnm_t  *enm = NULL;
174     t_enxframe   *fr;
175     real          lie;
176     double        lieaver = 0, lieav2 = 0;
177     output_env_t  oenv;
178
179     t_filenm      fnm[] = {
180         { efEDR, "-f",    "ener",     ffREAD   },
181         { efXVG, "-o",    "lie",      ffWRITE  }
182     };
183 #define NFILE asize(fnm)
184
185     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME,
186                            NFILE, fnm, NPA, pa, asize(desc), desc, 0, NULL, &oenv))
187     {
188         return 0;
189     }
190
191     fp = open_enx(ftp2fn(efEDR, NFILE, fnm), "r");
192     do_enxnms(fp, &nre, &enm);
193
194     ld = analyze_names(nre, enm, ligand);
195     snew(fr, 1);
196
197     out = xvgropen(ftp2fn(efXVG, NFILE, fnm), "LIE free energy estimate",
198                    "Time (ps)", "DGbind (kJ/mol)", oenv);
199     while (do_enx(fp, fr))
200     {
201         ct = check_times(fr->t);
202         if (ct == 0)
203         {
204             lie      = calc_lie(ld, fr->ener, lie_lj, lie_qq, fac_lj, fac_qq);
205             lieaver += lie;
206             lieav2  += lie*lie;
207             nframes++;
208             fprintf(out, "%10g  %10g\n", fr->t, lie);
209         }
210     }
211     close_enx(fp);
212     gmx_ffclose(out);
213     fprintf(stderr, "\n");
214
215     if (nframes > 0)
216     {
217         printf("DGbind = %.3f (%.3f)\n", lieaver/nframes,
218                sqrt(lieav2/nframes-sqr(lieaver/nframes)));
219     }
220
221     do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");
222
223     return 0;
224 }