Remove continuation from convert_tpr
[gromacs.git] / src / gromacs / gmxana / gmx_vanhove.cpp
blobdb3da07cff5eeb9089dbcedaf2d5d5407d39e1ce
1 /*
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,2015,2016, 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.
37 #include "gmxpre.h"
39 #include <cmath>
40 #include <cstdlib>
41 #include <cstring>
43 #include "gromacs/commandline/pargs.h"
44 #include "gromacs/commandline/viewit.h"
45 #include "gromacs/fileio/confio.h"
46 #include "gromacs/fileio/matio.h"
47 #include "gromacs/fileio/trxio.h"
48 #include "gromacs/fileio/xvgr.h"
49 #include "gromacs/gmxana/gmx_ana.h"
50 #include "gromacs/gmxana/gstat.h"
51 #include "gromacs/math/functions.h"
52 #include "gromacs/math/invertmatrix.h"
53 #include "gromacs/math/utilities.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/pbcutil/pbc.h"
56 #include "gromacs/topology/index.h"
57 #include "gromacs/topology/topology.h"
58 #include "gromacs/utility/arraysize.h"
59 #include "gromacs/utility/cstringutil.h"
60 #include "gromacs/utility/futil.h"
61 #include "gromacs/utility/gmxassert.h"
62 #include "gromacs/utility/smalloc.h"
65 int gmx_vanhove(int argc, char *argv[])
67 const char *desc[] = {
68 "[THISMODULE] computes the Van Hove correlation function.",
69 "The Van Hove G(r,t) is the probability that a particle that is at r[SUB]0[sub]",
70 "at time zero can be found at position r[SUB]0[sub]+r at time t.",
71 "[THISMODULE] determines G not for a vector r, but for the length of r.",
72 "Thus it gives the probability that a particle moves a distance of r",
73 "in time t.",
74 "Jumps across the periodic boundaries are removed.",
75 "Corrections are made for scaling due to isotropic",
76 "or anisotropic pressure coupling.",
77 "[PAR]",
78 "With option [TT]-om[tt] the whole matrix can be written as a function",
79 "of t and r or as a function of [SQRT]t[sqrt] and r (option [TT]-sqrt[tt]).",
80 "[PAR]",
81 "With option [TT]-or[tt] the Van Hove function is plotted for one",
82 "or more values of t. Option [TT]-nr[tt] sets the number of times,",
83 "option [TT]-fr[tt] the number spacing between the times.",
84 "The binwidth is set with option [TT]-rbin[tt]. The number of bins",
85 "is determined automatically.",
86 "[PAR]",
87 "With option [TT]-ot[tt] the integral up to a certain distance",
88 "(option [TT]-rt[tt]) is plotted as a function of time.",
89 "[PAR]",
90 "For all frames that are read the coordinates of the selected particles",
91 "are stored in memory. Therefore the program may use a lot of memory.",
92 "For options [TT]-om[tt] and [TT]-ot[tt] the program may be slow.",
93 "This is because the calculation scales as the number of frames times",
94 "[TT]-fm[tt] or [TT]-ft[tt].",
95 "Note that with the [TT]-dt[tt] option the memory usage and calculation",
96 "time can be reduced."
98 static int fmmax = 0, ftmax = 0, nlev = 81, nr = 1, fshift = 0;
99 static real sbin = 0, rmax = 2, rbin = 0.01, mmax = 0, rint = 0;
100 t_pargs pa[] = {
101 { "-sqrt", FALSE, etREAL, {&sbin},
102 "Use [SQRT]t[sqrt] on the matrix axis which binspacing # in [SQRT]ps[sqrt]" },
103 { "-fm", FALSE, etINT, {&fmmax},
104 "Number of frames in the matrix, 0 is plot all" },
105 { "-rmax", FALSE, etREAL, {&rmax},
106 "Maximum r in the matrix (nm)" },
107 { "-rbin", FALSE, etREAL, {&rbin},
108 "Binwidth in the matrix and for [TT]-or[tt] (nm)" },
109 { "-mmax", FALSE, etREAL, {&mmax},
110 "Maximum density in the matrix, 0 is calculate (1/nm)" },
111 { "-nlevels", FALSE, etINT, {&nlev},
112 "Number of levels in the matrix" },
113 { "-nr", FALSE, etINT, {&nr},
114 "Number of curves for the [TT]-or[tt] output" },
115 { "-fr", FALSE, etINT, {&fshift},
116 "Frame spacing for the [TT]-or[tt] output" },
117 { "-rt", FALSE, etREAL, {&rint},
118 "Integration limit for the [TT]-ot[tt] output (nm)" },
119 { "-ft", FALSE, etINT, {&ftmax},
120 "Number of frames in the [TT]-ot[tt] output, 0 is plot all" }
122 #define NPA asize(pa)
124 t_filenm fnm[] = {
125 { efTRX, NULL, NULL, ffREAD },
126 { efTPS, NULL, NULL, ffREAD },
127 { efNDX, NULL, NULL, ffOPTRD },
128 { efXPM, "-om", "vanhove", ffOPTWR },
129 { efXVG, "-or", "vanhove_r", ffOPTWR },
130 { efXVG, "-ot", "vanhove_t", ffOPTWR }
132 #define NFILE asize(fnm)
134 gmx_output_env_t *oenv;
135 const char *matfile, *otfile, *orfile;
136 t_topology top;
137 int ePBC;
138 matrix boxtop, box, *sbox, avbox, corr;
139 rvec *xtop, *x, **sx;
140 int isize, nalloc, nallocn;
141 t_trxstatus *status;
142 int *index;
143 char *grpname;
144 int nfr, f, ff, i, m, mat_nx = 0, nbin = 0, bin, mbin, fbin;
145 real *time, t, invbin = 0, rmax2 = 0, rint2 = 0, d2;
146 real invsbin = 0, matmax, normfac, dt, *tickx, *ticky;
147 char buf[STRLEN], **legend;
148 real **mat = NULL;
149 int *pt = NULL, **pr = NULL, *mcount = NULL, *tcount = NULL, *rcount = NULL;
150 FILE *fp;
151 t_rgb rlo = {1, 1, 1}, rhi = {0, 0, 0};
153 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME,
154 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
156 return 0;
159 matfile = opt2fn_null("-om", NFILE, fnm);
160 if (opt2parg_bSet("-fr", NPA, pa))
162 orfile = opt2fn("-or", NFILE, fnm);
164 else
166 orfile = opt2fn_null("-or", NFILE, fnm);
168 if (opt2parg_bSet("-rt", NPA, pa))
170 otfile = opt2fn("-ot", NFILE, fnm);
172 else
174 otfile = opt2fn_null("-ot", NFILE, fnm);
177 if (!matfile && !otfile && !orfile)
179 fprintf(stderr,
180 "For output set one (or more) of the output file options\n");
181 exit(0);
184 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &ePBC, &xtop, NULL, boxtop,
185 FALSE);
186 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
188 nalloc = 0;
189 time = NULL;
190 sbox = NULL;
191 sx = NULL;
192 clear_mat(avbox);
194 read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
195 nfr = 0;
198 if (nfr >= nalloc)
200 nalloc += 100;
201 srenew(time, nalloc);
202 srenew(sbox, nalloc);
203 srenew(sx, nalloc);
205 GMX_RELEASE_ASSERT(time != NULL, "Memory allocation failure; time array is NULL");
206 GMX_RELEASE_ASSERT(sbox != NULL, "Memory allocation failure; sbox array is NULL");
208 time[nfr] = t;
209 copy_mat(box, sbox[nfr]);
210 /* This assumes that the off-diagonal box elements
211 * are not affected by jumps across the periodic boundaries.
213 m_add(avbox, box, avbox);
214 snew(sx[nfr], isize);
215 for (i = 0; i < isize; i++)
217 copy_rvec(x[index[i]], sx[nfr][i]);
220 nfr++;
222 while (read_next_x(oenv, status, &t, x, box));
224 /* clean up */
225 sfree(x);
226 close_trj(status);
228 fprintf(stderr, "Read %d frames\n", nfr);
230 dt = (time[nfr-1] - time[0])/(nfr - 1);
231 /* Some ugly rounding to get nice nice times in the output */
232 dt = static_cast<int>((10000.0*dt + 0.5)/10000.0);
234 invbin = 1.0/rbin;
236 if (matfile)
238 if (fmmax <= 0 || fmmax >= nfr)
240 fmmax = nfr - 1;
242 snew(mcount, fmmax);
243 nbin = static_cast<int>(rmax*invbin + 0.5);
244 if (sbin == 0)
246 mat_nx = fmmax + 1;
248 else
250 invsbin = 1.0/sbin;
251 mat_nx = static_cast<int>(std::sqrt(fmmax*dt)*invsbin + 1);
253 snew(mat, mat_nx);
254 for (f = 0; f < mat_nx; f++)
256 snew(mat[f], nbin);
258 rmax2 = gmx::square(nbin*rbin);
259 /* Initialize time zero */
260 mat[0][0] = nfr*isize;
261 mcount[0] += nfr;
263 else
265 fmmax = 0;
268 if (orfile)
270 snew(pr, nr);
271 nalloc = 0;
272 snew(rcount, nr);
275 if (otfile)
277 if (ftmax <= 0)
279 ftmax = nfr - 1;
281 snew(tcount, ftmax);
282 snew(pt, nfr);
283 rint2 = rint*rint;
284 /* Initialize time zero */
285 pt[0] = nfr*isize;
286 tcount[0] += nfr;
288 else
290 ftmax = 0;
293 msmul(avbox, 1.0/nfr, avbox);
294 for (f = 0; f < nfr; f++)
296 if (f % 100 == 0)
298 fprintf(stderr, "\rProcessing frame %d", f);
299 fflush(stderr);
301 if (ePBC != epbcNONE)
303 /* Scale all the configuration to the average box */
304 gmx::invertBoxMatrix(sbox[f], corr);
305 mmul_ur0(avbox, corr, corr);
306 for (i = 0; i < isize; i++)
308 mvmul_ur0(corr, sx[f][i], sx[f][i]);
309 if (f > 0)
311 /* Correct for periodic jumps */
312 for (m = DIM-1; m >= 0; m--)
314 while (sx[f][i][m] - sx[f-1][i][m] > 0.5*avbox[m][m])
316 rvec_dec(sx[f][i], avbox[m]);
318 while (sx[f][i][m] - sx[f-1][i][m] <= -0.5*avbox[m][m])
320 rvec_inc(sx[f][i], avbox[m]);
326 for (ff = 0; ff < f; ff++)
328 fbin = f - ff;
329 if (fbin <= fmmax || fbin <= ftmax)
331 if (sbin == 0)
333 mbin = fbin;
335 else
337 mbin = static_cast<int>(std::sqrt(fbin*dt)*invsbin + 0.5);
339 for (i = 0; i < isize; i++)
341 d2 = distance2(sx[f][i], sx[ff][i]);
342 if (mbin < mat_nx && d2 < rmax2)
344 bin = static_cast<int>(std::sqrt(d2)*invbin + 0.5);
345 if (bin < nbin)
347 mat[mbin][bin] += 1;
350 if (fbin <= ftmax && d2 <= rint2)
352 pt[fbin]++;
355 if (matfile)
357 mcount[mbin]++;
359 if (otfile)
361 tcount[fbin]++;
365 if (orfile)
367 for (fbin = 0; fbin < nr; fbin++)
369 ff = f - (fbin + 1)*fshift;
370 if (ff >= 0)
372 for (i = 0; i < isize; i++)
374 d2 = distance2(sx[f][i], sx[ff][i]);
375 bin = static_cast<int>(std::sqrt(d2)*invbin + 0.5);
376 if (bin >= nalloc)
378 nallocn = 10*(bin/10) + 11;
379 for (m = 0; m < nr; m++)
381 srenew(pr[m], nallocn);
382 for (i = nalloc; i < nallocn; i++)
384 pr[m][i] = 0;
387 nalloc = nallocn;
389 pr[fbin][bin]++;
391 rcount[fbin]++;
396 fprintf(stderr, "\n");
398 if (matfile)
400 matmax = 0;
401 for (f = 0; f < mat_nx; f++)
403 normfac = 1.0/(mcount[f]*isize*rbin);
404 for (i = 0; i < nbin; i++)
406 mat[f][i] *= normfac;
407 if (mat[f][i] > matmax && (f != 0 || i != 0))
409 matmax = mat[f][i];
413 fprintf(stdout, "Value at (0,0): %.3f, maximum of the rest %.3f\n",
414 mat[0][0], matmax);
415 if (mmax > 0)
417 matmax = mmax;
419 snew(tickx, mat_nx);
420 for (f = 0; f < mat_nx; f++)
422 if (sbin == 0)
424 tickx[f] = f*dt;
426 else
428 tickx[f] = f*sbin;
431 snew(ticky, nbin+1);
432 for (i = 0; i <= nbin; i++)
434 ticky[i] = i*rbin;
436 fp = gmx_ffopen(matfile, "w");
437 write_xpm(fp, MAT_SPATIAL_Y, "Van Hove function", "G (1/nm)",
438 sbin == 0 ? "time (ps)" : "sqrt(time) (ps^1/2)", "r (nm)",
439 mat_nx, nbin, tickx, ticky, mat, 0, matmax, rlo, rhi, &nlev);
440 gmx_ffclose(fp);
443 if (orfile)
445 fp = xvgropen(orfile, "Van Hove function", "r (nm)", "G (nm\\S-1\\N)", oenv);
446 if (output_env_get_print_xvgr_codes(oenv))
448 fprintf(fp, "@ subtitle \"for particles in group %s\"\n", grpname);
450 snew(legend, nr);
451 for (fbin = 0; fbin < nr; fbin++)
453 sprintf(buf, "%g ps", (fbin + 1)*fshift*dt);
454 legend[fbin] = gmx_strdup(buf);
456 xvgr_legend(fp, nr, (const char**)legend, oenv);
457 for (i = 0; i < nalloc; i++)
459 fprintf(fp, "%g", i*rbin);
460 for (fbin = 0; fbin < nr; fbin++)
462 fprintf(fp, " %g", static_cast<real>(pr[fbin][i])/(rcount[fbin]*isize*rbin*(i == 0 ? 0.5 : 1.0)));
464 fprintf(fp, "\n");
466 xvgrclose(fp);
469 if (otfile)
471 sprintf(buf, "Probability of moving less than %g nm", rint);
472 fp = xvgropen(otfile, buf, "t (ps)", "", oenv);
473 if (output_env_get_print_xvgr_codes(oenv))
475 fprintf(fp, "@ subtitle \"for particles in group %s\"\n", grpname);
477 for (f = 0; f <= ftmax; f++)
479 fprintf(fp, "%g %g\n", f*dt, static_cast<real>(pt[f])/(tcount[f]*isize));
481 xvgrclose(fp);
484 do_view(oenv, matfile, NULL);
485 do_view(oenv, orfile, NULL);
486 do_view(oenv, otfile, NULL);
488 return 0;