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, 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/fileio/confio.h"
45 #include "gromacs/legacyheaders/copyrite.h"
46 #include "gromacs/utility/futil.h"
47 #include "gromacs/utility/fatalerror.h"
48 #include "gromacs/utility/smalloc.h"
49 #include "gromacs/utility/cstringutil.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/statistics/statistics.h"
52 #include "gromacs/commandline/pargs.h"
53 #include "gromacs/legacyheaders/typedefs.h"
54 #include "gromacs/fileio/xvgr.h"
55 #include "gromacs/legacyheaders/macros.h"
57 static const char *etitles[] = { "E-docked", "Free Energy" };
61 int index, cluster_id;
68 static t_pdbfile *read_pdbf(const char *fn)
77 get_stx_coordnum (fn, &natoms);
78 init_t_atoms(&(pdbf->atoms), natoms, FALSE);
79 snew(pdbf->x, natoms);
80 read_stx_conf(fn, buf, &pdbf->atoms, pdbf->x, NULL, &pdbf->ePBC, pdbf->box);
81 fp = gmx_ffopen(fn, "r");
84 ptr = fgets2(buf, 255, fp);
87 if (strstr(buf, "Intermolecular") != NULL)
89 ptr = strchr(buf, '=');
90 sscanf(ptr+1, "%lf", &e);
93 else if (strstr(buf, "Estimated Free") != NULL)
95 ptr = strchr(buf, '=');
96 sscanf(ptr+1, "%lf", &e);
107 static t_pdbfile **read_em_all(const char *fn, int *npdbf)
109 t_pdbfile **pdbf = 0;
111 char buf[256], name[256];
115 buf[strlen(buf)-4] = '\0';
121 sprintf(name, "%s_%d.pdb", buf, i+1);
122 if ((bExist = gmx_fexist(name)) == TRUE)
124 pdbf[i] = read_pdbf(name);
125 pdbf[i]->index = i+1;
130 srenew(pdbf, maxpdbf);
138 printf("Found %d pdb files\n", i);
143 static gmx_bool bFreeSort = FALSE;
145 static int pdbf_comp(const void *a, const void *b)
151 pa = *(t_pdbfile **)a;
152 pb = *(t_pdbfile **)b;
154 dc = pa->cluster_id - pb->cluster_id;
160 x = pa->efree - pb->efree;
164 x = pa->edocked - pb->edocked;
186 static void analyse_em_all(int npdb, t_pdbfile *pdbf[], const char *edocked,
187 const char *efree, const output_env_t oenv)
192 for (bFreeSort = FALSE; (bFreeSort <= TRUE); bFreeSort++)
194 qsort(pdbf, npdb, sizeof(pdbf[0]), pdbf_comp);
195 fp = xvgropen(bFreeSort ? efree : edocked,
196 etitles[bFreeSort], "()", "E (kJ/mol)", oenv);
197 for (i = 0; (i < npdb); i++)
199 fprintf(fp, "%12lf\n", bFreeSort ? pdbf[i]->efree : pdbf[i]->edocked);
205 static void clust_stat(FILE *fp, int start, int end, t_pdbfile *pdbf[])
211 ed = gmx_stats_init();
212 ef = gmx_stats_init();
213 for (i = start; (i < end); i++)
215 gmx_stats_add_point(ed, i-start, pdbf[i]->edocked, 0, 0);
216 gmx_stats_add_point(ef, i-start, pdbf[i]->efree, 0, 0);
218 gmx_stats_get_ase(ed, &aver, &sigma, NULL);
219 fprintf(fp, " <%12s> = %8.3f (+/- %6.3f)\n", etitles[FALSE], aver, sigma);
220 gmx_stats_get_ase(ef, &aver, &sigma, NULL);
221 fprintf(fp, " <%12s> = %8.3f (+/- %6.3f)\n", etitles[TRUE], aver, sigma);
228 static real rmsd_dist(t_pdbfile *pa, t_pdbfile *pb, gmx_bool bRMSD)
237 for (i = 0; (i < pa->atoms.nr); i++)
239 rvec_sub(pa->x[i], pb->x[i], dx);
240 rmsd += iprod(dx, dx);
242 rmsd = sqrt(rmsd/pa->atoms.nr);
249 for (i = 0; (i < pa->atoms.nr); i++)
251 rvec_inc(acm, pa->x[i]);
252 rvec_inc(bcm, pb->x[i]);
254 rvec_sub(acm, bcm, dx);
255 for (i = 0; (i < DIM); i++)
257 dx[i] /= pa->atoms.nr;
264 static void line(FILE *fp)
266 fprintf(fp, " - - - - - - - - - -\n");
269 static void cluster_em_all(FILE *fp, int npdb, t_pdbfile *pdbf[],
270 gmx_bool bFree, gmx_bool bRMSD, real cutoff)
277 qsort(pdbf, npdb, sizeof(pdbf[0]), pdbf_comp);
279 fprintf(fp, "Statistics over all structures:\n");
280 clust_stat(fp, 0, npdb, pdbf);
283 /* Index to first structure in a cluster */
287 for (i = 1; (i < npdb); i++)
289 for (j = 0; (j < ncluster); j++)
291 rmsd = rmsd_dist(pdbf[cndx[j]], pdbf[i], bRMSD);
294 /* Structure i is in cluster j */
295 pdbf[i]->cluster_id = pdbf[cndx[j]]->cluster_id;
301 /* New cluster! Cool! */
303 pdbf[i]->cluster_id = ncluster;
307 cndx[ncluster] = npdb;
308 qsort(pdbf, npdb, sizeof(pdbf[0]), pdbf_comp);
312 for (i = 1; (i < npdb); i++)
314 if (pdbf[i]->cluster_id != pdbf[i-1]->cluster_id)
322 gmx_fatal(FARGS, "Consistency error: j = %d, ncluster = %d", j, ncluster);
325 fprintf(fp, "I found %d clusters based on %s and %s with a %.3f nm cut-off\n",
326 ncluster, etitles[bFree], bRMSD ? "RMSD" : "distance", cutoff);
328 for (j = 0; (j < ncluster); j++)
330 fprintf(fp, "Cluster: %3d %s: %10.5f kJ/mol %3d elements\n",
332 bFree ? pdbf[cndx[j]]->efree : pdbf[cndx[j]]->edocked,
334 clust_stat(fp, cndx[j], cndx[j+1], pdbf);
335 for (k = cndx[j]; (k < cndx[j+1]); k++)
337 fprintf(fp, " %3d", pdbf[k]->index);
345 int gmx_anadock(int argc, char *argv[])
347 const char *desc[] = {
348 "[THISMODULE] analyses the results of an Autodock run and clusters the",
349 "structures together, based on distance or RMSD. The docked energy",
350 "and free energy estimates are analysed, and for each cluster the",
351 "energy statistics are printed.[PAR]",
352 "An alternative approach to this is to cluster the structures first",
353 "using [gmx-cluster] and then sort the clusters on either lowest",
354 "energy or average energy."
357 { efPDB, "-f", NULL, ffREAD },
358 { efXVG, "-od", "edocked", ffWRITE },
359 { efXVG, "-of", "efree", ffWRITE },
360 { efLOG, "-g", "anadock", ffWRITE }
363 #define NFILE asize(fnm)
364 static gmx_bool bFree = FALSE, bRMS = TRUE;
365 static real cutoff = 0.2;
367 { "-free", FALSE, etBOOL, {&bFree},
368 "Use Free energy estimate from autodock for sorting the classes" },
369 { "-rms", FALSE, etBOOL, {&bRMS},
370 "Cluster on RMS or distance" },
371 { "-cutoff", FALSE, etREAL, {&cutoff},
372 "Maximum RMSD/distance for belonging to the same cluster" }
374 #define NPA asize(pa)
377 t_pdbfile **pdbf = NULL;
380 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, NPA, pa, asize(desc), desc, 0,
386 fp = gmx_ffopen(opt2fn("-g", NFILE, fnm), "w");
387 please_cite(stdout, "Hetenyi2002b");
388 please_cite(fp, "Hetenyi2002b");
390 pdbf = read_em_all(opt2fn("-f", NFILE, fnm), &npdbf);
392 analyse_em_all(npdbf, pdbf, opt2fn("-od", NFILE, fnm), opt2fn("-of", NFILE, fnm),
395 cluster_em_all(fp, npdbf, pdbf, bFree, bRMS, cutoff);