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-2007, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015, 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.
41 #include "gromacs/commandline/pargs.h"
42 #include "gromacs/fileio/matio.h"
43 #include "gromacs/fileio/tpxio.h"
44 #include "gromacs/fileio/trxio.h"
45 #include "gromacs/fileio/xvgr.h"
46 #include "gromacs/gmxana/gmx_ana.h"
47 #include "gromacs/gmxana/gstat.h"
48 #include "gromacs/legacyheaders/macros.h"
49 #include "gromacs/legacyheaders/nrnb.h"
50 #include "gromacs/legacyheaders/typedefs.h"
51 #include "gromacs/math/units.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/pbcutil/pbc.h"
54 #include "gromacs/topology/index.h"
55 #include "gromacs/topology/mtop_util.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/futil.h"
58 #include "gromacs/utility/smalloc.h"
60 static void clust_size(const char *ndx, const char *trx, const char *xpm,
61 const char *xpmw, const char *ncl, const char *acl,
62 const char *mcl, const char *histo, const char *tempf,
63 const char *mcn, gmx_bool bMol, gmx_bool bPBC, const char *tpr,
64 real cut, int nskip, int nlevels,
65 t_rgb rmid, t_rgb rhi, int ndf,
66 const output_env_t oenv)
68 FILE *fp, *gp, *hp, *tp;
69 atom_id *index = NULL;
72 rvec *x = NULL, *v = NULL, dx;
76 gmx_bool bSame, bTPRwarn = TRUE;
80 gmx_mtop_t *mtop = NULL;
83 gmx_mtop_atomlookup_t alook;
85 int version, generation, ii, jj, nsame;
87 /* Cluster size distribution (matrix) */
88 real **cs_dist = NULL;
89 real tf, dx2, cut2, *t_x = NULL, *t_y, cmid, cmax, cav, ekin;
90 int i, j, k, ai, aj, ak, ci, cj, nframe, nclust, n_x, n_y, max_size = 0;
91 int *clust_index, *clust_size, max_clust_size, max_clust_ind, nav, nhisto;
92 t_rgb rlo = { 1.0, 1.0, 1.0 };
94 clear_trxframe(&fr, TRUE);
95 sprintf(timebuf, "Time (%s)", output_env_get_time_unit(oenv));
96 tf = output_env_get_time_factor(oenv);
97 fp = xvgropen(ncl, "Number of clusters", timebuf, "N", oenv);
98 gp = xvgropen(acl, "Average cluster size", timebuf, "#molecules", oenv);
99 hp = xvgropen(mcl, "Max cluster size", timebuf, "#molecules", oenv);
100 tp = xvgropen(tempf, "Temperature of largest cluster", timebuf, "T (K)",
103 if (!read_first_frame(oenv, &status, trx, &fr, TRX_NEED_X | TRX_READ_V))
114 read_tpxheader(tpr, &tpxh, TRUE, &version, &generation);
115 if (tpxh.natoms != natoms)
117 gmx_fatal(FARGS, "tpr (%d atoms) and trajectory (%d atoms) do not match!",
118 tpxh.natoms, natoms);
120 ePBC = read_tpx(tpr, NULL, NULL, &natoms, NULL, NULL, NULL, mtop);
128 tfac = ndf/(3.0*natoms);
135 printf("Using molecules rather than atoms. Not reading index file %s\n",
138 mols = &(mtop->mols);
140 /* Make dummy index */
143 for (i = 0; (i < nindex); i++)
147 gname = gmx_strdup("mols");
151 rd_index(ndx, 1, &nindex, &index, &gname);
154 alook = gmx_mtop_atomlookup_init(mtop);
156 snew(clust_index, nindex);
157 snew(clust_size, nindex);
162 for (i = 0; (i < nindex); i++)
170 if ((nskip == 0) || ((nskip > 0) && ((nframe % nskip) == 0)))
174 set_pbc(&pbc, ePBC, fr.box);
179 /* Put all atoms/molecules in their own cluster, with size 1 */
180 for (i = 0; (i < nindex); i++)
182 /* Cluster index is indexed with atom index number */
184 /* Cluster size is indexed with cluster number */
188 /* Loop over atoms */
189 for (i = 0; (i < nindex); i++)
194 /* Loop over atoms (only half a matrix) */
195 for (j = i+1; (j < nindex); j++)
199 /* If they are not in the same cluster already */
204 /* Compute distance */
208 for (ii = mols->index[ai]; !bSame && (ii < mols->index[ai+1]); ii++)
210 for (jj = mols->index[aj]; !bSame && (jj < mols->index[aj+1]); jj++)
214 pbc_dx(&pbc, x[ii], x[jj], dx);
218 rvec_sub(x[ii], x[jj], dx);
221 bSame = (dx2 < cut2);
229 pbc_dx(&pbc, x[ai], x[aj], dx);
233 rvec_sub(x[ai], x[aj], dx);
236 bSame = (dx2 < cut2);
238 /* If distance less than cut-off */
241 /* Merge clusters: check for all atoms whether they are in
242 * cluster cj and if so, put them in ci
244 for (k = 0; (k < nindex); k++)
246 if (clust_index[k] == cj)
248 if (clust_size[cj] <= 0)
250 gmx_fatal(FARGS, "negative cluster size %d for element %d",
264 t_x[n_x-1] = fr.time*tf;
265 srenew(cs_dist, n_x);
266 snew(cs_dist[n_x-1], nindex);
270 for (i = 0; (i < nindex); i++)
273 if (ci > max_clust_size)
281 cs_dist[n_x-1][ci-1] += 1.0;
282 max_size = max(max_size, ci);
290 fprintf(fp, "%14.6e %10d\n", fr.time, nclust);
293 fprintf(gp, "%14.6e %10.3f\n", fr.time, cav/nav);
295 fprintf(hp, "%14.6e %10d\n", fr.time, max_clust_size);
297 /* Analyse velocities, if present */
304 printf("You need a [REF].tpr[ref] file to analyse temperatures\n");
311 /* Loop over clusters and for each cluster compute 1/2 m v^2 */
312 if (max_clust_ind >= 0)
315 for (i = 0; (i < nindex); i++)
317 if (clust_index[i] == max_clust_ind)
320 gmx_mtop_atomnr_to_atom(alook, ai, &atom);
321 ekin += 0.5*atom->m*iprod(v[ai], v[ai]);
324 temp = (ekin*2.0)/(3.0*tfac*max_clust_size*BOLTZ);
325 fprintf(tp, "%10.3f %10.3f\n", fr.time, temp);
331 while (read_next_frame(oenv, status, &fr));
338 gmx_mtop_atomlookup_destroy(alook);
340 if (max_clust_ind >= 0)
342 fp = gmx_ffopen(mcn, "w");
343 fprintf(fp, "[ max_clust ]\n");
344 for (i = 0; (i < nindex); i++)
346 if (clust_index[i] == max_clust_ind)
350 for (j = mols->index[i]; (j < mols->index[i+1]); j++)
352 fprintf(fp, "%d\n", j+1);
357 fprintf(fp, "%d\n", index[i]+1);
364 /* Print the real distribution cluster-size/numer, averaged over the trajectory. */
365 fp = xvgropen(histo, "Cluster size distribution", "Cluster size", "()", oenv);
367 fprintf(fp, "%5d %8.3f\n", 0, 0.0);
368 for (j = 0; (j < max_size); j++)
371 for (i = 0; (i < n_x); i++)
373 nelem += cs_dist[i][j];
375 fprintf(fp, "%5d %8.3f\n", j+1, nelem/n_x);
376 nhisto += (int)((j+1)*nelem/n_x);
378 fprintf(fp, "%5d %8.3f\n", j+1, 0.0);
381 fprintf(stderr, "Total number of atoms in clusters = %d\n", nhisto);
383 /* Look for the smallest entry that is not zero
384 * This will make that zero is white, and not zero is coloured.
388 for (i = 0; (i < n_x); i++)
390 for (j = 0; (j < max_size); j++)
392 if ((cs_dist[i][j] > 0) && (cs_dist[i][j] < cmid))
394 cmid = cs_dist[i][j];
396 cmax = max(cs_dist[i][j], cmax);
399 fprintf(stderr, "cmid: %g, cmax: %g, max_size: %d\n", cmid, cmax, max_size);
401 fp = gmx_ffopen(xpm, "w");
402 write_xpm3(fp, 0, "Cluster size distribution", "# clusters", timebuf, "Size",
403 n_x, max_size, t_x, t_y, cs_dist, 0, cmid, cmax,
404 rlo, rmid, rhi, &nlevels);
408 for (i = 0; (i < n_x); i++)
410 for (j = 0; (j < max_size); j++)
412 cs_dist[i][j] *= (j+1);
413 if ((cs_dist[i][j] > 0) && (cs_dist[i][j] < cmid))
415 cmid = cs_dist[i][j];
417 cmax = max(cs_dist[i][j], cmax);
420 fprintf(stderr, "cmid: %g, cmax: %g, max_size: %d\n", cmid, cmax, max_size);
421 fp = gmx_ffopen(xpmw, "w");
422 write_xpm3(fp, 0, "Weighted cluster size distribution", "Fraction", timebuf,
423 "Size", n_x, max_size, t_x, t_y, cs_dist, 0, cmid, cmax,
424 rlo, rmid, rhi, &nlevels);
432 int gmx_clustsize(int argc, char *argv[])
434 const char *desc[] = {
435 "[THISMODULE] computes the size distributions of molecular/atomic clusters in",
436 "the gas phase. The output is given in the form of an [REF].xpm[ref] file.",
437 "The total number of clusters is written to an [REF].xvg[ref] file.[PAR]",
438 "When the [TT]-mol[tt] option is given clusters will be made out of",
439 "molecules rather than atoms, which allows clustering of large molecules.",
440 "In this case an index file would still contain atom numbers",
441 "or your calculation will die with a SEGV.[PAR]",
442 "When velocities are present in your trajectory, the temperature of",
443 "the largest cluster will be printed in a separate [REF].xvg[ref] file assuming",
444 "that the particles are free to move. If you are using constraints,",
445 "please correct the temperature. For instance water simulated with SHAKE",
446 "or SETTLE will yield a temperature that is 1.5 times too low. You can",
447 "compensate for this with the [TT]-ndf[tt] option. Remember to take the removal",
448 "of center of mass motion into account.[PAR]",
449 "The [TT]-mc[tt] option will produce an index file containing the",
450 "atom numbers of the largest cluster."
453 static real cutoff = 0.35;
454 static int nskip = 0;
455 static int nlevels = 20;
457 static gmx_bool bMol = FALSE;
458 static gmx_bool bPBC = TRUE;
459 static rvec rlo = { 1.0, 1.0, 0.0 };
460 static rvec rhi = { 0.0, 0.0, 1.0 };
465 { "-cut", FALSE, etREAL, {&cutoff},
466 "Largest distance (nm) to be considered in a cluster" },
467 { "-mol", FALSE, etBOOL, {&bMol},
468 "Cluster molecules rather than atoms (needs [REF].tpr[ref] file)" },
469 { "-pbc", FALSE, etBOOL, {&bPBC},
470 "Use periodic boundary conditions" },
471 { "-nskip", FALSE, etINT, {&nskip},
472 "Number of frames to skip between writing" },
473 { "-nlevels", FALSE, etINT, {&nlevels},
474 "Number of levels of grey in [REF].xpm[ref] output" },
475 { "-ndf", FALSE, etINT, {&ndf},
476 "Number of degrees of freedom of the entire system for temperature calculation. If not set, the number of atoms times three is used." },
477 { "-rgblo", FALSE, etRVEC, {rlo},
478 "RGB values for the color of the lowest occupied cluster size" },
479 { "-rgbhi", FALSE, etRVEC, {rhi},
480 "RGB values for the color of the highest occupied cluster size" }
482 #define NPA asize(pa)
483 const char *fnNDX, *fnTPR;
488 { efTRX, "-f", NULL, ffREAD },
489 { efTPR, NULL, NULL, ffOPTRD },
490 { efNDX, NULL, NULL, ffOPTRD },
491 { efXPM, "-o", "csize", ffWRITE },
492 { efXPM, "-ow", "csizew", ffWRITE },
493 { efXVG, "-nc", "nclust", ffWRITE },
494 { efXVG, "-mc", "maxclust", ffWRITE },
495 { efXVG, "-ac", "avclust", ffWRITE },
496 { efXVG, "-hc", "histo-clust", ffWRITE },
497 { efXVG, "-temp", "temp", ffOPTWR },
498 { efNDX, "-mcn", "maxclust", ffOPTWR }
500 #define NFILE asize(fnm)
502 if (!parse_common_args(&argc, argv,
503 PCA_CAN_VIEW | PCA_CAN_TIME | PCA_TIME_UNIT,
504 NFILE, fnm, NPA, pa, asize(desc), desc, 0, NULL, &oenv))
509 fnNDX = ftp2fn_null(efNDX, NFILE, fnm);
510 rgblo.r = rlo[XX], rgblo.g = rlo[YY], rgblo.b = rlo[ZZ];
511 rgbhi.r = rhi[XX], rgbhi.g = rhi[YY], rgbhi.b = rhi[ZZ];
513 fnTPR = ftp2fn_null(efTPR, NFILE, fnm);
516 gmx_fatal(FARGS, "You need a tpr file for the -mol option");
519 clust_size(fnNDX, ftp2fn(efTRX, NFILE, fnm), opt2fn("-o", NFILE, fnm),
520 opt2fn("-ow", NFILE, fnm),
521 opt2fn("-nc", NFILE, fnm), opt2fn("-ac", NFILE, fnm),
522 opt2fn("-mc", NFILE, fnm), opt2fn("-hc", NFILE, fnm),
523 opt2fn("-temp", NFILE, fnm), opt2fn("-mcn", NFILE, fnm),
525 cutoff, nskip, nlevels, rgblo, rgbhi, ndf, oenv);