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.
45 #include "gromacs/utility/smalloc.h"
47 #include "gromacs/commandline/pargs.h"
48 #include "gromacs/math/utilities.h"
49 #include "gromacs/utility/futil.h"
52 #include "gromacs/fileio/xvgr.h"
55 #include "gromacs/fileio/tpxio.h"
56 #include "gromacs/fileio/trxio.h"
58 #include "gromacs/fileio/matio.h"
62 int gmx_vanhove(int argc, char *argv[])
64 const char *desc[] = {
65 "[THISMODULE] computes the Van Hove correlation function.",
66 "The Van Hove G(r,t) is the probability that a particle that is at r[SUB]0[sub]",
67 "at time zero can be found at position r[SUB]0[sub]+r at time t.",
68 "[THISMODULE] determines G not for a vector r, but for the length of r.",
69 "Thus it gives the probability that a particle moves a distance of r",
71 "Jumps across the periodic boundaries are removed.",
72 "Corrections are made for scaling due to isotropic",
73 "or anisotropic pressure coupling.",
75 "With option [TT]-om[tt] the whole matrix can be written as a function",
76 "of t and r or as a function of [SQRT]t[sqrt] and r (option [TT]-sqrt[tt]).",
78 "With option [TT]-or[tt] the Van Hove function is plotted for one",
79 "or more values of t. Option [TT]-nr[tt] sets the number of times,",
80 "option [TT]-fr[tt] the number spacing between the times.",
81 "The binwidth is set with option [TT]-rbin[tt]. The number of bins",
82 "is determined automatically.",
84 "With option [TT]-ot[tt] the integral up to a certain distance",
85 "(option [TT]-rt[tt]) is plotted as a function of time.",
87 "For all frames that are read the coordinates of the selected particles",
88 "are stored in memory. Therefore the program may use a lot of memory.",
89 "For options [TT]-om[tt] and [TT]-ot[tt] the program may be slow.",
90 "This is because the calculation scales as the number of frames times",
91 "[TT]-fm[tt] or [TT]-ft[tt].",
92 "Note that with the [TT]-dt[tt] option the memory usage and calculation",
93 "time can be reduced."
95 static int fmmax = 0, ftmax = 0, nlev = 81, nr = 1, fshift = 0;
96 static real sbin = 0, rmax = 2, rbin = 0.01, mmax = 0, rint = 0;
98 { "-sqrt", FALSE, etREAL, {&sbin},
99 "Use [SQRT]t[sqrt] on the matrix axis which binspacing # in [SQRT]ps[sqrt]" },
100 { "-fm", FALSE, etINT, {&fmmax},
101 "Number of frames in the matrix, 0 is plot all" },
102 { "-rmax", FALSE, etREAL, {&rmax},
103 "Maximum r in the matrix (nm)" },
104 { "-rbin", FALSE, etREAL, {&rbin},
105 "Binwidth in the matrix and for [TT]-or[tt] (nm)" },
106 { "-mmax", FALSE, etREAL, {&mmax},
107 "Maximum density in the matrix, 0 is calculate (1/nm)" },
108 { "-nlevels", FALSE, etINT, {&nlev},
109 "Number of levels in the matrix" },
110 { "-nr", FALSE, etINT, {&nr},
111 "Number of curves for the [TT]-or[tt] output" },
112 { "-fr", FALSE, etINT, {&fshift},
113 "Frame spacing for the [TT]-or[tt] output" },
114 { "-rt", FALSE, etREAL, {&rint},
115 "Integration limit for the [TT]-ot[tt] output (nm)" },
116 { "-ft", FALSE, etINT, {&ftmax},
117 "Number of frames in the [TT]-ot[tt] output, 0 is plot all" }
119 #define NPA asize(pa)
122 { efTRX, NULL, NULL, ffREAD },
123 { efTPS, NULL, NULL, ffREAD },
124 { efNDX, NULL, NULL, ffOPTRD },
125 { efXPM, "-om", "vanhove", ffOPTWR },
126 { efXVG, "-or", "vanhove_r", ffOPTWR },
127 { efXVG, "-ot", "vanhove_t", ffOPTWR }
129 #define NFILE asize(fnm)
132 const char *matfile, *otfile, *orfile;
136 matrix boxtop, box, *sbox, avbox, corr;
137 rvec *xtop, *x, **sx;
138 int isize, nalloc, nallocn, natom;
142 int nfr, f, ff, i, m, mat_nx = 0, nbin = 0, bin, mbin, fbin;
143 real *time, t, invbin = 0, rmax2 = 0, rint2 = 0, d2;
144 real invsbin = 0, matmax, normfac, dt, *tickx, *ticky;
145 char buf[STRLEN], **legend;
147 int *pt = NULL, **pr = NULL, *mcount = NULL, *tcount = NULL, *rcount = NULL;
149 t_rgb rlo = {1, 1, 1}, rhi = {0, 0, 0};
151 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
152 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
157 matfile = opt2fn_null("-om", NFILE, fnm);
158 if (opt2parg_bSet("-fr", NPA, pa))
160 orfile = opt2fn("-or", NFILE, fnm);
164 orfile = opt2fn_null("-or", NFILE, fnm);
166 if (opt2parg_bSet("-rt", NPA, pa))
168 otfile = opt2fn("-ot", NFILE, fnm);
172 otfile = opt2fn_null("-ot", NFILE, fnm);
175 if (!matfile && !otfile && !orfile)
178 "For output set one (or more) of the output file options\n");
182 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, boxtop,
184 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
192 natom = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
199 srenew(time, nalloc);
200 srenew(sbox, nalloc);
205 copy_mat(box, sbox[nfr]);
206 /* This assumes that the off-diagonal box elements
207 * are not affected by jumps across the periodic boundaries.
209 m_add(avbox, box, avbox);
210 snew(sx[nfr], isize);
211 for (i = 0; i < isize; i++)
213 copy_rvec(x[index[i]], sx[nfr][i]);
218 while (read_next_x(oenv, status, &t, x, box));
224 fprintf(stderr, "Read %d frames\n", nfr);
226 dt = (time[nfr-1] - time[0])/(nfr - 1);
227 /* Some ugly rounding to get nice nice times in the output */
228 dt = (int)(10000.0*dt + 0.5)/10000.0;
234 if (fmmax <= 0 || fmmax >= nfr)
239 nbin = (int)(rmax*invbin + 0.5);
247 mat_nx = sqrt(fmmax*dt)*invsbin + 1;
250 for (f = 0; f < mat_nx; f++)
254 rmax2 = sqr(nbin*rbin);
255 /* Initialize time zero */
256 mat[0][0] = nfr*isize;
280 /* Initialize time zero */
289 msmul(avbox, 1.0/nfr, avbox);
290 for (f = 0; f < nfr; f++)
294 fprintf(stderr, "\rProcessing frame %d", f);
296 /* Scale all the configuration to the average box */
297 m_inv_ur0(sbox[f], corr);
298 mmul_ur0(avbox, corr, corr);
299 for (i = 0; i < isize; i++)
301 mvmul_ur0(corr, sx[f][i], sx[f][i]);
304 /* Correct for periodic jumps */
305 for (m = DIM-1; m >= 0; m--)
307 while (sx[f][i][m] - sx[f-1][i][m] > 0.5*avbox[m][m])
309 rvec_dec(sx[f][i], avbox[m]);
311 while (sx[f][i][m] - sx[f-1][i][m] <= -0.5*avbox[m][m])
313 rvec_inc(sx[f][i], avbox[m]);
318 for (ff = 0; ff < f; ff++)
321 if (fbin <= fmmax || fbin <= ftmax)
329 mbin = (int)(sqrt(fbin*dt)*invsbin + 0.5);
331 for (i = 0; i < isize; i++)
333 d2 = distance2(sx[f][i], sx[ff][i]);
334 if (mbin < mat_nx && d2 < rmax2)
336 bin = (int)(sqrt(d2)*invbin + 0.5);
342 if (fbin <= ftmax && d2 <= rint2)
359 for (fbin = 0; fbin < nr; fbin++)
361 ff = f - (fbin + 1)*fshift;
364 for (i = 0; i < isize; i++)
366 d2 = distance2(sx[f][i], sx[ff][i]);
367 bin = (int)(sqrt(d2)*invbin + 0.5);
370 nallocn = 10*(bin/10) + 11;
371 for (m = 0; m < nr; m++)
373 srenew(pr[m], nallocn);
374 for (i = nalloc; i < nallocn; i++)
388 fprintf(stderr, "\n");
393 for (f = 0; f < mat_nx; f++)
395 normfac = 1.0/(mcount[f]*isize*rbin);
396 for (i = 0; i < nbin; i++)
398 mat[f][i] *= normfac;
399 if (mat[f][i] > matmax && (f != 0 || i != 0))
405 fprintf(stdout, "Value at (0,0): %.3f, maximum of the rest %.3f\n",
412 for (f = 0; f < mat_nx; f++)
424 for (i = 0; i <= nbin; i++)
428 fp = gmx_ffopen(matfile, "w");
429 write_xpm(fp, MAT_SPATIAL_Y, "Van Hove function", "G (1/nm)",
430 sbin == 0 ? "time (ps)" : "sqrt(time) (ps^1/2)", "r (nm)",
431 mat_nx, nbin, tickx, ticky, mat, 0, matmax, rlo, rhi, &nlev);
437 fp = xvgropen(orfile, "Van Hove function", "r (nm)", "G (nm\\S-1\\N)", oenv);
438 fprintf(fp, "@ subtitle \"for particles in group %s\"\n", grpname);
440 for (fbin = 0; fbin < nr; fbin++)
442 sprintf(buf, "%g ps", (fbin + 1)*fshift*dt);
443 legend[fbin] = strdup(buf);
445 xvgr_legend(fp, nr, (const char**)legend, oenv);
446 for (i = 0; i < nalloc; i++)
448 fprintf(fp, "%g", i*rbin);
449 for (fbin = 0; fbin < nr; fbin++)
452 (real)pr[fbin][i]/(rcount[fbin]*isize*rbin*(i == 0 ? 0.5 : 1)));
461 sprintf(buf, "Probability of moving less than %g nm", rint);
462 fp = xvgropen(otfile, buf, "t (ps)", "", oenv);
463 fprintf(fp, "@ subtitle \"for particles in group %s\"\n", grpname);
464 for (f = 0; f <= ftmax; f++)
466 fprintf(fp, "%g %g\n", f*dt, (real)pt[f]/(tcount[f]*isize));
471 do_view(oenv, matfile, NULL);
472 do_view(oenv, orfile, NULL);
473 do_view(oenv, otfile, NULL);