Renamed entropy.* to thermochemistry.*
[gromacs.git] / src / gromacs / gmxana / gmx_anadock.cpp
blobd9f1fcede9e01d0aa26183898aaef3ba7a110035
1 /*
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) 2012,2013,2014,2015,2017, 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.
37 #include "gmxpre.h"
39 #include <cmath>
40 #include <cstdlib>
41 #include <cstring>
43 #include "gromacs/commandline/pargs.h"
44 #include "gromacs/fileio/confio.h"
45 #include "gromacs/fileio/xvgr.h"
46 #include "gromacs/gmxana/gmx_ana.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/statistics/statistics.h"
49 #include "gromacs/topology/topology.h"
50 #include "gromacs/utility/arraysize.h"
51 #include "gromacs/utility/cstringutil.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/futil.h"
54 #include "gromacs/utility/pleasecite.h"
55 #include "gromacs/utility/smalloc.h"
57 static const char *etitles[] = { "E-docked", "Free Energy" };
59 typedef struct {
60 real edocked, efree;
61 int index, cluster_id;
62 t_atoms atoms;
63 rvec *x;
64 int ePBC;
65 matrix box;
66 } t_pdbfile;
68 static t_pdbfile *read_pdbf(const char *fn)
70 t_pdbfile *pdbf;
71 double e;
72 FILE *fp;
74 snew(pdbf, 1);
75 t_topology top;
76 read_tps_conf(fn, &top, &pdbf->ePBC, &pdbf->x, nullptr, pdbf->box, FALSE);
77 pdbf->atoms = top.atoms;
78 fp = gmx_ffopen(fn, "r");
79 char buf[256], *ptr;
80 while ((ptr = fgets2(buf, 255, fp)) != nullptr)
82 if (std::strstr(buf, "Intermolecular") != nullptr)
84 ptr = std::strchr(buf, '=');
85 sscanf(ptr+1, "%lf", &e);
86 pdbf->edocked = e;
88 else if (std::strstr(buf, "Estimated Free") != nullptr)
90 ptr = std::strchr(buf, '=');
91 sscanf(ptr+1, "%lf", &e);
92 pdbf->efree = e;
95 gmx_ffclose(fp);
97 return pdbf;
100 static t_pdbfile **read_em_all(const char *fn, int *npdbf)
102 t_pdbfile **pdbf = nullptr;
103 int i, maxpdbf;
104 char buf[256], name[256];
105 gmx_bool bExist;
107 std::strcpy(buf, fn);
108 buf[std::strlen(buf)-4] = '\0';
109 maxpdbf = 100;
110 snew(pdbf, maxpdbf);
111 i = 0;
114 sprintf(name, "%s_%d.pdb", buf, i+1);
115 if ((bExist = gmx_fexist(name)) == TRUE)
117 pdbf[i] = read_pdbf(name);
118 pdbf[i]->index = i+1;
119 i++;
120 if (i >= maxpdbf)
122 maxpdbf += 100;
123 srenew(pdbf, maxpdbf);
127 while (bExist);
129 *npdbf = i;
131 printf("Found %d pdb files\n", i);
133 return pdbf;
136 static gmx_bool bFreeSort = FALSE;
138 static int pdbf_comp(const void *a, const void *b)
140 t_pdbfile *pa, *pb;
141 real x;
142 int dc;
144 pa = *(t_pdbfile **)a;
145 pb = *(t_pdbfile **)b;
147 dc = pa->cluster_id - pb->cluster_id;
149 if (dc == 0)
151 if (bFreeSort)
153 x = pa->efree - pb->efree;
155 else
157 x = pa->edocked - pb->edocked;
160 if (x < 0)
162 return -1;
164 else if (x > 0)
166 return 1;
168 else
170 return 0;
173 else
175 return dc;
179 static void analyse_em_all(int npdb, t_pdbfile *pdbf[], const char *edocked,
180 const char *efree, const gmx_output_env_t *oenv)
182 FILE *fp;
183 int i;
185 for (bFreeSort = FALSE; (bFreeSort <= TRUE); bFreeSort++)
187 qsort(pdbf, npdb, sizeof(pdbf[0]), pdbf_comp);
188 fp = xvgropen(bFreeSort ? efree : edocked,
189 etitles[bFreeSort], "()", "E (kJ/mol)", oenv);
190 for (i = 0; (i < npdb); i++)
192 fprintf(fp, "%12lf\n", bFreeSort ? pdbf[i]->efree : pdbf[i]->edocked);
194 xvgrclose(fp);
198 static void clust_stat(FILE *fp, int start, int end, t_pdbfile *pdbf[])
200 int i;
201 gmx_stats_t ed, ef;
202 real aver, sigma;
204 ed = gmx_stats_init();
205 ef = gmx_stats_init();
206 for (i = start; (i < end); i++)
208 gmx_stats_add_point(ed, i-start, pdbf[i]->edocked, 0, 0);
209 gmx_stats_add_point(ef, i-start, pdbf[i]->efree, 0, 0);
211 gmx_stats_get_ase(ed, &aver, &sigma, nullptr);
212 fprintf(fp, " <%12s> = %8.3f (+/- %6.3f)\n", etitles[FALSE], aver, sigma);
213 gmx_stats_get_ase(ef, &aver, &sigma, nullptr);
214 fprintf(fp, " <%12s> = %8.3f (+/- %6.3f)\n", etitles[TRUE], aver, sigma);
215 gmx_stats_free(ed);
216 gmx_stats_free(ef);
219 static real rmsd_dist(t_pdbfile *pa, t_pdbfile *pb, gmx_bool bRMSD)
221 int i;
222 real rmsd;
223 rvec acm, bcm, dx;
225 if (bRMSD)
227 rmsd = 0;
228 for (i = 0; (i < pa->atoms.nr); i++)
230 rvec_sub(pa->x[i], pb->x[i], dx);
231 rmsd += iprod(dx, dx);
233 rmsd = std::sqrt(rmsd/pa->atoms.nr);
235 else
237 clear_rvec(acm);
238 clear_rvec(bcm);
239 for (i = 0; (i < pa->atoms.nr); i++)
241 rvec_inc(acm, pa->x[i]);
242 rvec_inc(bcm, pb->x[i]);
244 rvec_sub(acm, bcm, dx);
245 for (i = 0; (i < DIM); i++)
247 dx[i] /= pa->atoms.nr;
249 rmsd = norm(dx);
251 return rmsd;
254 static void line(FILE *fp)
256 fprintf(fp, " - - - - - - - - - -\n");
259 static void cluster_em_all(FILE *fp, int npdb, t_pdbfile *pdbf[],
260 gmx_bool bFree, gmx_bool bRMSD, real cutoff)
262 int i, j, k;
263 int *cndx, ncluster;
264 real rmsd;
266 bFreeSort = bFree;
267 qsort(pdbf, npdb, sizeof(pdbf[0]), pdbf_comp);
269 fprintf(fp, "Statistics over all structures:\n");
270 clust_stat(fp, 0, npdb, pdbf);
271 line(fp);
273 /* Index to first structure in a cluster */
274 snew(cndx, npdb);
275 ncluster = 1;
277 for (i = 1; (i < npdb); i++)
279 for (j = 0; (j < ncluster); j++)
281 rmsd = rmsd_dist(pdbf[cndx[j]], pdbf[i], bRMSD);
282 if (rmsd <= cutoff)
284 /* Structure i is in cluster j */
285 pdbf[i]->cluster_id = pdbf[cndx[j]]->cluster_id;
286 break;
289 if (j == ncluster)
291 /* New cluster! Cool! */
292 cndx[ncluster] = i;
293 pdbf[i]->cluster_id = ncluster;
294 ncluster++;
297 cndx[ncluster] = npdb;
298 qsort(pdbf, npdb, sizeof(pdbf[0]), pdbf_comp);
300 j = 0;
301 cndx[j++] = 0;
302 for (i = 1; (i < npdb); i++)
304 if (pdbf[i]->cluster_id != pdbf[i-1]->cluster_id)
306 cndx[j++] = i;
309 cndx[j] = npdb;
310 if (j != ncluster)
312 gmx_fatal(FARGS, "Consistency error: j = %d, ncluster = %d", j, ncluster);
315 fprintf(fp, "I found %d clusters based on %s and %s with a %.3f nm cut-off\n",
316 ncluster, etitles[bFree], bRMSD ? "RMSD" : "distance", cutoff);
317 line(fp);
318 for (j = 0; (j < ncluster); j++)
320 fprintf(fp, "Cluster: %3d %s: %10.5f kJ/mol %3d elements\n",
321 j, etitles[bFree],
322 bFree ? pdbf[cndx[j]]->efree : pdbf[cndx[j]]->edocked,
323 cndx[j+1]-cndx[j]);
324 clust_stat(fp, cndx[j], cndx[j+1], pdbf);
325 for (k = cndx[j]; (k < cndx[j+1]); k++)
327 fprintf(fp, " %3d", pdbf[k]->index);
329 fprintf(fp, "\n");
330 line(fp);
332 sfree(cndx);
335 int gmx_anadock(int argc, char *argv[])
337 const char *desc[] = {
338 "[THISMODULE] analyses the results of an Autodock run and clusters the",
339 "structures together, based on distance or RMSD. The docked energy",
340 "and free energy estimates are analysed, and for each cluster the",
341 "energy statistics are printed.[PAR]",
342 "An alternative approach to this is to cluster the structures first",
343 "using [gmx-cluster] and then sort the clusters on either lowest",
344 "energy or average energy."
346 t_filenm fnm[] = {
347 { efPDB, "-f", nullptr, ffREAD },
348 { efXVG, "-od", "edocked", ffWRITE },
349 { efXVG, "-of", "efree", ffWRITE },
350 { efLOG, "-g", "anadock", ffWRITE }
352 gmx_output_env_t *oenv;
353 #define NFILE asize(fnm)
354 static gmx_bool bFree = FALSE, bRMS = TRUE;
355 static real cutoff = 0.2;
356 t_pargs pa[] = {
357 { "-free", FALSE, etBOOL, {&bFree},
358 "Use Free energy estimate from autodock for sorting the classes" },
359 { "-rms", FALSE, etBOOL, {&bRMS},
360 "Cluster on RMS or distance" },
361 { "-cutoff", FALSE, etREAL, {&cutoff},
362 "Maximum RMSD/distance for belonging to the same cluster" }
364 #define NPA asize(pa)
366 FILE *fp;
367 t_pdbfile **pdbf = nullptr;
368 int npdbf;
370 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, NPA, pa, asize(desc), desc, 0,
371 nullptr, &oenv))
373 return 0;
376 fp = gmx_ffopen(opt2fn("-g", NFILE, fnm), "w");
377 please_cite(stdout, "Hetenyi2002b");
378 please_cite(fp, "Hetenyi2002b");
380 pdbf = read_em_all(opt2fn("-f", NFILE, fnm), &npdbf);
382 analyse_em_all(npdbf, pdbf, opt2fn("-od", NFILE, fnm), opt2fn("-of", NFILE, fnm),
383 oenv);
385 cluster_em_all(fp, npdbf, pdbf, bFree, bRMS, cutoff);
387 gmx_ffclose(fp);
389 return 0;