Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_lie.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
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.
14
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.
19  *
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.
26  *
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.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Green Red Orange Magenta Azure Cyan Skyblue
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <stdio.h>
40 #include <stdlib.h>
41 #include <math.h>
42 #include <string.h>
43
44 #include "gromacs/commandline/pargs.h"
45 #include "sysstuff.h"
46 #include "typedefs.h"
47 #include "smalloc.h"
48 #include "macros.h"
49 #include "gmx_fatal.h"
50 #include "vec.h"
51 #include "gromacs/fileio/futil.h"
52 #include "txtdump.h"
53 #include "gromacs/fileio/enxio.h"
54 #include "gstat.h"
55 #include "xvgr.h"
56 #include "gmx_ana.h"
57 #include "gromacs/fileio/trxio.h"
58
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. One needs an energy file with the following components:",
145         "Coul (A-B) LJ-SR (A-B) etc."
146     };
147     static real        lie_lj = 0, lie_qq = 0, fac_lj = 0.181, fac_qq = 0.5;
148     static const char *ligand = "none";
149     t_pargs            pa[]   = {
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" }
160     };
161 #define NPA asize(pa)
162
163     FILE         *out;
164     int           nre, nframes = 0, ct = 0;
165     ener_file_t   fp;
166     gmx_bool      bCont;
167     t_liedata    *ld;
168     gmx_enxnm_t  *enm = NULL;
169     t_enxframe   *fr;
170     real          lie;
171     double        lieaver = 0, lieav2 = 0;
172     output_env_t  oenv;
173
174     t_filenm      fnm[] = {
175         { efEDR, "-f",    "ener",     ffREAD   },
176         { efXVG, "-o",    "lie",      ffWRITE  }
177     };
178 #define NFILE asize(fnm)
179
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))
182     {
183         return 0;
184     }
185
186     fp = open_enx(ftp2fn(efEDR, NFILE, fnm), "r");
187     do_enxnms(fp, &nre, &enm);
188
189     ld = analyze_names(nre, enm, ligand);
190     snew(fr, 1);
191     out = xvgropen(ftp2fn(efXVG, NFILE, fnm), "LIE free energy estimate",
192                    "Time (ps)", "DGbind (kJ/mol)", oenv);
193     while (do_enx(fp, fr))
194     {
195         ct = check_times(fr->t);
196         if (ct == 0)
197         {
198             lie      = calc_lie(ld, fr->ener, lie_lj, lie_qq, fac_lj, fac_qq);
199             lieaver += lie;
200             lieav2  += lie*lie;
201             nframes++;
202             fprintf(out, "%10g  %10g\n", fr->t, lie);
203         }
204     }
205     close_enx(fp);
206     ffclose(out);
207     fprintf(stderr, "\n");
208
209     if (nframes > 0)
210     {
211         printf("DGbind = %.3f (%.3f)\n", lieaver/nframes,
212                sqrt(lieav2/nframes-sqr(lieaver/nframes)));
213     }
214
215     do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");
216
217     return 0;
218 }