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, 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/legacyheaders/typedefs.h"
44 #include "gromacs/legacyheaders/macros.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/pbcutil/pbc.h"
47 #include "gromacs/commandline/pargs.h"
48 #include "gromacs/fileio/xvgr.h"
49 #include "gromacs/utility/futil.h"
50 #include "gromacs/fileio/tpxio.h"
51 #include "gromacs/fileio/trxio.h"
52 #include "gromacs/topology/index.h"
53 #include "gromacs/utility/smalloc.h"
54 #include "gromacs/legacyheaders/nrnb.h"
55 #include "gromacs/math/units.h"
56 #include "gromacs/legacyheaders/coulomb.h"
57 #include "gromacs/legacyheaders/pme.h"
59 #include "gromacs/fileio/matio.h"
60 #include "gromacs/topology/mtop_util.h"
63 #include "gromacs/utility/fatalerror.h"
65 static void clust_size(const char *ndx, const char *trx, const char *xpm,
66 const char *xpmw, const char *ncl, const char *acl,
67 const char *mcl, const char *histo, const char *tempf,
68 const char *mcn, gmx_bool bMol, gmx_bool bPBC, const char *tpr,
69 real cut, int nskip, int nlevels,
70 t_rgb rmid, t_rgb rhi, int ndf,
71 const output_env_t oenv)
73 FILE *fp, *gp, *hp, *tp;
74 atom_id *index = NULL;
77 rvec *x = NULL, *v = NULL, dx;
81 gmx_bool bSame, bTPRwarn = TRUE;
85 gmx_mtop_t *mtop = NULL;
88 gmx_mtop_atomlookup_t alook;
90 int version, generation, ii, jj, nsame;
92 /* Cluster size distribution (matrix) */
93 real **cs_dist = NULL;
94 real tf, dx2, cut2, *t_x = NULL, *t_y, cmid, cmax, cav, ekin;
95 int i, j, k, ai, aj, ak, ci, cj, nframe, nclust, n_x, n_y, max_size = 0;
96 int *clust_index, *clust_size, max_clust_size, max_clust_ind, nav, nhisto;
97 t_rgb rlo = { 1.0, 1.0, 1.0 };
99 clear_trxframe(&fr, TRUE);
100 sprintf(timebuf, "Time (%s)", output_env_get_time_unit(oenv));
101 tf = output_env_get_time_factor(oenv);
102 fp = xvgropen(ncl, "Number of clusters", timebuf, "N", oenv);
103 gp = xvgropen(acl, "Average cluster size", timebuf, "#molecules", oenv);
104 hp = xvgropen(mcl, "Max cluster size", timebuf, "#molecules", oenv);
105 tp = xvgropen(tempf, "Temperature of largest cluster", timebuf, "T (K)",
108 if (!read_first_frame(oenv, &status, trx, &fr, TRX_NEED_X | TRX_READ_V))
119 read_tpxheader(tpr, &tpxh, TRUE, &version, &generation);
120 if (tpxh.natoms != natoms)
122 gmx_fatal(FARGS, "tpr (%d atoms) and trajectory (%d atoms) do not match!",
123 tpxh.natoms, natoms);
125 ePBC = read_tpx(tpr, NULL, NULL, &natoms, NULL, NULL, NULL, mtop);
133 tfac = ndf/(3.0*natoms);
140 printf("Using molecules rather than atoms. Not reading index file %s\n",
143 mols = &(mtop->mols);
145 /* Make dummy index */
148 for (i = 0; (i < nindex); i++)
152 gname = gmx_strdup("mols");
156 rd_index(ndx, 1, &nindex, &index, &gname);
159 alook = gmx_mtop_atomlookup_init(mtop);
161 snew(clust_index, nindex);
162 snew(clust_size, nindex);
167 for (i = 0; (i < nindex); i++)
175 if ((nskip == 0) || ((nskip > 0) && ((nframe % nskip) == 0)))
179 set_pbc(&pbc, ePBC, fr.box);
184 /* Put all atoms/molecules in their own cluster, with size 1 */
185 for (i = 0; (i < nindex); i++)
187 /* Cluster index is indexed with atom index number */
189 /* Cluster size is indexed with cluster number */
193 /* Loop over atoms */
194 for (i = 0; (i < nindex); i++)
199 /* Loop over atoms (only half a matrix) */
200 for (j = i+1; (j < nindex); j++)
204 /* If they are not in the same cluster already */
209 /* Compute distance */
213 for (ii = mols->index[ai]; !bSame && (ii < mols->index[ai+1]); ii++)
215 for (jj = mols->index[aj]; !bSame && (jj < mols->index[aj+1]); jj++)
219 pbc_dx(&pbc, x[ii], x[jj], dx);
223 rvec_sub(x[ii], x[jj], dx);
226 bSame = (dx2 < cut2);
234 pbc_dx(&pbc, x[ai], x[aj], dx);
238 rvec_sub(x[ai], x[aj], dx);
241 bSame = (dx2 < cut2);
243 /* If distance less than cut-off */
246 /* Merge clusters: check for all atoms whether they are in
247 * cluster cj and if so, put them in ci
249 for (k = 0; (k < nindex); k++)
251 if (clust_index[k] == cj)
253 if (clust_size[cj] <= 0)
255 gmx_fatal(FARGS, "negative cluster size %d for element %d",
269 t_x[n_x-1] = fr.time*tf;
270 srenew(cs_dist, n_x);
271 snew(cs_dist[n_x-1], nindex);
275 for (i = 0; (i < nindex); i++)
278 if (ci > max_clust_size)
286 cs_dist[n_x-1][ci-1] += 1.0;
287 max_size = max(max_size, ci);
295 fprintf(fp, "%14.6e %10d\n", fr.time, nclust);
298 fprintf(gp, "%14.6e %10.3f\n", fr.time, cav/nav);
300 fprintf(hp, "%14.6e %10d\n", fr.time, max_clust_size);
302 /* Analyse velocities, if present */
309 printf("You need a [TT].tpr[tt] file to analyse temperatures\n");
316 /* Loop over clusters and for each cluster compute 1/2 m v^2 */
317 if (max_clust_ind >= 0)
320 for (i = 0; (i < nindex); i++)
322 if (clust_index[i] == max_clust_ind)
325 gmx_mtop_atomnr_to_atom(alook, ai, &atom);
326 ekin += 0.5*atom->m*iprod(v[ai], v[ai]);
329 temp = (ekin*2.0)/(3.0*tfac*max_clust_size*BOLTZ);
330 fprintf(tp, "%10.3f %10.3f\n", fr.time, temp);
336 while (read_next_frame(oenv, status, &fr));
343 gmx_mtop_atomlookup_destroy(alook);
345 if (max_clust_ind >= 0)
347 fp = gmx_ffopen(mcn, "w");
348 fprintf(fp, "[ max_clust ]\n");
349 for (i = 0; (i < nindex); i++)
351 if (clust_index[i] == max_clust_ind)
355 for (j = mols->index[i]; (j < mols->index[i+1]); j++)
357 fprintf(fp, "%d\n", j+1);
362 fprintf(fp, "%d\n", index[i]+1);
369 /* Print the real distribution cluster-size/numer, averaged over the trajectory. */
370 fp = xvgropen(histo, "Cluster size distribution", "Cluster size", "()", oenv);
372 fprintf(fp, "%5d %8.3f\n", 0, 0.0);
373 for (j = 0; (j < max_size); j++)
376 for (i = 0; (i < n_x); i++)
378 nelem += cs_dist[i][j];
380 fprintf(fp, "%5d %8.3f\n", j+1, nelem/n_x);
381 nhisto += (int)((j+1)*nelem/n_x);
383 fprintf(fp, "%5d %8.3f\n", j+1, 0.0);
386 fprintf(stderr, "Total number of atoms in clusters = %d\n", nhisto);
388 /* Look for the smallest entry that is not zero
389 * This will make that zero is white, and not zero is coloured.
393 for (i = 0; (i < n_x); i++)
395 for (j = 0; (j < max_size); j++)
397 if ((cs_dist[i][j] > 0) && (cs_dist[i][j] < cmid))
399 cmid = cs_dist[i][j];
401 cmax = max(cs_dist[i][j], cmax);
404 fprintf(stderr, "cmid: %g, cmax: %g, max_size: %d\n", cmid, cmax, max_size);
406 fp = gmx_ffopen(xpm, "w");
407 write_xpm3(fp, 0, "Cluster size distribution", "# clusters", timebuf, "Size",
408 n_x, max_size, t_x, t_y, cs_dist, 0, cmid, cmax,
409 rlo, rmid, rhi, &nlevels);
413 for (i = 0; (i < n_x); i++)
415 for (j = 0; (j < max_size); j++)
417 cs_dist[i][j] *= (j+1);
418 if ((cs_dist[i][j] > 0) && (cs_dist[i][j] < cmid))
420 cmid = cs_dist[i][j];
422 cmax = max(cs_dist[i][j], cmax);
425 fprintf(stderr, "cmid: %g, cmax: %g, max_size: %d\n", cmid, cmax, max_size);
426 fp = gmx_ffopen(xpmw, "w");
427 write_xpm3(fp, 0, "Weighted cluster size distribution", "Fraction", timebuf,
428 "Size", n_x, max_size, t_x, t_y, cs_dist, 0, cmid, cmax,
429 rlo, rmid, rhi, &nlevels);
437 int gmx_clustsize(int argc, char *argv[])
439 const char *desc[] = {
440 "[THISMODULE] computes the size distributions of molecular/atomic clusters in",
441 "the gas phase. The output is given in the form of an [TT].xpm[tt] file.",
442 "The total number of clusters is written to an [TT].xvg[tt] file.[PAR]",
443 "When the [TT]-mol[tt] option is given clusters will be made out of",
444 "molecules rather than atoms, which allows clustering of large molecules.",
445 "In this case an index file would still contain atom numbers",
446 "or your calculation will die with a SEGV.[PAR]",
447 "When velocities are present in your trajectory, the temperature of",
448 "the largest cluster will be printed in a separate [TT].xvg[tt] file assuming",
449 "that the particles are free to move. If you are using constraints,",
450 "please correct the temperature. For instance water simulated with SHAKE",
451 "or SETTLE will yield a temperature that is 1.5 times too low. You can",
452 "compensate for this with the [TT]-ndf[tt] option. Remember to take the removal",
453 "of center of mass motion into account.[PAR]",
454 "The [TT]-mc[tt] option will produce an index file containing the",
455 "atom numbers of the largest cluster."
458 static real cutoff = 0.35;
459 static int nskip = 0;
460 static int nlevels = 20;
462 static gmx_bool bMol = FALSE;
463 static gmx_bool bPBC = TRUE;
464 static rvec rlo = { 1.0, 1.0, 0.0 };
465 static rvec rhi = { 0.0, 0.0, 1.0 };
470 { "-cut", FALSE, etREAL, {&cutoff},
471 "Largest distance (nm) to be considered in a cluster" },
472 { "-mol", FALSE, etBOOL, {&bMol},
473 "Cluster molecules rather than atoms (needs [TT].tpr[tt] file)" },
474 { "-pbc", FALSE, etBOOL, {&bPBC},
475 "Use periodic boundary conditions" },
476 { "-nskip", FALSE, etINT, {&nskip},
477 "Number of frames to skip between writing" },
478 { "-nlevels", FALSE, etINT, {&nlevels},
479 "Number of levels of grey in [TT].xpm[tt] output" },
480 { "-ndf", FALSE, etINT, {&ndf},
481 "Number of degrees of freedom of the entire system for temperature calculation. If not set, the number of atoms times three is used." },
482 { "-rgblo", FALSE, etRVEC, {rlo},
483 "RGB values for the color of the lowest occupied cluster size" },
484 { "-rgbhi", FALSE, etRVEC, {rhi},
485 "RGB values for the color of the highest occupied cluster size" }
487 #define NPA asize(pa)
488 const char *fnNDX, *fnTPR;
493 { efTRX, "-f", NULL, ffREAD },
494 { efTPR, NULL, NULL, ffOPTRD },
495 { efNDX, NULL, NULL, ffOPTRD },
496 { efXPM, "-o", "csize", ffWRITE },
497 { efXPM, "-ow", "csizew", ffWRITE },
498 { efXVG, "-nc", "nclust", ffWRITE },
499 { efXVG, "-mc", "maxclust", ffWRITE },
500 { efXVG, "-ac", "avclust", ffWRITE },
501 { efXVG, "-hc", "histo-clust", ffWRITE },
502 { efXVG, "-temp", "temp", ffOPTWR },
503 { efNDX, "-mcn", "maxclust", ffOPTWR }
505 #define NFILE asize(fnm)
507 if (!parse_common_args(&argc, argv,
508 PCA_CAN_VIEW | PCA_CAN_TIME | PCA_TIME_UNIT,
509 NFILE, fnm, NPA, pa, asize(desc), desc, 0, NULL, &oenv))
514 fnNDX = ftp2fn_null(efNDX, NFILE, fnm);
515 rgblo.r = rlo[XX], rgblo.g = rlo[YY], rgblo.b = rlo[ZZ];
516 rgbhi.r = rhi[XX], rgbhi.g = rhi[YY], rgbhi.b = rhi[ZZ];
518 fnTPR = ftp2fn_null(efTPR, NFILE, fnm);
521 gmx_fatal(FARGS, "You need a tpr file for the -mol option");
524 clust_size(fnNDX, ftp2fn(efTRX, NFILE, fnm), opt2fn("-o", NFILE, fnm),
525 opt2fn("-ow", NFILE, fnm),
526 opt2fn("-nc", NFILE, fnm), opt2fn("-ac", NFILE, fnm),
527 opt2fn("-mc", NFILE, fnm), opt2fn("-hc", NFILE, fnm),
528 opt2fn("-temp", NFILE, fnm), opt2fn("-mcn", NFILE, fnm),
530 cutoff, nskip, nlevels, rgblo, rgbhi, ndf, oenv);