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.
42 #include "gromacs/fileio/confio.h"
44 #include "gromacs/utility/futil.h"
45 #include "gromacs/utility/fatalerror.h"
46 #include "gromacs/utility/smalloc.h"
47 #include "gromacs/utility/cstringutil.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/statistics/statistics.h"
50 #include "gromacs/commandline/pargs.h"
52 #include "gromacs/fileio/xvgr.h"
55 static const char *etitles[] = { "E-docked", "Free Energy" };
59 int index, cluster_id;
66 static t_pdbfile *read_pdbf(const char *fn)
75 get_stx_coordnum (fn, &natoms);
76 init_t_atoms(&(pdbf->atoms), natoms, FALSE);
77 snew(pdbf->x, natoms);
78 read_stx_conf(fn, buf, &pdbf->atoms, pdbf->x, NULL, &pdbf->ePBC, pdbf->box);
79 fp = gmx_ffopen(fn, "r");
82 ptr = fgets2(buf, 255, fp);
85 if (strstr(buf, "Intermolecular") != NULL)
87 ptr = strchr(buf, '=');
88 sscanf(ptr+1, "%lf", &e);
91 else if (strstr(buf, "Estimated Free") != NULL)
93 ptr = strchr(buf, '=');
94 sscanf(ptr+1, "%lf", &e);
105 static t_pdbfile **read_em_all(const char *fn, int *npdbf)
107 t_pdbfile **pdbf = 0;
109 char buf[256], name[256];
113 buf[strlen(buf)-4] = '\0';
119 sprintf(name, "%s_%d.pdb", buf, i+1);
120 if ((bExist = gmx_fexist(name)) == TRUE)
122 pdbf[i] = read_pdbf(name);
123 pdbf[i]->index = i+1;
128 srenew(pdbf, maxpdbf);
136 printf("Found %d pdb files\n", i);
141 static gmx_bool bFreeSort = FALSE;
143 static int pdbf_comp(const void *a, const void *b)
149 pa = *(t_pdbfile **)a;
150 pb = *(t_pdbfile **)b;
152 dc = pa->cluster_id - pb->cluster_id;
158 x = pa->efree - pb->efree;
162 x = pa->edocked - pb->edocked;
184 static void analyse_em_all(int npdb, t_pdbfile *pdbf[], const char *edocked,
185 const char *efree, const output_env_t oenv)
190 for (bFreeSort = FALSE; (bFreeSort <= TRUE); bFreeSort++)
192 qsort(pdbf, npdb, sizeof(pdbf[0]), pdbf_comp);
193 fp = xvgropen(bFreeSort ? efree : edocked,
194 etitles[bFreeSort], "()", "E (kJ/mol)", oenv);
195 for (i = 0; (i < npdb); i++)
197 fprintf(fp, "%12lf\n", bFreeSort ? pdbf[i]->efree : pdbf[i]->edocked);
203 static void clust_stat(FILE *fp, int start, int end, t_pdbfile *pdbf[])
209 ed = gmx_stats_init();
210 ef = gmx_stats_init();
211 for (i = start; (i < end); i++)
213 gmx_stats_add_point(ed, i-start, pdbf[i]->edocked, 0, 0);
214 gmx_stats_add_point(ef, i-start, pdbf[i]->efree, 0, 0);
216 gmx_stats_get_ase(ed, &aver, &sigma, NULL);
217 fprintf(fp, " <%12s> = %8.3f (+/- %6.3f)\n", etitles[FALSE], aver, sigma);
218 gmx_stats_get_ase(ef, &aver, &sigma, NULL);
219 fprintf(fp, " <%12s> = %8.3f (+/- %6.3f)\n", etitles[TRUE], aver, sigma);
226 static real rmsd_dist(t_pdbfile *pa, t_pdbfile *pb, gmx_bool bRMSD)
235 for (i = 0; (i < pa->atoms.nr); i++)
237 rvec_sub(pa->x[i], pb->x[i], dx);
238 rmsd += iprod(dx, dx);
240 rmsd = sqrt(rmsd/pa->atoms.nr);
247 for (i = 0; (i < pa->atoms.nr); i++)
249 rvec_inc(acm, pa->x[i]);
250 rvec_inc(bcm, pb->x[i]);
252 rvec_sub(acm, bcm, dx);
253 for (i = 0; (i < DIM); i++)
255 dx[i] /= pa->atoms.nr;
262 static void line(FILE *fp)
264 fprintf(fp, " - - - - - - - - - -\n");
267 static void cluster_em_all(FILE *fp, int npdb, t_pdbfile *pdbf[],
268 gmx_bool bFree, gmx_bool bRMSD, real cutoff)
275 qsort(pdbf, npdb, sizeof(pdbf[0]), pdbf_comp);
277 fprintf(fp, "Statistics over all structures:\n");
278 clust_stat(fp, 0, npdb, pdbf);
281 /* Index to first structure in a cluster */
285 for (i = 1; (i < npdb); i++)
287 for (j = 0; (j < ncluster); j++)
289 rmsd = rmsd_dist(pdbf[cndx[j]], pdbf[i], bRMSD);
292 /* Structure i is in cluster j */
293 pdbf[i]->cluster_id = pdbf[cndx[j]]->cluster_id;
299 /* New cluster! Cool! */
301 pdbf[i]->cluster_id = ncluster;
305 cndx[ncluster] = npdb;
306 qsort(pdbf, npdb, sizeof(pdbf[0]), pdbf_comp);
310 for (i = 1; (i < npdb); i++)
312 if (pdbf[i]->cluster_id != pdbf[i-1]->cluster_id)
320 gmx_fatal(FARGS, "Consistency error: j = %d, ncluster = %d", j, ncluster);
323 fprintf(fp, "I found %d clusters based on %s and %s with a %.3f nm cut-off\n",
324 ncluster, etitles[bFree], bRMSD ? "RMSD" : "distance", cutoff);
326 for (j = 0; (j < ncluster); j++)
328 fprintf(fp, "Cluster: %3d %s: %10.5f kJ/mol %3d elements\n",
330 bFree ? pdbf[cndx[j]]->efree : pdbf[cndx[j]]->edocked,
332 clust_stat(fp, cndx[j], cndx[j+1], pdbf);
333 for (k = cndx[j]; (k < cndx[j+1]); k++)
335 fprintf(fp, " %3d", pdbf[k]->index);
343 int gmx_anadock(int argc, char *argv[])
345 const char *desc[] = {
346 "[THISMODULE] analyses the results of an Autodock run and clusters the",
347 "structures together, based on distance or RMSD. The docked energy",
348 "and free energy estimates are analysed, and for each cluster the",
349 "energy statistics are printed.[PAR]",
350 "An alternative approach to this is to cluster the structures first",
351 "using [gmx-cluster] and then sort the clusters on either lowest",
352 "energy or average energy."
355 { efPDB, "-f", NULL, ffREAD },
356 { efXVG, "-od", "edocked", ffWRITE },
357 { efXVG, "-of", "efree", ffWRITE },
358 { efLOG, "-g", "anadock", ffWRITE }
361 #define NFILE asize(fnm)
362 static gmx_bool bFree = FALSE, bRMS = TRUE;
363 static real cutoff = 0.2;
365 { "-free", FALSE, etBOOL, {&bFree},
366 "Use Free energy estimate from autodock for sorting the classes" },
367 { "-rms", FALSE, etBOOL, {&bRMS},
368 "Cluster on RMS or distance" },
369 { "-cutoff", FALSE, etREAL, {&cutoff},
370 "Maximum RMSD/distance for belonging to the same cluster" }
372 #define NPA asize(pa)
375 t_pdbfile **pdbf = NULL;
378 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, NPA, pa, asize(desc), desc, 0,
384 fp = gmx_ffopen(opt2fn("-g", NFILE, fnm), "w");
385 please_cite(stdout, "Hetenyi2002b");
386 please_cite(fp, "Hetenyi2002b");
388 pdbf = read_em_all(opt2fn("-f", NFILE, fnm), &npdbf);
390 analyse_em_all(npdbf, pdbf, opt2fn("-od", NFILE, fnm), opt2fn("-of", NFILE, fnm),
393 cluster_em_all(fp, npdbf, pdbf, bFree, bRMS, cutoff);