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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
52 #include "gmx_fatal.h"
70 static t_liedata
*analyze_names(int nre
, gmx_enxnm_t
*names
, const char *ligand
)
76 /* Skip until we come to pressure */
77 for (i
= 0; (i
< F_NRE
); i
++)
79 if (strcmp(names
[i
].name
, interaction_function
[F_PRES
].longname
) == 0)
85 /* Now real analysis: find components of energies */
86 sprintf(self
, "%s-%s", ligand
, ligand
);
88 for (; (i
< nre
); i
++)
90 if ((strstr(names
[i
].name
, ligand
) != NULL
) &&
91 (strstr(names
[i
].name
, self
) == NULL
))
93 if (strstr(names
[i
].name
, "LJ") != NULL
)
96 srenew(ld
->lj
, ld
->nlj
);
97 ld
->lj
[ld
->nlj
-1] = i
;
99 else if (strstr(names
[i
].name
, "Coul") != NULL
)
102 srenew(ld
->qq
, ld
->nqq
);
103 ld
->qq
[ld
->nqq
-1] = i
;
107 printf("Using the following energy terms:\n");
109 for (i
= 0; (i
< ld
->nlj
); i
++)
111 printf(" %12s", names
[ld
->lj
[i
]].name
);
114 for (i
= 0; (i
< ld
->nqq
); i
++)
116 printf(" %12s", names
[ld
->qq
[i
]].name
);
123 real
calc_lie(t_liedata
*ld
, t_energy ee
[], real lie_lj
, real lie_qq
,
124 real fac_lj
, real fac_qq
)
130 for (i
= 0; (i
< ld
->nlj
); i
++)
132 lj_tot
+= ee
[ld
->lj
[i
]].e
;
135 for (i
= 0; (i
< ld
->nqq
); i
++)
137 qq_tot
+= ee
[ld
->qq
[i
]].e
;
140 /* And now the great LIE formula: */
141 return fac_lj
*(lj_tot
-lie_lj
)+fac_qq
*(qq_tot
-lie_qq
);
144 int gmx_lie(int argc
, char *argv
[])
146 const char *desc
[] = {
147 "[TT]g_lie[tt] computes a free energy estimate based on an energy analysis",
148 "from nonbonded energies. One needs an energy file with the following components:",
149 "Coul-(A-B) LJ-SR (A-B) etc.[PAR]",
150 "To utilize [TT]g_lie[tt] correctly, two simulations are required: one with the",
151 "molecule of interest bound to its receptor and one with the molecule in water.",
152 "Both need to utilize [TT]energygrps[tt] such that Coul-SR(A-B), LJ-SR(A-B), etc. terms",
153 "are written to the [TT].edr[tt] file. Values from the molecule-in-water simulation",
154 "are necessary for supplying suitable values for -Elj and -Eqq."
156 static real lie_lj
= 0, lie_qq
= 0, fac_lj
= 0.181, fac_qq
= 0.5;
157 static const char *ligand
= "none";
159 { "-Elj", FALSE
, etREAL
, {&lie_lj
},
160 "Lennard-Jones interaction between ligand and solvent" },
161 { "-Eqq", FALSE
, etREAL
, {&lie_qq
},
162 "Coulomb interaction between ligand and solvent" },
163 { "-Clj", FALSE
, etREAL
, {&fac_lj
},
164 "Factor in the LIE equation for Lennard-Jones component of energy" },
165 { "-Cqq", FALSE
, etREAL
, {&fac_qq
},
166 "Factor in the LIE equation for Coulomb component of energy" },
167 { "-ligand", FALSE
, etSTR
, {&ligand
},
168 "Name of the ligand in the energy file" }
170 #define NPA asize(pa)
173 int nre
, nframes
= 0, ct
= 0;
177 gmx_enxnm_t
*enm
= NULL
;
180 double lieaver
= 0, lieav2
= 0;
184 { efEDR
, "-f", "ener", ffREAD
},
185 { efXVG
, "-o", "lie", ffWRITE
}
187 #define NFILE asize(fnm)
189 CopyRight(stderr
, argv
[0]);
190 parse_common_args(&argc
, argv
, PCA_CAN_VIEW
| PCA_CAN_TIME
| PCA_BE_NICE
,
191 NFILE
, fnm
, NPA
, pa
, asize(desc
), desc
, 0, NULL
, &oenv
);
193 fp
= open_enx(ftp2fn(efEDR
, NFILE
, fnm
), "r");
194 do_enxnms(fp
, &nre
, &enm
);
196 ld
= analyze_names(nre
, enm
, ligand
);
199 out
= xvgropen(ftp2fn(efXVG
, NFILE
, fnm
), "LIE free energy estimate",
200 "Time (ps)", "DGbind (kJ/mol)", oenv
);
201 while (do_enx(fp
, fr
))
203 ct
= check_times(fr
->t
);
206 lie
= calc_lie(ld
, fr
->ener
, lie_lj
, lie_qq
, fac_lj
, fac_qq
);
210 fprintf(out
, "%10g %10g\n", fr
->t
, lie
);
215 fprintf(stderr
, "\n");
219 printf("DGbind = %.3f (%.3f)\n", lieaver
/nframes
,
220 sqrt(lieav2
/nframes
-sqr(lieaver
/nframes
)));
223 do_view(oenv
, ftp2fn(efXVG
, NFILE
, fnm
), "-nxy");