Fix problems with intel-mpi
[gromacs.git] / src / tools / gmx_lie.c
blob6b7117ab0d1d4371cb93fe7ca37eb2eea4769242
1 /*
3 * This source code is part of
5 * G R O M A C S
7 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * Green Red Orange Magenta Azure Cyan Skyblue
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
39 #include <stdio.h>
40 #include <stdlib.h>
41 #include <math.h>
42 #include <string.h>
44 #include "statutil.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 "copyrite.h"
52 #include "futil.h"
53 #include "statutil.h"
54 #include "txtdump.h"
55 #include "enxio.h"
56 #include "gstat.h"
57 #include "xvgr.h"
58 #include "gmx_ana.h"
61 typedef struct {
62 int nlj, nqq;
63 int *lj;
64 int *qq;
65 } t_liedata;
67 static t_liedata *analyze_names(int nre, gmx_enxnm_t *names, const char *ligand)
69 int i;
70 t_liedata *ld;
71 char self[256];
73 /* Skip until we come to pressure */
74 for (i = 0; (i < F_NRE); i++)
76 if (strcmp(names[i].name, interaction_function[F_PRES].longname) == 0)
78 break;
82 /* Now real analysis: find components of energies */
83 sprintf(self, "%s-%s", ligand, ligand);
84 snew(ld, 1);
85 for (; (i < nre); i++)
87 if ((strstr(names[i].name, ligand) != NULL) &&
88 (strstr(names[i].name, self) == NULL))
90 if (strstr(names[i].name, "LJ") != NULL)
92 ld->nlj++;
93 srenew(ld->lj, ld->nlj);
94 ld->lj[ld->nlj-1] = i;
96 else if (strstr(names[i].name, "Coul") != NULL)
98 ld->nqq++;
99 srenew(ld->qq, ld->nqq);
100 ld->qq[ld->nqq-1] = i;
104 printf("Using the following energy terms:\n");
105 printf("LJ: ");
106 for (i = 0; (i < ld->nlj); i++)
108 printf(" %12s", names[ld->lj[i]].name);
110 printf("\nCoul:");
111 for (i = 0; (i < ld->nqq); i++)
113 printf(" %12s", names[ld->qq[i]].name);
115 printf("\n");
117 return ld;
120 real calc_lie(t_liedata *ld, t_energy ee[], real lie_lj, real lie_qq,
121 real fac_lj, real fac_qq)
123 int i;
124 real lj_tot, qq_tot;
126 lj_tot = 0;
127 for (i = 0; (i < ld->nlj); i++)
129 lj_tot += ee[ld->lj[i]].e;
131 qq_tot = 0;
132 for (i = 0; (i < ld->nqq); i++)
134 qq_tot += ee[ld->qq[i]].e;
137 /* And now the great LIE formula: */
138 return fac_lj*(lj_tot-lie_lj)+fac_qq*(qq_tot-lie_qq);
141 int gmx_lie(int argc, char *argv[])
143 const char *desc[] = {
144 "[TT]g_lie[tt] computes a free energy estimate based on an energy analysis",
145 "from. One needs an energy file with the following components:",
146 "Coul (A-B) LJ-SR (A-B) etc."
148 static real lie_lj = 0, lie_qq = 0, fac_lj = 0.181, fac_qq = 0.5;
149 static const char *ligand = "none";
150 t_pargs pa[] = {
151 { "-Elj", FALSE, etREAL, {&lie_lj},
152 "Lennard-Jones interaction between ligand and solvent" },
153 { "-Eqq", FALSE, etREAL, {&lie_qq},
154 "Coulomb interaction between ligand and solvent" },
155 { "-Clj", FALSE, etREAL, {&fac_lj},
156 "Factor in the LIE equation for Lennard-Jones component of energy" },
157 { "-Cqq", FALSE, etREAL, {&fac_qq},
158 "Factor in the LIE equation for Coulomb component of energy" },
159 { "-ligand", FALSE, etSTR, {&ligand},
160 "Name of the ligand in the energy file" }
162 #define NPA asize(pa)
164 FILE *out;
165 int nre, nframes = 0, ct = 0;
166 ener_file_t fp;
167 gmx_bool bCont;
168 t_liedata *ld;
169 gmx_enxnm_t *enm = NULL;
170 t_enxframe *fr;
171 real lie;
172 double lieaver = 0, lieav2 = 0;
173 output_env_t oenv;
175 t_filenm fnm[] = {
176 { efEDR, "-f", "ener", ffREAD },
177 { efXVG, "-o", "lie", ffWRITE }
179 #define NFILE asize(fnm)
181 parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
182 NFILE, fnm, NPA, pa, asize(desc), desc, 0, NULL, &oenv);
184 fp = open_enx(ftp2fn(efEDR, NFILE, fnm), "r");
185 do_enxnms(fp, &nre, &enm);
187 ld = analyze_names(nre, enm, ligand);
188 snew(fr, 1);
189 out = xvgropen(ftp2fn(efXVG, NFILE, fnm), "LIE free energy estimate",
190 "Time (ps)", "DGbind (kJ/mol)", oenv);
191 while (do_enx(fp, fr))
193 ct = check_times(fr->t);
194 if (ct == 0)
196 lie = calc_lie(ld, fr->ener, lie_lj, lie_qq, fac_lj, fac_qq);
197 lieaver += lie;
198 lieav2 += lie*lie;
199 nframes++;
200 fprintf(out, "%10g %10g\n", fr->t, lie);
203 close_enx(fp);
204 ffclose(out);
205 fprintf(stderr, "\n");
207 if (nframes > 0)
209 printf("DGbind = %.3f (%.3f)\n", lieaver/nframes,
210 sqrt(lieav2/nframes-sqr(lieaver/nframes)));
213 do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");
215 thanx(stderr);
217 return 0;