3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
42 #include "gmx_fatal.h"
46 #include "gmx_statistics.h"
52 static const char *etitles[] = { "E-docked", "Free Energy" };
56 int index, cluster_id;
63 static t_pdbfile *read_pdbf(const char *fn)
72 get_stx_coordnum (fn, &natoms);
73 init_t_atoms(&(pdbf->atoms), natoms, FALSE);
74 snew(pdbf->x, natoms);
75 read_stx_conf(fn, buf, &pdbf->atoms, pdbf->x, NULL, &pdbf->ePBC, pdbf->box);
79 ptr = fgets2(buf, 255, fp);
82 if (strstr(buf, "Intermolecular") != NULL)
84 ptr = strchr(buf, '=');
85 sscanf(ptr+1, "%lf", &e);
88 else if (strstr(buf, "Estimated Free") != NULL)
90 ptr = strchr(buf, '=');
91 sscanf(ptr+1, "%lf", &e);
102 static t_pdbfile **read_em_all(const char *fn, int *npdbf)
104 t_pdbfile **pdbf = 0;
106 char buf[256], name[256];
110 buf[strlen(buf)-4] = '\0';
116 sprintf(name, "%s_%d.pdb", buf, i+1);
117 if ((bExist = gmx_fexist(name)) == TRUE)
119 pdbf[i] = read_pdbf(name);
120 pdbf[i]->index = i+1;
125 srenew(pdbf, maxpdbf);
133 printf("Found %d pdb files\n", i);
138 static gmx_bool bFreeSort = FALSE;
140 static int pdbf_comp(const void *a, const void *b)
146 pa = *(t_pdbfile **)a;
147 pb = *(t_pdbfile **)b;
149 dc = pa->cluster_id - pb->cluster_id;
155 x = pa->efree - pb->efree;
159 x = pa->edocked - pb->edocked;
181 static void analyse_em_all(int npdb, t_pdbfile *pdbf[], const char *edocked,
182 const char *efree, const output_env_t oenv)
187 for (bFreeSort = FALSE; (bFreeSort <= TRUE); bFreeSort++)
189 qsort(pdbf, npdb, sizeof(pdbf[0]), pdbf_comp);
190 fp = xvgropen(bFreeSort ? efree : edocked,
191 etitles[bFreeSort], "()", "E (kJ/mol)", oenv);
192 for (i = 0; (i < npdb); i++)
194 fprintf(fp, "%12lf\n", bFreeSort ? pdbf[i]->efree : pdbf[i]->edocked);
200 static void clust_stat(FILE *fp, int start, int end, t_pdbfile *pdbf[])
206 ed = gmx_stats_init();
207 ef = gmx_stats_init();
208 for (i = start; (i < end); i++)
210 gmx_stats_add_point(ed, i-start, pdbf[i]->edocked, 0, 0);
211 gmx_stats_add_point(ef, i-start, pdbf[i]->efree, 0, 0);
213 gmx_stats_get_ase(ed, &aver, &sigma, NULL);
214 fprintf(fp, " <%12s> = %8.3f (+/- %6.3f)\n", etitles[FALSE], aver, sigma);
215 gmx_stats_get_ase(ef, &aver, &sigma, NULL);
216 fprintf(fp, " <%12s> = %8.3f (+/- %6.3f)\n", etitles[TRUE], aver, sigma);
223 static real rmsd_dist(t_pdbfile *pa, t_pdbfile *pb, gmx_bool bRMSD)
232 for (i = 0; (i < pa->atoms.nr); i++)
234 rvec_sub(pa->x[i], pb->x[i], dx);
235 rmsd += iprod(dx, dx);
237 rmsd = sqrt(rmsd/pa->atoms.nr);
244 for (i = 0; (i < pa->atoms.nr); i++)
246 rvec_inc(acm, pa->x[i]);
247 rvec_inc(bcm, pb->x[i]);
249 rvec_sub(acm, bcm, dx);
250 for (i = 0; (i < DIM); i++)
252 dx[i] /= pa->atoms.nr;
259 static void line(FILE *fp)
261 fprintf(fp, " - - - - - - - - - -\n");
264 static void cluster_em_all(FILE *fp, int npdb, t_pdbfile *pdbf[],
265 const char *pdbout, gmx_bool bFree, gmx_bool bRMSD, real cutoff)
272 qsort(pdbf, npdb, sizeof(pdbf[0]), pdbf_comp);
274 fprintf(fp, "Statistics over all structures:\n");
275 clust_stat(fp, 0, npdb, pdbf);
278 /* Index to first structure in a cluster */
282 for (i = 1; (i < npdb); i++)
284 for (j = 0; (j < ncluster); j++)
286 rmsd = rmsd_dist(pdbf[cndx[j]], pdbf[i], bRMSD);
289 /* Structure i is in cluster j */
290 pdbf[i]->cluster_id = pdbf[cndx[j]]->cluster_id;
296 /* New cluster! Cool! */
298 pdbf[i]->cluster_id = ncluster;
302 cndx[ncluster] = npdb;
303 qsort(pdbf, npdb, sizeof(pdbf[0]), pdbf_comp);
307 for (i = 1; (i < npdb); i++)
309 if (pdbf[i]->cluster_id != pdbf[i-1]->cluster_id)
317 gmx_fatal(FARGS, "Consistency error: j = %d, ncluster = %d", j, ncluster);
320 fprintf(fp, "I found %d clusters based on %s and %s with a %.3f nm cut-off\n",
321 ncluster, etitles[bFree], bRMSD ? "RMSD" : "distance", cutoff);
323 for (j = 0; (j < ncluster); j++)
325 fprintf(fp, "Cluster: %3d %s: %10.5f kJ/mol %3d elements\n",
327 bFree ? pdbf[cndx[j]]->efree : pdbf[cndx[j]]->edocked,
329 clust_stat(fp, cndx[j], cndx[j+1], pdbf);
330 for (k = cndx[j]; (k < cndx[j+1]); k++)
332 fprintf(fp, " %3d", pdbf[k]->index);
340 int gmx_anadock(int argc, char *argv[])
342 const char *desc[] = {
343 "[TT]g_anadock[tt] analyses the results of an Autodock run and clusters the",
344 "structures together, based on distance or RMSD. The docked energy",
345 "and free energy estimates are analysed, and for each cluster the",
346 "energy statistics are printed.[PAR]",
347 "An alternative approach to this is to cluster the structures first",
348 "using [TT]g_cluster[tt] and then sort the clusters on either lowest",
349 "energy or average energy."
352 { efPDB, "-f", NULL, ffREAD },
353 { efPDB, "-ox", "cluster", ffWRITE },
354 { efXVG, "-od", "edocked", ffWRITE },
355 { efXVG, "-of", "efree", ffWRITE },
356 { efLOG, "-g", "anadock", ffWRITE }
359 #define NFILE asize(fnm)
360 static gmx_bool bFree = FALSE, bRMS = TRUE;
361 static real cutoff = 0.2;
363 { "-free", FALSE, etBOOL, {&bFree},
364 "Use Free energy estimate from autodock for sorting the classes" },
365 { "-rms", FALSE, etBOOL, {&bRMS},
366 "Cluster on RMS or distance" },
367 { "-cutoff", FALSE, etREAL, {&cutoff},
368 "Maximum RMSD/distance for belonging to the same cluster" }
370 #define NPA asize(pa)
373 t_pdbfile **pdbf = NULL;
376 CopyRight(stderr, argv[0]);
377 parse_common_args(&argc, argv, 0, NFILE, fnm, NPA, pa, asize(desc), desc, 0,
380 fp = ffopen(opt2fn("-g", NFILE, fnm), "w");
381 please_cite(stdout, "Hetenyi2002b");
382 please_cite(fp, "Hetenyi2002b");
384 pdbf = read_em_all(opt2fn("-f", NFILE, fnm), &npdbf);
386 analyse_em_all(npdbf, pdbf, opt2fn("-od", NFILE, fnm), opt2fn("-of", NFILE, fnm),
389 cluster_em_all(fp, npdbf, pdbf, opt2fn("-ox", NFILE, fnm),
390 bFree, bRMS, cutoff);