3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
60 int gmx_vanhove(int argc, char *argv[])
62 const char *desc[] = {
63 "[TT]g_vanhove[tt] computes the Van Hove correlation function.",
64 "The Van Hove G(r,t) is the probability that a particle that is at r[SUB]0[sub]",
65 "at time zero can be found at position r[SUB]0[sub]+r at time t.",
66 "[TT]g_vanhove[tt] determines G not for a vector r, but for the length of r.",
67 "Thus it gives the probability that a particle moves a distance of r",
69 "Jumps across the periodic boundaries are removed.",
70 "Corrections are made for scaling due to isotropic",
71 "or anisotropic pressure coupling.",
73 "With option [TT]-om[tt] the whole matrix can be written as a function",
74 "of t and r or as a function of [SQRT]t[sqrt] and r (option [TT]-sqrt[tt]).",
76 "With option [TT]-or[tt] the Van Hove function is plotted for one",
77 "or more values of t. Option [TT]-nr[tt] sets the number of times,",
78 "option [TT]-fr[tt] the number spacing between the times.",
79 "The binwidth is set with option [TT]-rbin[tt]. The number of bins",
80 "is determined automatically.",
82 "With option [TT]-ot[tt] the integral up to a certain distance",
83 "(option [TT]-rt[tt]) is plotted as a function of time.",
85 "For all frames that are read the coordinates of the selected particles",
86 "are stored in memory. Therefore the program may use a lot of memory.",
87 "For options [TT]-om[tt] and [TT]-ot[tt] the program may be slow.",
88 "This is because the calculation scales as the number of frames times",
89 "[TT]-fm[tt] or [TT]-ft[tt].",
90 "Note that with the [TT]-dt[tt] option the memory usage and calculation",
91 "time can be reduced."
93 static int fmmax = 0, ftmax = 0, nlev = 81, nr = 1, fshift = 0;
94 static real sbin = 0, rmax = 2, rbin = 0.01, mmax = 0, rint = 0;
96 { "-sqrt", FALSE, etREAL, {&sbin},
97 "Use [SQRT]t[sqrt] on the matrix axis which binspacing # in [SQRT]ps[sqrt]" },
98 { "-fm", FALSE, etINT, {&fmmax},
99 "Number of frames in the matrix, 0 is plot all" },
100 { "-rmax", FALSE, etREAL, {&rmax},
101 "Maximum r in the matrix (nm)" },
102 { "-rbin", FALSE, etREAL, {&rbin},
103 "Binwidth in the matrix and for [TT]-or[tt] (nm)" },
104 { "-mmax", FALSE, etREAL, {&mmax},
105 "Maximum density in the matrix, 0 is calculate (1/nm)" },
106 { "-nlevels", FALSE, etINT, {&nlev},
107 "Number of levels in the matrix" },
108 { "-nr", FALSE, etINT, {&nr},
109 "Number of curves for the [TT]-or[tt] output" },
110 { "-fr", FALSE, etINT, {&fshift},
111 "Frame spacing for the [TT]-or[tt] output" },
112 { "-rt", FALSE, etREAL, {&rint},
113 "Integration limit for the [TT]-ot[tt] output (nm)" },
114 { "-ft", FALSE, etINT, {&ftmax},
115 "Number of frames in the [TT]-ot[tt] output, 0 is plot all" }
117 #define NPA asize(pa)
120 { efTRX, NULL, NULL, ffREAD },
121 { efTPS, NULL, NULL, ffREAD },
122 { efNDX, NULL, NULL, ffOPTRD },
123 { efXPM, "-om", "vanhove", ffOPTWR },
124 { efXVG, "-or", "vanhove_r", ffOPTWR },
125 { efXVG, "-ot", "vanhove_t", ffOPTWR }
127 #define NFILE asize(fnm)
130 const char *matfile, *otfile, *orfile;
134 matrix boxtop, box, *sbox, avbox, corr;
135 rvec *xtop, *x, **sx;
136 int isize, nalloc, nallocn, natom;
140 int nfr, f, ff, i, m, mat_nx = 0, nbin = 0, bin, mbin, fbin;
141 real *time, t, invbin = 0, rmax2 = 0, rint2 = 0, d2;
142 real invsbin = 0, matmax, normfac, dt, *tickx, *ticky;
143 char buf[STRLEN], **legend;
145 int *pt = NULL, **pr = NULL, *mcount = NULL, *tcount = NULL, *rcount = NULL;
147 t_rgb rlo = {1, 1, 1}, rhi = {0, 0, 0};
149 parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
150 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
152 matfile = opt2fn_null("-om", NFILE, fnm);
153 if (opt2parg_bSet("-fr", NPA, pa))
155 orfile = opt2fn("-or", NFILE, fnm);
159 orfile = opt2fn_null("-or", NFILE, fnm);
161 if (opt2parg_bSet("-rt", NPA, pa))
163 otfile = opt2fn("-ot", NFILE, fnm);
167 otfile = opt2fn_null("-ot", NFILE, fnm);
170 if (!matfile && !otfile && !orfile)
173 "For output set one (or more) of the output file options\n");
177 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, boxtop,
179 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
187 natom = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
194 srenew(time, nalloc);
195 srenew(sbox, nalloc);
200 copy_mat(box, sbox[nfr]);
201 /* This assumes that the off-diagonal box elements
202 * are not affected by jumps across the periodic boundaries.
204 m_add(avbox, box, avbox);
205 snew(sx[nfr], isize);
206 for (i = 0; i < isize; i++)
208 copy_rvec(x[index[i]], sx[nfr][i]);
213 while (read_next_x(oenv, status, &t, natom, x, box));
219 fprintf(stderr, "Read %d frames\n", nfr);
221 dt = (time[nfr-1] - time[0])/(nfr - 1);
222 /* Some ugly rounding to get nice nice times in the output */
223 dt = (int)(10000.0*dt + 0.5)/10000.0;
229 if (fmmax <= 0 || fmmax >= nfr)
234 nbin = (int)(rmax*invbin + 0.5);
242 mat_nx = sqrt(fmmax*dt)*invsbin + 1;
245 for (f = 0; f < mat_nx; f++)
249 rmax2 = sqr(nbin*rbin);
250 /* Initialize time zero */
251 mat[0][0] = nfr*isize;
275 /* Initialize time zero */
284 msmul(avbox, 1.0/nfr, avbox);
285 for (f = 0; f < nfr; f++)
289 fprintf(stderr, "\rProcessing frame %d", f);
291 /* Scale all the configuration to the average box */
292 m_inv_ur0(sbox[f], corr);
293 mmul_ur0(avbox, corr, corr);
294 for (i = 0; i < isize; i++)
296 mvmul_ur0(corr, sx[f][i], sx[f][i]);
299 /* Correct for periodic jumps */
300 for (m = DIM-1; m >= 0; m--)
302 while (sx[f][i][m] - sx[f-1][i][m] > 0.5*avbox[m][m])
304 rvec_dec(sx[f][i], avbox[m]);
306 while (sx[f][i][m] - sx[f-1][i][m] <= -0.5*avbox[m][m])
308 rvec_inc(sx[f][i], avbox[m]);
313 for (ff = 0; ff < f; ff++)
316 if (fbin <= fmmax || fbin <= ftmax)
324 mbin = (int)(sqrt(fbin*dt)*invsbin + 0.5);
326 for (i = 0; i < isize; i++)
328 d2 = distance2(sx[f][i], sx[ff][i]);
329 if (mbin < mat_nx && d2 < rmax2)
331 bin = (int)(sqrt(d2)*invbin + 0.5);
337 if (fbin <= ftmax && d2 <= rint2)
354 for (fbin = 0; fbin < nr; fbin++)
356 ff = f - (fbin + 1)*fshift;
359 for (i = 0; i < isize; i++)
361 d2 = distance2(sx[f][i], sx[ff][i]);
362 bin = (int)(sqrt(d2)*invbin);
365 nallocn = 10*(bin/10) + 11;
366 for (m = 0; m < nr; m++)
368 srenew(pr[m], nallocn);
369 for (i = nalloc; i < nallocn; i++)
383 fprintf(stderr, "\n");
388 for (f = 0; f < mat_nx; f++)
390 normfac = 1.0/(mcount[f]*isize*rbin);
391 for (i = 0; i < nbin; i++)
393 mat[f][i] *= normfac;
394 if (mat[f][i] > matmax && (f != 0 || i != 0))
400 fprintf(stdout, "Value at (0,0): %.3f, maximum of the rest %.3f\n",
407 for (f = 0; f < mat_nx; f++)
419 for (i = 0; i <= nbin; i++)
423 fp = ffopen(matfile, "w");
424 write_xpm(fp, MAT_SPATIAL_Y, "Van Hove function", "G (1/nm)",
425 sbin == 0 ? "time (ps)" : "sqrt(time) (ps^1/2)", "r (nm)",
426 mat_nx, nbin, tickx, ticky, mat, 0, matmax, rlo, rhi, &nlev);
432 fp = xvgropen(orfile, "Van Hove function", "r (nm)", "G (nm\\S-1\\N)", oenv);
433 fprintf(fp, "@ subtitle \"for particles in group %s\"\n", grpname);
435 for (fbin = 0; fbin < nr; fbin++)
437 sprintf(buf, "%g ps", (fbin + 1)*fshift*dt);
438 legend[fbin] = strdup(buf);
440 xvgr_legend(fp, nr, (const char**)legend, oenv);
441 for (i = 0; i < nalloc; i++)
443 fprintf(fp, "%g", i*rbin);
444 for (fbin = 0; fbin < nr; fbin++)
447 (real)pr[fbin][i]/(rcount[fbin]*isize*rbin*(i == 0 ? 0.5 : 1)));
456 sprintf(buf, "Probability of moving less than %g nm", rint);
457 fp = xvgropen(otfile, buf, "t (ps)", "", oenv);
458 fprintf(fp, "@ subtitle \"for particles in group %s\"\n", grpname);
459 for (f = 0; f <= ftmax; f++)
461 fprintf(fp, "%g %g\n", f*dt, (real)pt[f]/(tcount[f]*isize));
466 do_view(oenv, matfile, NULL);
467 do_view(oenv, orfile, NULL);
468 do_view(oenv, otfile, NULL);