Remove sysstuff.h
[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 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <stdio.h>
42 #include <stdlib.h>
43 #include <math.h>
44 #include <string.h>
45
46 #include "gromacs/commandline/pargs.h"
47 #include "typedefs.h"
48 #include "gromacs/utility/smalloc.h"
49 #include "macros.h"
50 #include "vec.h"
51 #include "gromacs/utility/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 typedef struct {
60     int  nlj, nqq;
61     int *lj;
62     int *qq;
63 } t_liedata;
64
65 static t_liedata *analyze_names(int nre, gmx_enxnm_t *names, const char *ligand)
66 {
67     int        i;
68     t_liedata *ld;
69     char       self[256];
70
71     /* Skip until we come to pressure */
72     for (i = 0; (i < F_NRE); i++)
73     {
74         if (strcmp(names[i].name, interaction_function[F_PRES].longname) == 0)
75         {
76             break;
77         }
78     }
79
80     /* Now real analysis: find components of energies */
81     sprintf(self, "%s-%s", ligand, ligand);
82     snew(ld, 1);
83     for (; (i < nre); i++)
84     {
85         if ((strstr(names[i].name, ligand) != NULL) &&
86             (strstr(names[i].name, self) == NULL))
87         {
88             if (strstr(names[i].name, "LJ") != NULL)
89             {
90                 ld->nlj++;
91                 srenew(ld->lj, ld->nlj);
92                 ld->lj[ld->nlj-1] = i;
93             }
94             else if (strstr(names[i].name, "Coul") != NULL)
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 real calc_lie(t_liedata *ld, t_energy ee[], real lie_lj, real lie_qq,
119               real fac_lj, real fac_qq)
120 {
121     int  i;
122     real lj_tot, qq_tot;
123
124     lj_tot = 0;
125     for (i = 0; (i < ld->nlj); i++)
126     {
127         lj_tot += ee[ld->lj[i]].e;
128     }
129     qq_tot = 0;
130     for (i = 0; (i < ld->nqq); i++)
131     {
132         qq_tot += ee[ld->qq[i]].e;
133     }
134
135     /* And now the great LIE formula: */
136     return fac_lj*(lj_tot-lie_lj)+fac_qq*(qq_tot-lie_qq);
137 }
138
139 int gmx_lie(int argc, char *argv[])
140 {
141     const char        *desc[] = {
142         "[THISMODULE] computes a free energy estimate based on an energy analysis",
143         "from. One needs an energy file with the following components:",
144         "Coul (A-B) LJ-SR (A-B) etc."
145     };
146     static real        lie_lj = 0, lie_qq = 0, fac_lj = 0.181, fac_qq = 0.5;
147     static const char *ligand = "none";
148     t_pargs            pa[]   = {
149         { "-Elj",  FALSE, etREAL, {&lie_lj},
150           "Lennard-Jones interaction between ligand and solvent" },
151         { "-Eqq",  FALSE, etREAL, {&lie_qq},
152           "Coulomb interaction between ligand and solvent" },
153         { "-Clj",  FALSE, etREAL, {&fac_lj},
154           "Factor in the LIE equation for Lennard-Jones component of energy" },
155         { "-Cqq",  FALSE, etREAL, {&fac_qq},
156           "Factor in the LIE equation for Coulomb component of energy" },
157         { "-ligand",  FALSE, etSTR, {&ligand},
158           "Name of the ligand in the energy file" }
159     };
160 #define NPA asize(pa)
161
162     FILE         *out;
163     int           nre, nframes = 0, ct = 0;
164     ener_file_t   fp;
165     gmx_bool      bCont;
166     t_liedata    *ld;
167     gmx_enxnm_t  *enm = NULL;
168     t_enxframe   *fr;
169     real          lie;
170     double        lieaver = 0, lieav2 = 0;
171     output_env_t  oenv;
172
173     t_filenm      fnm[] = {
174         { efEDR, "-f",    "ener",     ffREAD   },
175         { efXVG, "-o",    "lie",      ffWRITE  }
176     };
177 #define NFILE asize(fnm)
178
179     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
180                            NFILE, fnm, NPA, pa, asize(desc), desc, 0, NULL, &oenv))
181     {
182         return 0;
183     }
184
185     fp = open_enx(ftp2fn(efEDR, NFILE, fnm), "r");
186     do_enxnms(fp, &nre, &enm);
187
188     ld = analyze_names(nre, enm, ligand);
189     snew(fr, 1);
190     out = xvgropen(ftp2fn(efXVG, NFILE, fnm), "LIE free energy estimate",
191                    "Time (ps)", "DGbind (kJ/mol)", oenv);
192     while (do_enx(fp, fr))
193     {
194         ct = check_times(fr->t);
195         if (ct == 0)
196         {
197             lie      = calc_lie(ld, fr->ener, lie_lj, lie_qq, fac_lj, fac_qq);
198             lieaver += lie;
199             lieav2  += lie*lie;
200             nframes++;
201             fprintf(out, "%10g  %10g\n", fr->t, lie);
202         }
203     }
204     close_enx(fp);
205     gmx_ffclose(out);
206     fprintf(stderr, "\n");
207
208     if (nframes > 0)
209     {
210         printf("DGbind = %.3f (%.3f)\n", lieaver/nframes,
211                sqrt(lieav2/nframes-sqr(lieaver/nframes)));
212     }
213
214     do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");
215
216     return 0;
217 }