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.
51 #include "gromacs/commandline/pargs.h"
53 #include "gromacs/fileio/futil.h"
54 #include "gromacs/fileio/tpxio.h"
55 #include "gromacs/fileio/trxio.h"
64 #include "gromacs/fileio/matio.h"
65 #include "mtop_util.h"
69 static void clust_size(const char *ndx, const char *trx, const char *xpm,
70 const char *xpmw, const char *ncl, const char *acl,
71 const char *mcl, const char *histo, const char *tempf,
72 const char *mcn, gmx_bool bMol, gmx_bool bPBC, const char *tpr,
73 real cut, int nskip, int nlevels,
74 t_rgb rmid, t_rgb rhi, int ndf,
75 const output_env_t oenv)
77 FILE *fp, *gp, *hp, *tp;
78 atom_id *index = NULL;
81 rvec *x = NULL, *v = NULL, dx;
85 gmx_bool bSame, bTPRwarn = TRUE;
89 gmx_mtop_t *mtop = NULL;
92 gmx_mtop_atomlookup_t alook;
94 int version, generation, ii, jj, nsame;
96 /* Cluster size distribution (matrix) */
97 real **cs_dist = NULL;
98 real tf, dx2, cut2, *t_x = NULL, *t_y, cmid, cmax, cav, ekin;
99 int i, j, k, ai, aj, ak, ci, cj, nframe, nclust, n_x, n_y, max_size = 0;
100 int *clust_index, *clust_size, max_clust_size, max_clust_ind, nav, nhisto;
101 t_rgb rlo = { 1.0, 1.0, 1.0 };
103 clear_trxframe(&fr, TRUE);
104 sprintf(timebuf, "Time (%s)", output_env_get_time_unit(oenv));
105 tf = output_env_get_time_factor(oenv);
106 fp = xvgropen(ncl, "Number of clusters", timebuf, "N", oenv);
107 gp = xvgropen(acl, "Average cluster size", timebuf, "#molecules", oenv);
108 hp = xvgropen(mcl, "Max cluster size", timebuf, "#molecules", oenv);
109 tp = xvgropen(tempf, "Temperature of largest cluster", timebuf, "T (K)",
112 if (!read_first_frame(oenv, &status, trx, &fr, TRX_NEED_X | TRX_READ_V))
123 read_tpxheader(tpr, &tpxh, TRUE, &version, &generation);
124 if (tpxh.natoms != natoms)
126 gmx_fatal(FARGS, "tpr (%d atoms) and trajectory (%d atoms) do not match!",
127 tpxh.natoms, natoms);
129 ePBC = read_tpx(tpr, NULL, NULL, &natoms, NULL, NULL, NULL, mtop);
137 tfac = ndf/(3.0*natoms);
144 printf("Using molecules rather than atoms. Not reading index file %s\n",
147 mols = &(mtop->mols);
149 /* Make dummy index */
152 for (i = 0; (i < nindex); i++)
156 gname = strdup("mols");
160 rd_index(ndx, 1, &nindex, &index, &gname);
163 alook = gmx_mtop_atomlookup_init(mtop);
165 snew(clust_index, nindex);
166 snew(clust_size, nindex);
171 for (i = 0; (i < nindex); i++)
179 if ((nskip == 0) || ((nskip > 0) && ((nframe % nskip) == 0)))
183 set_pbc(&pbc, ePBC, fr.box);
188 /* Put all atoms/molecules in their own cluster, with size 1 */
189 for (i = 0; (i < nindex); i++)
191 /* Cluster index is indexed with atom index number */
193 /* Cluster size is indexed with cluster number */
197 /* Loop over atoms */
198 for (i = 0; (i < nindex); i++)
203 /* Loop over atoms (only half a matrix) */
204 for (j = i+1; (j < nindex); j++)
208 /* If they are not in the same cluster already */
213 /* Compute distance */
217 for (ii = mols->index[ai]; !bSame && (ii < mols->index[ai+1]); ii++)
219 for (jj = mols->index[aj]; !bSame && (jj < mols->index[aj+1]); jj++)
223 pbc_dx(&pbc, x[ii], x[jj], dx);
227 rvec_sub(x[ii], x[jj], dx);
230 bSame = (dx2 < cut2);
238 pbc_dx(&pbc, x[ai], x[aj], dx);
242 rvec_sub(x[ai], x[aj], dx);
245 bSame = (dx2 < cut2);
247 /* If distance less than cut-off */
250 /* Merge clusters: check for all atoms whether they are in
251 * cluster cj and if so, put them in ci
253 for (k = 0; (k < nindex); k++)
255 if (clust_index[k] == cj)
257 if (clust_size[cj] <= 0)
259 gmx_fatal(FARGS, "negative cluster size %d for element %d",
273 t_x[n_x-1] = fr.time*tf;
274 srenew(cs_dist, n_x);
275 snew(cs_dist[n_x-1], nindex);
279 for (i = 0; (i < nindex); i++)
282 if (ci > max_clust_size)
290 cs_dist[n_x-1][ci-1] += 1.0;
291 max_size = max(max_size, ci);
299 fprintf(fp, "%14.6e %10d\n", fr.time, nclust);
302 fprintf(gp, "%14.6e %10.3f\n", fr.time, cav/nav);
304 fprintf(hp, "%14.6e %10d\n", fr.time, max_clust_size);
306 /* Analyse velocities, if present */
313 printf("You need a [TT].tpr[tt] file to analyse temperatures\n");
320 /* Loop over clusters and for each cluster compute 1/2 m v^2 */
321 if (max_clust_ind >= 0)
324 for (i = 0; (i < nindex); i++)
326 if (clust_index[i] == max_clust_ind)
329 gmx_mtop_atomnr_to_atom(alook, ai, &atom);
330 ekin += 0.5*atom->m*iprod(v[ai], v[ai]);
333 temp = (ekin*2.0)/(3.0*tfac*max_clust_size*BOLTZ);
334 fprintf(tp, "%10.3f %10.3f\n", fr.time, temp);
340 while (read_next_frame(oenv, status, &fr));
347 gmx_mtop_atomlookup_destroy(alook);
349 if (max_clust_ind >= 0)
351 fp = gmx_ffopen(mcn, "w");
352 fprintf(fp, "[ max_clust ]\n");
353 for (i = 0; (i < nindex); i++)
355 if (clust_index[i] == max_clust_ind)
359 for (j = mols->index[i]; (j < mols->index[i+1]); j++)
361 fprintf(fp, "%d\n", j+1);
366 fprintf(fp, "%d\n", index[i]+1);
373 /* Print the real distribution cluster-size/numer, averaged over the trajectory. */
374 fp = xvgropen(histo, "Cluster size distribution", "Cluster size", "()", oenv);
376 fprintf(fp, "%5d %8.3f\n", 0, 0.0);
377 for (j = 0; (j < max_size); j++)
380 for (i = 0; (i < n_x); i++)
382 nelem += cs_dist[i][j];
384 fprintf(fp, "%5d %8.3f\n", j+1, nelem/n_x);
385 nhisto += (int)((j+1)*nelem/n_x);
387 fprintf(fp, "%5d %8.3f\n", j+1, 0.0);
390 fprintf(stderr, "Total number of atoms in clusters = %d\n", nhisto);
392 /* Look for the smallest entry that is not zero
393 * This will make that zero is white, and not zero is coloured.
397 for (i = 0; (i < n_x); i++)
399 for (j = 0; (j < max_size); j++)
401 if ((cs_dist[i][j] > 0) && (cs_dist[i][j] < cmid))
403 cmid = cs_dist[i][j];
405 cmax = max(cs_dist[i][j], cmax);
408 fprintf(stderr, "cmid: %g, cmax: %g, max_size: %d\n", cmid, cmax, max_size);
410 fp = gmx_ffopen(xpm, "w");
411 write_xpm3(fp, 0, "Cluster size distribution", "# clusters", timebuf, "Size",
412 n_x, max_size, t_x, t_y, cs_dist, 0, cmid, cmax,
413 rlo, rmid, rhi, &nlevels);
417 for (i = 0; (i < n_x); i++)
419 for (j = 0; (j < max_size); j++)
421 cs_dist[i][j] *= (j+1);
422 if ((cs_dist[i][j] > 0) && (cs_dist[i][j] < cmid))
424 cmid = cs_dist[i][j];
426 cmax = max(cs_dist[i][j], cmax);
429 fprintf(stderr, "cmid: %g, cmax: %g, max_size: %d\n", cmid, cmax, max_size);
430 fp = gmx_ffopen(xpmw, "w");
431 write_xpm3(fp, 0, "Weighted cluster size distribution", "Fraction", timebuf,
432 "Size", n_x, max_size, t_x, t_y, cs_dist, 0, cmid, cmax,
433 rlo, rmid, rhi, &nlevels);
441 int gmx_clustsize(int argc, char *argv[])
443 const char *desc[] = {
444 "[THISMODULE] computes the size distributions of molecular/atomic clusters in",
445 "the gas phase. The output is given in the form of an [TT].xpm[tt] file.",
446 "The total number of clusters is written to an [TT].xvg[tt] file.[PAR]",
447 "When the [TT]-mol[tt] option is given clusters will be made out of",
448 "molecules rather than atoms, which allows clustering of large molecules.",
449 "In this case an index file would still contain atom numbers",
450 "or your calculation will die with a SEGV.[PAR]",
451 "When velocities are present in your trajectory, the temperature of",
452 "the largest cluster will be printed in a separate [TT].xvg[tt] file assuming",
453 "that the particles are free to move. If you are using constraints,",
454 "please correct the temperature. For instance water simulated with SHAKE",
455 "or SETTLE will yield a temperature that is 1.5 times too low. You can",
456 "compensate for this with the [TT]-ndf[tt] option. Remember to take the removal",
457 "of center of mass motion into account.[PAR]",
458 "The [TT]-mc[tt] option will produce an index file containing the",
459 "atom numbers of the largest cluster."
462 static real cutoff = 0.35;
463 static int nskip = 0;
464 static int nlevels = 20;
466 static gmx_bool bMol = FALSE;
467 static gmx_bool bPBC = TRUE;
468 static rvec rlo = { 1.0, 1.0, 0.0 };
469 static rvec rhi = { 0.0, 0.0, 1.0 };
474 { "-cut", FALSE, etREAL, {&cutoff},
475 "Largest distance (nm) to be considered in a cluster" },
476 { "-mol", FALSE, etBOOL, {&bMol},
477 "Cluster molecules rather than atoms (needs [TT].tpr[tt] file)" },
478 { "-pbc", FALSE, etBOOL, {&bPBC},
479 "Use periodic boundary conditions" },
480 { "-nskip", FALSE, etINT, {&nskip},
481 "Number of frames to skip between writing" },
482 { "-nlevels", FALSE, etINT, {&nlevels},
483 "Number of levels of grey in [TT].xpm[tt] output" },
484 { "-ndf", FALSE, etINT, {&ndf},
485 "Number of degrees of freedom of the entire system for temperature calculation. If not set, the number of atoms times three is used." },
486 { "-rgblo", FALSE, etRVEC, {rlo},
487 "RGB values for the color of the lowest occupied cluster size" },
488 { "-rgbhi", FALSE, etRVEC, {rhi},
489 "RGB values for the color of the highest occupied cluster size" }
491 #define NPA asize(pa)
492 const char *fnNDX, *fnTPR;
497 { efTRX, "-f", NULL, ffREAD },
498 { efTPR, NULL, NULL, ffOPTRD },
499 { efNDX, NULL, NULL, ffOPTRD },
500 { efXPM, "-o", "csize", ffWRITE },
501 { efXPM, "-ow", "csizew", ffWRITE },
502 { efXVG, "-nc", "nclust", ffWRITE },
503 { efXVG, "-mc", "maxclust", ffWRITE },
504 { efXVG, "-ac", "avclust", ffWRITE },
505 { efXVG, "-hc", "histo-clust", ffWRITE },
506 { efXVG, "-temp", "temp", ffOPTWR },
507 { efNDX, "-mcn", "maxclust", ffOPTWR }
509 #define NFILE asize(fnm)
511 if (!parse_common_args(&argc, argv,
512 PCA_CAN_VIEW | PCA_CAN_TIME | PCA_TIME_UNIT | PCA_BE_NICE,
513 NFILE, fnm, NPA, pa, asize(desc), desc, 0, NULL, &oenv))
518 fnNDX = ftp2fn_null(efNDX, NFILE, fnm);
519 rgblo.r = rlo[XX], rgblo.g = rlo[YY], rgblo.b = rlo[ZZ];
520 rgbhi.r = rhi[XX], rgbhi.g = rhi[YY], rgbhi.b = rhi[ZZ];
522 fnTPR = ftp2fn_null(efTPR, NFILE, fnm);
525 gmx_fatal(FARGS, "You need a tpr file for the -mol option");
528 clust_size(fnNDX, ftp2fn(efTRX, NFILE, fnm), opt2fn("-o", NFILE, fnm),
529 opt2fn("-ow", NFILE, fnm),
530 opt2fn("-nc", NFILE, fnm), opt2fn("-ac", NFILE, fnm),
531 opt2fn("-mc", NFILE, fnm), opt2fn("-hc", NFILE, fnm),
532 opt2fn("-temp", NFILE, fnm), opt2fn("-mcn", NFILE, fnm),
534 cutoff, nskip, nlevels, rgblo, rgbhi, ndf, oenv);