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) 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.
46 #include "gromacs/utility/smalloc.h"
47 #include "gromacs/legacyheaders/macros.h"
48 #include "gromacs/commandline/pargs.h"
49 #include "gromacs/math/utilities.h"
50 #include "gromacs/utility/futil.h"
51 #include "gromacs/topology/index.h"
52 #include "gromacs/legacyheaders/typedefs.h"
53 #include "gromacs/fileio/xvgr.h"
54 #include "gromacs/legacyheaders/viewit.h"
56 #include "gromacs/fileio/tpxio.h"
57 #include "gromacs/fileio/trxio.h"
58 #include "gromacs/math/vec.h"
59 #include "gromacs/fileio/matio.h"
63 int gmx_vanhove(int argc, char *argv[])
65 const char *desc[] = {
66 "[THISMODULE] computes the Van Hove correlation function.",
67 "The Van Hove G(r,t) is the probability that a particle that is at r[SUB]0[sub]",
68 "at time zero can be found at position r[SUB]0[sub]+r at time t.",
69 "[THISMODULE] determines G not for a vector r, but for the length of r.",
70 "Thus it gives the probability that a particle moves a distance of r",
72 "Jumps across the periodic boundaries are removed.",
73 "Corrections are made for scaling due to isotropic",
74 "or anisotropic pressure coupling.",
76 "With option [TT]-om[tt] the whole matrix can be written as a function",
77 "of t and r or as a function of [SQRT]t[sqrt] and r (option [TT]-sqrt[tt]).",
79 "With option [TT]-or[tt] the Van Hove function is plotted for one",
80 "or more values of t. Option [TT]-nr[tt] sets the number of times,",
81 "option [TT]-fr[tt] the number spacing between the times.",
82 "The binwidth is set with option [TT]-rbin[tt]. The number of bins",
83 "is determined automatically.",
85 "With option [TT]-ot[tt] the integral up to a certain distance",
86 "(option [TT]-rt[tt]) is plotted as a function of time.",
88 "For all frames that are read the coordinates of the selected particles",
89 "are stored in memory. Therefore the program may use a lot of memory.",
90 "For options [TT]-om[tt] and [TT]-ot[tt] the program may be slow.",
91 "This is because the calculation scales as the number of frames times",
92 "[TT]-fm[tt] or [TT]-ft[tt].",
93 "Note that with the [TT]-dt[tt] option the memory usage and calculation",
94 "time can be reduced."
96 static int fmmax = 0, ftmax = 0, nlev = 81, nr = 1, fshift = 0;
97 static real sbin = 0, rmax = 2, rbin = 0.01, mmax = 0, rint = 0;
99 { "-sqrt", FALSE, etREAL, {&sbin},
100 "Use [SQRT]t[sqrt] on the matrix axis which binspacing # in [SQRT]ps[sqrt]" },
101 { "-fm", FALSE, etINT, {&fmmax},
102 "Number of frames in the matrix, 0 is plot all" },
103 { "-rmax", FALSE, etREAL, {&rmax},
104 "Maximum r in the matrix (nm)" },
105 { "-rbin", FALSE, etREAL, {&rbin},
106 "Binwidth in the matrix and for [TT]-or[tt] (nm)" },
107 { "-mmax", FALSE, etREAL, {&mmax},
108 "Maximum density in the matrix, 0 is calculate (1/nm)" },
109 { "-nlevels", FALSE, etINT, {&nlev},
110 "Number of levels in the matrix" },
111 { "-nr", FALSE, etINT, {&nr},
112 "Number of curves for the [TT]-or[tt] output" },
113 { "-fr", FALSE, etINT, {&fshift},
114 "Frame spacing for the [TT]-or[tt] output" },
115 { "-rt", FALSE, etREAL, {&rint},
116 "Integration limit for the [TT]-ot[tt] output (nm)" },
117 { "-ft", FALSE, etINT, {&ftmax},
118 "Number of frames in the [TT]-ot[tt] output, 0 is plot all" }
120 #define NPA asize(pa)
123 { efTRX, NULL, NULL, ffREAD },
124 { efTPS, NULL, NULL, ffREAD },
125 { efNDX, NULL, NULL, ffOPTRD },
126 { efXPM, "-om", "vanhove", ffOPTWR },
127 { efXVG, "-or", "vanhove_r", ffOPTWR },
128 { efXVG, "-ot", "vanhove_t", ffOPTWR }
130 #define NFILE asize(fnm)
133 const char *matfile, *otfile, *orfile;
137 matrix boxtop, box, *sbox, avbox, corr;
138 rvec *xtop, *x, **sx;
139 int isize, nalloc, nallocn, natom;
143 int nfr, f, ff, i, m, mat_nx = 0, nbin = 0, bin, mbin, fbin;
144 real *time, t, invbin = 0, rmax2 = 0, rint2 = 0, d2;
145 real invsbin = 0, matmax, normfac, dt, *tickx, *ticky;
146 char buf[STRLEN], **legend;
148 int *pt = NULL, **pr = NULL, *mcount = NULL, *tcount = NULL, *rcount = NULL;
150 t_rgb rlo = {1, 1, 1}, rhi = {0, 0, 0};
152 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
153 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
158 matfile = opt2fn_null("-om", NFILE, fnm);
159 if (opt2parg_bSet("-fr", NPA, pa))
161 orfile = opt2fn("-or", NFILE, fnm);
165 orfile = opt2fn_null("-or", NFILE, fnm);
167 if (opt2parg_bSet("-rt", NPA, pa))
169 otfile = opt2fn("-ot", NFILE, fnm);
173 otfile = opt2fn_null("-ot", NFILE, fnm);
176 if (!matfile && !otfile && !orfile)
179 "For output set one (or more) of the output file options\n");
183 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, boxtop,
185 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
193 natom = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
200 srenew(time, nalloc);
201 srenew(sbox, nalloc);
204 assert(time != NULL); assert(sbox != NULL);
207 copy_mat(box, sbox[nfr]);
208 /* This assumes that the off-diagonal box elements
209 * are not affected by jumps across the periodic boundaries.
211 m_add(avbox, box, avbox);
212 snew(sx[nfr], isize);
213 for (i = 0; i < isize; i++)
215 copy_rvec(x[index[i]], sx[nfr][i]);
220 while (read_next_x(oenv, status, &t, x, box));
226 fprintf(stderr, "Read %d frames\n", nfr);
228 dt = (time[nfr-1] - time[0])/(nfr - 1);
229 /* Some ugly rounding to get nice nice times in the output */
230 dt = (int)(10000.0*dt + 0.5)/10000.0;
236 if (fmmax <= 0 || fmmax >= nfr)
241 nbin = (int)(rmax*invbin + 0.5);
249 mat_nx = sqrt(fmmax*dt)*invsbin + 1;
252 for (f = 0; f < mat_nx; f++)
256 rmax2 = sqr(nbin*rbin);
257 /* Initialize time zero */
258 mat[0][0] = nfr*isize;
282 /* Initialize time zero */
291 msmul(avbox, 1.0/nfr, avbox);
292 for (f = 0; f < nfr; f++)
296 fprintf(stderr, "\rProcessing frame %d", f);
298 /* Scale all the configuration to the average box */
299 m_inv_ur0(sbox[f], corr);
300 mmul_ur0(avbox, corr, corr);
301 for (i = 0; i < isize; i++)
303 mvmul_ur0(corr, sx[f][i], sx[f][i]);
306 /* Correct for periodic jumps */
307 for (m = DIM-1; m >= 0; m--)
309 while (sx[f][i][m] - sx[f-1][i][m] > 0.5*avbox[m][m])
311 rvec_dec(sx[f][i], avbox[m]);
313 while (sx[f][i][m] - sx[f-1][i][m] <= -0.5*avbox[m][m])
315 rvec_inc(sx[f][i], avbox[m]);
320 for (ff = 0; ff < f; ff++)
323 if (fbin <= fmmax || fbin <= ftmax)
331 mbin = (int)(sqrt(fbin*dt)*invsbin + 0.5);
333 for (i = 0; i < isize; i++)
335 d2 = distance2(sx[f][i], sx[ff][i]);
336 if (mbin < mat_nx && d2 < rmax2)
338 bin = (int)(sqrt(d2)*invbin + 0.5);
344 if (fbin <= ftmax && d2 <= rint2)
361 for (fbin = 0; fbin < nr; fbin++)
363 ff = f - (fbin + 1)*fshift;
366 for (i = 0; i < isize; i++)
368 d2 = distance2(sx[f][i], sx[ff][i]);
369 bin = (int)(sqrt(d2)*invbin + 0.5);
372 nallocn = 10*(bin/10) + 11;
373 for (m = 0; m < nr; m++)
375 srenew(pr[m], nallocn);
376 for (i = nalloc; i < nallocn; i++)
390 fprintf(stderr, "\n");
395 for (f = 0; f < mat_nx; f++)
397 normfac = 1.0/(mcount[f]*isize*rbin);
398 for (i = 0; i < nbin; i++)
400 mat[f][i] *= normfac;
401 if (mat[f][i] > matmax && (f != 0 || i != 0))
407 fprintf(stdout, "Value at (0,0): %.3f, maximum of the rest %.3f\n",
414 for (f = 0; f < mat_nx; f++)
426 for (i = 0; i <= nbin; i++)
430 fp = gmx_ffopen(matfile, "w");
431 write_xpm(fp, MAT_SPATIAL_Y, "Van Hove function", "G (1/nm)",
432 sbin == 0 ? "time (ps)" : "sqrt(time) (ps^1/2)", "r (nm)",
433 mat_nx, nbin, tickx, ticky, mat, 0, matmax, rlo, rhi, &nlev);
439 fp = xvgropen(orfile, "Van Hove function", "r (nm)", "G (nm\\S-1\\N)", oenv);
440 if (output_env_get_print_xvgr_codes(oenv))
442 fprintf(fp, "@ subtitle \"for particles in group %s\"\n", grpname);
445 for (fbin = 0; fbin < nr; fbin++)
447 sprintf(buf, "%g ps", (fbin + 1)*fshift*dt);
448 legend[fbin] = gmx_strdup(buf);
450 xvgr_legend(fp, nr, (const char**)legend, oenv);
451 for (i = 0; i < nalloc; i++)
453 fprintf(fp, "%g", i*rbin);
454 for (fbin = 0; fbin < nr; fbin++)
457 (real)pr[fbin][i]/(rcount[fbin]*isize*rbin*(i == 0 ? 0.5 : 1)));
466 sprintf(buf, "Probability of moving less than %g nm", rint);
467 fp = xvgropen(otfile, buf, "t (ps)", "", oenv);
468 if (output_env_get_print_xvgr_codes(oenv))
470 fprintf(fp, "@ subtitle \"for particles in group %s\"\n", grpname);
472 for (f = 0; f <= ftmax; f++)
474 fprintf(fp, "%g %g\n", f*dt, (real)pt[f]/(tcount[f]*isize));
479 do_view(oenv, matfile, NULL);
480 do_view(oenv, orfile, NULL);
481 do_view(oenv, otfile, NULL);