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,2018, 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.
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" };
61 int index, cluster_id;
68 static t_pdbfile *read_pdbf(const char *fn)
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");
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);
88 else if (std::strstr(buf, "Estimated Free") != nullptr)
90 ptr = std::strchr(buf, '=');
91 sscanf(ptr+1, "%lf", &e);
100 static t_pdbfile **read_em_all(const char *fn, int *npdbf)
102 t_pdbfile **pdbf = nullptr;
104 char buf[256], name[256];
107 std::strcpy(buf, fn);
108 buf[std::strlen(buf)-4] = '\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;
123 srenew(pdbf, maxpdbf);
131 printf("Found %d pdb files\n", i);
136 static gmx_bool bFreeSort = FALSE;
138 static int pdbf_comp(const void *a, const void *b)
144 pa = *(t_pdbfile **)a;
145 pb = *(t_pdbfile **)b;
147 dc = pa->cluster_id - pb->cluster_id;
153 x = pa->efree - pb->efree;
157 x = pa->edocked - pb->edocked;
179 static void analyse_em_all(int npdb, t_pdbfile *pdbf[], const char *edocked,
180 const char *efree, const gmx_output_env_t *oenv)
185 for (int freeSort = 0; freeSort < 2; freeSort++)
187 bFreeSort = (freeSort == 1);
188 qsort(pdbf, npdb, sizeof(pdbf[0]), pdbf_comp);
189 fp = xvgropen(bFreeSort ? efree : edocked,
190 etitles[bFreeSort], "()", "E (kJ/mol)", oenv);
191 for (i = 0; (i < npdb); i++)
193 fprintf(fp, "%12lf\n", bFreeSort ? pdbf[i]->efree : pdbf[i]->edocked);
199 static void clust_stat(FILE *fp, int start, int end, t_pdbfile *pdbf[])
205 ed = gmx_stats_init();
206 ef = gmx_stats_init();
207 for (i = start; (i < end); i++)
209 gmx_stats_add_point(ed, i-start, pdbf[i]->edocked, 0, 0);
210 gmx_stats_add_point(ef, i-start, pdbf[i]->efree, 0, 0);
212 gmx_stats_get_ase(ed, &aver, &sigma, nullptr);
213 fprintf(fp, " <%12s> = %8.3f (+/- %6.3f)\n", etitles[FALSE], aver, sigma);
214 gmx_stats_get_ase(ef, &aver, &sigma, nullptr);
215 fprintf(fp, " <%12s> = %8.3f (+/- %6.3f)\n", etitles[TRUE], aver, sigma);
220 static real rmsd_dist(t_pdbfile *pa, t_pdbfile *pb, gmx_bool bRMSD)
229 for (i = 0; (i < pa->atoms.nr); i++)
231 rvec_sub(pa->x[i], pb->x[i], dx);
232 rmsd += iprod(dx, dx);
234 rmsd = std::sqrt(rmsd/pa->atoms.nr);
240 for (i = 0; (i < pa->atoms.nr); i++)
242 rvec_inc(acm, pa->x[i]);
243 rvec_inc(bcm, pb->x[i]);
245 rvec_sub(acm, bcm, dx);
246 for (i = 0; (i < DIM); i++)
248 dx[i] /= pa->atoms.nr;
255 static void line(FILE *fp)
257 fprintf(fp, " - - - - - - - - - -\n");
260 static void cluster_em_all(FILE *fp, int npdb, t_pdbfile *pdbf[],
261 gmx_bool bFree, gmx_bool bRMSD, real cutoff)
268 qsort(pdbf, npdb, sizeof(pdbf[0]), pdbf_comp);
270 fprintf(fp, "Statistics over all structures:\n");
271 clust_stat(fp, 0, npdb, pdbf);
274 /* Index to first structure in a cluster */
278 for (i = 1; (i < npdb); i++)
280 for (j = 0; (j < ncluster); j++)
282 rmsd = rmsd_dist(pdbf[cndx[j]], pdbf[i], bRMSD);
285 /* Structure i is in cluster j */
286 pdbf[i]->cluster_id = pdbf[cndx[j]]->cluster_id;
292 /* New cluster! Cool! */
294 pdbf[i]->cluster_id = ncluster;
298 cndx[ncluster] = npdb;
299 qsort(pdbf, npdb, sizeof(pdbf[0]), pdbf_comp);
303 for (i = 1; (i < npdb); i++)
305 if (pdbf[i]->cluster_id != pdbf[i-1]->cluster_id)
313 gmx_fatal(FARGS, "Consistency error: j = %d, ncluster = %d", j, ncluster);
316 fprintf(fp, "I found %d clusters based on %s and %s with a %.3f nm cut-off\n",
317 ncluster, etitles[bFree], bRMSD ? "RMSD" : "distance", cutoff);
319 for (j = 0; (j < ncluster); j++)
321 fprintf(fp, "Cluster: %3d %s: %10.5f kJ/mol %3d elements\n",
323 bFree ? pdbf[cndx[j]]->efree : pdbf[cndx[j]]->edocked,
325 clust_stat(fp, cndx[j], cndx[j+1], pdbf);
326 for (k = cndx[j]; (k < cndx[j+1]); k++)
328 fprintf(fp, " %3d", pdbf[k]->index);
336 int gmx_anadock(int argc, char *argv[])
338 const char *desc[] = {
339 "[THISMODULE] analyses the results of an Autodock run and clusters the",
340 "structures together, based on distance or RMSD. The docked energy",
341 "and free energy estimates are analysed, and for each cluster the",
342 "energy statistics are printed.[PAR]",
343 "An alternative approach to this is to cluster the structures first",
344 "using [gmx-cluster] and then sort the clusters on either lowest",
345 "energy or average energy."
348 { efPDB, "-f", nullptr, ffREAD },
349 { efXVG, "-od", "edocked", ffWRITE },
350 { efXVG, "-of", "efree", ffWRITE },
351 { efLOG, "-g", "anadock", ffWRITE }
353 gmx_output_env_t *oenv;
354 #define NFILE asize(fnm)
355 static gmx_bool bFree = FALSE, bRMS = TRUE;
356 static real cutoff = 0.2;
358 { "-free", FALSE, etBOOL, {&bFree},
359 "Use Free energy estimate from autodock for sorting the classes" },
360 { "-rms", FALSE, etBOOL, {&bRMS},
361 "Cluster on RMS or distance" },
362 { "-cutoff", FALSE, etREAL, {&cutoff},
363 "Maximum RMSD/distance for belonging to the same cluster" }
365 #define NPA asize(pa)
368 t_pdbfile **pdbf = nullptr;
371 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, NPA, pa, asize(desc), desc, 0,
377 fp = gmx_ffopen(opt2fn("-g", NFILE, fnm), "w");
378 please_cite(stdout, "Hetenyi2002b");
379 please_cite(fp, "Hetenyi2002b");
381 pdbf = read_em_all(opt2fn("-f", NFILE, fnm), &npdbf);
383 analyse_em_all(npdbf, pdbf, opt2fn("-od", NFILE, fnm), opt2fn("-of", NFILE, fnm),
386 cluster_em_all(fp, npdbf, pdbf, bFree, bRMS, cutoff);