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.
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/coulomb.h"
49 #include "gromacs/legacyheaders/macros.h"
50 #include "gromacs/legacyheaders/nrnb.h"
51 #include "gromacs/legacyheaders/pme.h"
52 #include "gromacs/legacyheaders/typedefs.h"
53 #include "gromacs/math/units.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/pbcutil/pbc.h"
56 #include "gromacs/topology/index.h"
57 #include "gromacs/topology/mtop_util.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/futil.h"
60 #include "gromacs/utility/smalloc.h"
62 static void clust_size(const char *ndx, const char *trx, const char *xpm,
63 const char *xpmw, const char *ncl, const char *acl,
64 const char *mcl, const char *histo, const char *tempf,
65 const char *mcn, gmx_bool bMol, gmx_bool bPBC, const char *tpr,
66 real cut, int nskip, int nlevels,
67 t_rgb rmid, t_rgb rhi, int ndf,
68 const output_env_t oenv)
70 FILE *fp, *gp, *hp, *tp;
71 atom_id *index = NULL;
74 rvec *x = NULL, *v = NULL, dx;
78 gmx_bool bSame, bTPRwarn = TRUE;
82 gmx_mtop_t *mtop = NULL;
85 gmx_mtop_atomlookup_t alook;
87 int version, generation, ii, jj, nsame;
89 /* Cluster size distribution (matrix) */
90 real **cs_dist = NULL;
91 real tf, dx2, cut2, *t_x = NULL, *t_y, cmid, cmax, cav, ekin;
92 int i, j, k, ai, aj, ak, ci, cj, nframe, nclust, n_x, n_y, max_size = 0;
93 int *clust_index, *clust_size, max_clust_size, max_clust_ind, nav, nhisto;
94 t_rgb rlo = { 1.0, 1.0, 1.0 };
96 clear_trxframe(&fr, TRUE);
97 sprintf(timebuf, "Time (%s)", output_env_get_time_unit(oenv));
98 tf = output_env_get_time_factor(oenv);
99 fp = xvgropen(ncl, "Number of clusters", timebuf, "N", oenv);
100 gp = xvgropen(acl, "Average cluster size", timebuf, "#molecules", oenv);
101 hp = xvgropen(mcl, "Max cluster size", timebuf, "#molecules", oenv);
102 tp = xvgropen(tempf, "Temperature of largest cluster", timebuf, "T (K)",
105 if (!read_first_frame(oenv, &status, trx, &fr, TRX_NEED_X | TRX_READ_V))
116 read_tpxheader(tpr, &tpxh, TRUE, &version, &generation);
117 if (tpxh.natoms != natoms)
119 gmx_fatal(FARGS, "tpr (%d atoms) and trajectory (%d atoms) do not match!",
120 tpxh.natoms, natoms);
122 ePBC = read_tpx(tpr, NULL, NULL, &natoms, NULL, NULL, NULL, mtop);
130 tfac = ndf/(3.0*natoms);
137 printf("Using molecules rather than atoms. Not reading index file %s\n",
140 mols = &(mtop->mols);
142 /* Make dummy index */
145 for (i = 0; (i < nindex); i++)
149 gname = gmx_strdup("mols");
153 rd_index(ndx, 1, &nindex, &index, &gname);
156 alook = gmx_mtop_atomlookup_init(mtop);
158 snew(clust_index, nindex);
159 snew(clust_size, nindex);
164 for (i = 0; (i < nindex); i++)
172 if ((nskip == 0) || ((nskip > 0) && ((nframe % nskip) == 0)))
176 set_pbc(&pbc, ePBC, fr.box);
181 /* Put all atoms/molecules in their own cluster, with size 1 */
182 for (i = 0; (i < nindex); i++)
184 /* Cluster index is indexed with atom index number */
186 /* Cluster size is indexed with cluster number */
190 /* Loop over atoms */
191 for (i = 0; (i < nindex); i++)
196 /* Loop over atoms (only half a matrix) */
197 for (j = i+1; (j < nindex); j++)
201 /* If they are not in the same cluster already */
206 /* Compute distance */
210 for (ii = mols->index[ai]; !bSame && (ii < mols->index[ai+1]); ii++)
212 for (jj = mols->index[aj]; !bSame && (jj < mols->index[aj+1]); jj++)
216 pbc_dx(&pbc, x[ii], x[jj], dx);
220 rvec_sub(x[ii], x[jj], dx);
223 bSame = (dx2 < cut2);
231 pbc_dx(&pbc, x[ai], x[aj], dx);
235 rvec_sub(x[ai], x[aj], dx);
238 bSame = (dx2 < cut2);
240 /* If distance less than cut-off */
243 /* Merge clusters: check for all atoms whether they are in
244 * cluster cj and if so, put them in ci
246 for (k = 0; (k < nindex); k++)
248 if (clust_index[k] == cj)
250 if (clust_size[cj] <= 0)
252 gmx_fatal(FARGS, "negative cluster size %d for element %d",
266 t_x[n_x-1] = fr.time*tf;
267 srenew(cs_dist, n_x);
268 snew(cs_dist[n_x-1], nindex);
272 for (i = 0; (i < nindex); i++)
275 if (ci > max_clust_size)
283 cs_dist[n_x-1][ci-1] += 1.0;
284 max_size = max(max_size, ci);
292 fprintf(fp, "%14.6e %10d\n", fr.time, nclust);
295 fprintf(gp, "%14.6e %10.3f\n", fr.time, cav/nav);
297 fprintf(hp, "%14.6e %10d\n", fr.time, max_clust_size);
299 /* Analyse velocities, if present */
306 printf("You need a [TT].tpr[tt] file to analyse temperatures\n");
313 /* Loop over clusters and for each cluster compute 1/2 m v^2 */
314 if (max_clust_ind >= 0)
317 for (i = 0; (i < nindex); i++)
319 if (clust_index[i] == max_clust_ind)
322 gmx_mtop_atomnr_to_atom(alook, ai, &atom);
323 ekin += 0.5*atom->m*iprod(v[ai], v[ai]);
326 temp = (ekin*2.0)/(3.0*tfac*max_clust_size*BOLTZ);
327 fprintf(tp, "%10.3f %10.3f\n", fr.time, temp);
333 while (read_next_frame(oenv, status, &fr));
340 gmx_mtop_atomlookup_destroy(alook);
342 if (max_clust_ind >= 0)
344 fp = gmx_ffopen(mcn, "w");
345 fprintf(fp, "[ max_clust ]\n");
346 for (i = 0; (i < nindex); i++)
348 if (clust_index[i] == max_clust_ind)
352 for (j = mols->index[i]; (j < mols->index[i+1]); j++)
354 fprintf(fp, "%d\n", j+1);
359 fprintf(fp, "%d\n", index[i]+1);
366 /* Print the real distribution cluster-size/numer, averaged over the trajectory. */
367 fp = xvgropen(histo, "Cluster size distribution", "Cluster size", "()", oenv);
369 fprintf(fp, "%5d %8.3f\n", 0, 0.0);
370 for (j = 0; (j < max_size); j++)
373 for (i = 0; (i < n_x); i++)
375 nelem += cs_dist[i][j];
377 fprintf(fp, "%5d %8.3f\n", j+1, nelem/n_x);
378 nhisto += (int)((j+1)*nelem/n_x);
380 fprintf(fp, "%5d %8.3f\n", j+1, 0.0);
383 fprintf(stderr, "Total number of atoms in clusters = %d\n", nhisto);
385 /* Look for the smallest entry that is not zero
386 * This will make that zero is white, and not zero is coloured.
390 for (i = 0; (i < n_x); i++)
392 for (j = 0; (j < max_size); j++)
394 if ((cs_dist[i][j] > 0) && (cs_dist[i][j] < cmid))
396 cmid = cs_dist[i][j];
398 cmax = max(cs_dist[i][j], cmax);
401 fprintf(stderr, "cmid: %g, cmax: %g, max_size: %d\n", cmid, cmax, max_size);
403 fp = gmx_ffopen(xpm, "w");
404 write_xpm3(fp, 0, "Cluster size distribution", "# clusters", timebuf, "Size",
405 n_x, max_size, t_x, t_y, cs_dist, 0, cmid, cmax,
406 rlo, rmid, rhi, &nlevels);
410 for (i = 0; (i < n_x); i++)
412 for (j = 0; (j < max_size); j++)
414 cs_dist[i][j] *= (j+1);
415 if ((cs_dist[i][j] > 0) && (cs_dist[i][j] < cmid))
417 cmid = cs_dist[i][j];
419 cmax = max(cs_dist[i][j], cmax);
422 fprintf(stderr, "cmid: %g, cmax: %g, max_size: %d\n", cmid, cmax, max_size);
423 fp = gmx_ffopen(xpmw, "w");
424 write_xpm3(fp, 0, "Weighted cluster size distribution", "Fraction", timebuf,
425 "Size", n_x, max_size, t_x, t_y, cs_dist, 0, cmid, cmax,
426 rlo, rmid, rhi, &nlevels);
434 int gmx_clustsize(int argc, char *argv[])
436 const char *desc[] = {
437 "[THISMODULE] computes the size distributions of molecular/atomic clusters in",
438 "the gas phase. The output is given in the form of an [TT].xpm[tt] file.",
439 "The total number of clusters is written to an [TT].xvg[tt] file.[PAR]",
440 "When the [TT]-mol[tt] option is given clusters will be made out of",
441 "molecules rather than atoms, which allows clustering of large molecules.",
442 "In this case an index file would still contain atom numbers",
443 "or your calculation will die with a SEGV.[PAR]",
444 "When velocities are present in your trajectory, the temperature of",
445 "the largest cluster will be printed in a separate [TT].xvg[tt] file assuming",
446 "that the particles are free to move. If you are using constraints,",
447 "please correct the temperature. For instance water simulated with SHAKE",
448 "or SETTLE will yield a temperature that is 1.5 times too low. You can",
449 "compensate for this with the [TT]-ndf[tt] option. Remember to take the removal",
450 "of center of mass motion into account.[PAR]",
451 "The [TT]-mc[tt] option will produce an index file containing the",
452 "atom numbers of the largest cluster."
455 static real cutoff = 0.35;
456 static int nskip = 0;
457 static int nlevels = 20;
459 static gmx_bool bMol = FALSE;
460 static gmx_bool bPBC = TRUE;
461 static rvec rlo = { 1.0, 1.0, 0.0 };
462 static rvec rhi = { 0.0, 0.0, 1.0 };
467 { "-cut", FALSE, etREAL, {&cutoff},
468 "Largest distance (nm) to be considered in a cluster" },
469 { "-mol", FALSE, etBOOL, {&bMol},
470 "Cluster molecules rather than atoms (needs [TT].tpr[tt] file)" },
471 { "-pbc", FALSE, etBOOL, {&bPBC},
472 "Use periodic boundary conditions" },
473 { "-nskip", FALSE, etINT, {&nskip},
474 "Number of frames to skip between writing" },
475 { "-nlevels", FALSE, etINT, {&nlevels},
476 "Number of levels of grey in [TT].xpm[tt] output" },
477 { "-ndf", FALSE, etINT, {&ndf},
478 "Number of degrees of freedom of the entire system for temperature calculation. If not set, the number of atoms times three is used." },
479 { "-rgblo", FALSE, etRVEC, {rlo},
480 "RGB values for the color of the lowest occupied cluster size" },
481 { "-rgbhi", FALSE, etRVEC, {rhi},
482 "RGB values for the color of the highest occupied cluster size" }
484 #define NPA asize(pa)
485 const char *fnNDX, *fnTPR;
490 { efTRX, "-f", NULL, ffREAD },
491 { efTPR, NULL, NULL, ffOPTRD },
492 { efNDX, NULL, NULL, ffOPTRD },
493 { efXPM, "-o", "csize", ffWRITE },
494 { efXPM, "-ow", "csizew", ffWRITE },
495 { efXVG, "-nc", "nclust", ffWRITE },
496 { efXVG, "-mc", "maxclust", ffWRITE },
497 { efXVG, "-ac", "avclust", ffWRITE },
498 { efXVG, "-hc", "histo-clust", ffWRITE },
499 { efXVG, "-temp", "temp", ffOPTWR },
500 { efNDX, "-mcn", "maxclust", ffOPTWR }
502 #define NFILE asize(fnm)
504 if (!parse_common_args(&argc, argv,
505 PCA_CAN_VIEW | PCA_CAN_TIME | PCA_TIME_UNIT,
506 NFILE, fnm, NPA, pa, asize(desc), desc, 0, NULL, &oenv))
511 fnNDX = ftp2fn_null(efNDX, NFILE, fnm);
512 rgblo.r = rlo[XX], rgblo.g = rlo[YY], rgblo.b = rlo[ZZ];
513 rgbhi.r = rhi[XX], rgbhi.g = rhi[YY], rgbhi.b = rhi[ZZ];
515 fnTPR = ftp2fn_null(efTPR, NFILE, fnm);
518 gmx_fatal(FARGS, "You need a tpr file for the -mol option");
521 clust_size(fnNDX, ftp2fn(efTRX, NFILE, fnm), opt2fn("-o", NFILE, fnm),
522 opt2fn("-ow", NFILE, fnm),
523 opt2fn("-nc", NFILE, fnm), opt2fn("-ac", NFILE, fnm),
524 opt2fn("-mc", NFILE, fnm), opt2fn("-hc", NFILE, fnm),
525 opt2fn("-temp", NFILE, fnm), opt2fn("-mcn", NFILE, fnm),
527 cutoff, nskip, nlevels, rgblo, rgbhi, ndf, oenv);