2 * This file is part of the GROMACS molecular simulation package.
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, 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.
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.
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.
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.
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.
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.
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/utility/arraysize.h"
55 #include "gromacs/utility/futil.h"
56 #include "gromacs/utility/smalloc.h"
64 static t_liedata
*analyze_names(int nre
, gmx_enxnm_t
*names
, const char *ligand
)
70 /* Skip until we come to pressure */
71 for (i
= 0; (i
< F_NRE
); i
++)
73 if (std::strcmp(names
[i
].name
, interaction_function
[F_PRES
].longname
) == 0)
79 /* Now real analysis: find components of energies */
80 sprintf(self
, "%s-%s", ligand
, ligand
);
82 for (; (i
< nre
); i
++)
84 if ((std::strstr(names
[i
].name
, ligand
) != NULL
) &&
85 (std::strstr(names
[i
].name
, self
) == NULL
))
87 if (std::strstr(names
[i
].name
, "LJ") != NULL
)
90 srenew(ld
->lj
, ld
->nlj
);
91 ld
->lj
[ld
->nlj
-1] = i
;
93 else if (std::strstr(names
[i
].name
, "Coul") != NULL
)
96 srenew(ld
->qq
, ld
->nqq
);
97 ld
->qq
[ld
->nqq
-1] = i
;
101 printf("Using the following energy terms:\n");
103 for (i
= 0; (i
< ld
->nlj
); i
++)
105 printf(" %12s", names
[ld
->lj
[i
]].name
);
108 for (i
= 0; (i
< ld
->nqq
); i
++)
110 printf(" %12s", names
[ld
->qq
[i
]].name
);
117 real
calc_lie(t_liedata
*ld
, t_energy ee
[], real lie_lj
, real lie_qq
,
118 real fac_lj
, real fac_qq
)
124 for (i
= 0; (i
< ld
->nlj
); i
++)
126 lj_tot
+= ee
[ld
->lj
[i
]].e
;
129 for (i
= 0; (i
< ld
->nqq
); i
++)
131 qq_tot
+= ee
[ld
->qq
[i
]].e
;
134 /* And now the great LIE formula: */
135 return fac_lj
*(lj_tot
-lie_lj
)+fac_qq
*(qq_tot
-lie_qq
);
138 int gmx_lie(int argc
, char *argv
[])
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."
150 static real lie_lj
= 0, lie_qq
= 0, fac_lj
= 0.181, fac_qq
= 0.5;
151 static const char *ligand
= "none";
153 { "-Elj", FALSE
, etREAL
, {&lie_lj
},
154 "Lennard-Jones interaction between ligand and solvent" },
155 { "-Eqq", FALSE
, etREAL
, {&lie_qq
},
156 "Coulomb interaction between ligand and solvent" },
157 { "-Clj", FALSE
, etREAL
, {&fac_lj
},
158 "Factor in the LIE equation for Lennard-Jones component of energy" },
159 { "-Cqq", FALSE
, etREAL
, {&fac_qq
},
160 "Factor in the LIE equation for Coulomb component of energy" },
161 { "-ligand", FALSE
, etSTR
, {&ligand
},
162 "Name of the ligand in the energy file" }
164 #define NPA asize(pa)
167 int nre
, nframes
= 0, ct
= 0;
170 gmx_enxnm_t
*enm
= NULL
;
173 double lieaver
= 0, lieav2
= 0;
174 gmx_output_env_t
*oenv
;
177 { efEDR
, "-f", "ener", ffREAD
},
178 { efXVG
, "-o", "lie", ffWRITE
}
180 #define NFILE asize(fnm)
182 if (!parse_common_args(&argc
, argv
, PCA_CAN_VIEW
| PCA_CAN_TIME
,
183 NFILE
, fnm
, NPA
, pa
, asize(desc
), desc
, 0, NULL
, &oenv
))
188 fp
= open_enx(ftp2fn(efEDR
, NFILE
, fnm
), "r");
189 do_enxnms(fp
, &nre
, &enm
);
191 ld
= analyze_names(nre
, enm
, ligand
);
194 out
= xvgropen(ftp2fn(efXVG
, NFILE
, fnm
), "LIE free energy estimate",
195 "Time (ps)", "DGbind (kJ/mol)", oenv
);
196 while (do_enx(fp
, fr
))
198 ct
= check_times(fr
->t
);
201 lie
= calc_lie(ld
, fr
->ener
, lie_lj
, lie_qq
, fac_lj
, fac_qq
);
205 fprintf(out
, "%10g %10g\n", fr
->t
, lie
);
210 fprintf(stderr
, "\n");
214 printf("DGbind = %.3f (%.3f)\n", lieaver
/nframes
,
215 std::sqrt(lieav2
/nframes
-gmx::square(lieaver
/nframes
)));
218 do_view(oenv
, ftp2fn(efXVG
, NFILE
, fnm
), "-nxy");