7dae7ea80b1f45d5041725ebbe5c417cf9b462b9
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_vanhove.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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,2017,2018,2019, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include <cmath>
40 #include <cstdlib>
41 #include <cstring>
42
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"
63
64
65 int gmx_vanhove(int argc, char* argv[])
66 {
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."
97     };
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",
102           FALSE,
103           etREAL,
104           { &sbin },
105           "Use [SQRT]t[sqrt] on the matrix axis which binspacing # in [SQRT]ps[sqrt]" },
106         { "-fm", FALSE, etINT, { &fmmax }, "Number of frames in the matrix, 0 is plot all" },
107         { "-rmax", FALSE, etREAL, { &rmax }, "Maximum r in the matrix (nm)" },
108         { "-rbin", FALSE, etREAL, { &rbin }, "Binwidth in the matrix and for [TT]-or[tt] (nm)" },
109         { "-mmax",
110           FALSE,
111           etREAL,
112           { &mmax },
113           "Maximum density in the matrix, 0 is calculate (1/nm)" },
114         { "-nlevels", FALSE, etINT, { &nlev }, "Number of levels in the matrix" },
115         { "-nr", FALSE, etINT, { &nr }, "Number of curves for the [TT]-or[tt] output" },
116         { "-fr", FALSE, etINT, { &fshift }, "Frame spacing for the [TT]-or[tt] output" },
117         { "-rt", FALSE, etREAL, { &rint }, "Integration limit for the [TT]-ot[tt] output (nm)" },
118         { "-ft",
119           FALSE,
120           etINT,
121           { &ftmax },
122           "Number of frames in the [TT]-ot[tt] output, 0 is plot all" }
123     };
124 #define NPA asize(pa)
125
126     t_filenm fnm[] = {
127         { efTRX, nullptr, nullptr, ffREAD },    { efTPS, nullptr, nullptr, ffREAD },
128         { efNDX, nullptr, nullptr, ffOPTRD },   { efXPM, "-om", "vanhove", ffOPTWR },
129         { efXVG, "-or", "vanhove_r", ffOPTWR }, { efXVG, "-ot", "vanhove_t", ffOPTWR }
130     };
131 #define NFILE asize(fnm)
132
133     gmx_output_env_t* oenv;
134     const char *      matfile, *otfile, *orfile;
135     t_topology        top;
136     int               ePBC;
137     matrix            boxtop, box, *sbox, avbox, corr;
138     rvec *            xtop, *x, **sx;
139     int               isize, nalloc, nallocn;
140     t_trxstatus*      status;
141     int*              index;
142     char*             grpname;
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;
147     real**            mat = nullptr;
148     int * pt = nullptr, **pr = nullptr, *mcount = nullptr, *tcount = nullptr, *rcount = nullptr;
149     FILE* fp;
150     t_rgb rlo = { 1, 1, 1 }, rhi = { 0, 0, 0 };
151
152     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME, NFILE, fnm, asize(pa), pa,
153                            asize(desc), desc, 0, nullptr, &oenv))
154     {
155         return 0;
156     }
157
158     matfile = opt2fn_null("-om", NFILE, fnm);
159     if (opt2parg_bSet("-fr", NPA, pa))
160     {
161         orfile = opt2fn("-or", NFILE, fnm);
162     }
163     else
164     {
165         orfile = opt2fn_null("-or", NFILE, fnm);
166     }
167     if (opt2parg_bSet("-rt", NPA, pa))
168     {
169         otfile = opt2fn("-ot", NFILE, fnm);
170     }
171     else
172     {
173         otfile = opt2fn_null("-ot", NFILE, fnm);
174     }
175
176     if (!matfile && !otfile && !orfile)
177     {
178         fprintf(stderr, "For output set one (or more) of the output file options\n");
179         exit(0);
180     }
181
182     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &ePBC, &xtop, nullptr, boxtop, FALSE);
183     get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
184
185     nalloc = 0;
186     time   = nullptr;
187     sbox   = nullptr;
188     sx     = nullptr;
189     clear_mat(avbox);
190
191     read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
192     nfr = 0;
193     do
194     {
195         if (nfr >= nalloc)
196         {
197             nalloc += 100;
198             srenew(time, nalloc);
199             srenew(sbox, nalloc);
200             srenew(sx, nalloc);
201         }
202         GMX_RELEASE_ASSERT(time != nullptr, "Memory allocation failure; time array is NULL");
203         GMX_RELEASE_ASSERT(sbox != nullptr, "Memory allocation failure; sbox array is NULL");
204
205         time[nfr] = t;
206         copy_mat(box, sbox[nfr]);
207         /* This assumes that the off-diagonal box elements
208          * are not affected by jumps across the periodic boundaries.
209          */
210         m_add(avbox, box, avbox);
211         snew(sx[nfr], isize);
212         for (i = 0; i < isize; i++)
213         {
214             copy_rvec(x[index[i]], sx[nfr][i]);
215         }
216
217         nfr++;
218     } while (read_next_x(oenv, status, &t, x, box));
219
220     /* clean up */
221     sfree(x);
222     close_trx(status);
223
224     fprintf(stderr, "Read %d frames\n", nfr);
225
226     dt = (time[nfr - 1] - time[0]) / (nfr - 1);
227     /* Some ugly rounding to get nice nice times in the output */
228     dt = std::round(10000.0 * dt) / 10000.0;
229
230     invbin = 1.0 / rbin;
231
232     if (matfile)
233     {
234         if (fmmax <= 0 || fmmax >= nfr)
235         {
236             fmmax = nfr - 1;
237         }
238         snew(mcount, fmmax);
239         nbin = gmx::roundToInt(rmax * invbin);
240         if (sbin == 0)
241         {
242             mat_nx = fmmax + 1;
243         }
244         else
245         {
246             invsbin = 1.0 / sbin;
247             mat_nx  = static_cast<int>(std::sqrt(fmmax * dt) * invsbin + 1);
248         }
249         snew(mat, mat_nx);
250         for (f = 0; f < mat_nx; f++)
251         {
252             snew(mat[f], nbin);
253         }
254         rmax2 = gmx::square(nbin * rbin);
255         /* Initialize time zero */
256         mat[0][0] = nfr * isize;
257         mcount[0] += nfr;
258     }
259     else
260     {
261         fmmax = 0;
262     }
263
264     if (orfile)
265     {
266         snew(pr, nr);
267         nalloc = 0;
268         snew(rcount, nr);
269     }
270
271     if (otfile)
272     {
273         if (ftmax <= 0)
274         {
275             ftmax = nfr - 1;
276         }
277         snew(tcount, ftmax);
278         snew(pt, nfr);
279         rint2 = rint * rint;
280         /* Initialize time zero */
281         pt[0] = nfr * isize;
282         tcount[0] += nfr;
283     }
284     else
285     {
286         ftmax = 0;
287     }
288
289     msmul(avbox, 1.0 / nfr, avbox);
290     for (f = 0; f < nfr; f++)
291     {
292         if (f % 100 == 0)
293         {
294             fprintf(stderr, "\rProcessing frame %d", f);
295             fflush(stderr);
296         }
297         if (ePBC != epbcNONE)
298         {
299             /* Scale all the configuration to the average box */
300             gmx::invertBoxMatrix(sbox[f], corr);
301             mmul_ur0(avbox, corr, corr);
302             for (i = 0; i < isize; i++)
303             {
304                 mvmul_ur0(corr, sx[f][i], sx[f][i]);
305                 if (f > 0)
306                 {
307                     /* Correct for periodic jumps */
308                     for (m = DIM - 1; m >= 0; m--)
309                     {
310                         while (sx[f][i][m] - sx[f - 1][i][m] > 0.5 * avbox[m][m])
311                         {
312                             rvec_dec(sx[f][i], avbox[m]);
313                         }
314                         while (sx[f][i][m] - sx[f - 1][i][m] <= -0.5 * avbox[m][m])
315                         {
316                             rvec_inc(sx[f][i], avbox[m]);
317                         }
318                     }
319                 }
320             }
321         }
322         for (ff = 0; ff < f; ff++)
323         {
324             fbin = f - ff;
325             if (fbin <= fmmax || fbin <= ftmax)
326             {
327                 if (sbin == 0)
328                 {
329                     mbin = fbin;
330                 }
331                 else
332                 {
333                     mbin = gmx::roundToInt(std::sqrt(fbin * dt) * invsbin);
334                 }
335                 for (i = 0; i < isize; i++)
336                 {
337                     d2 = distance2(sx[f][i], sx[ff][i]);
338                     if (mbin < mat_nx && d2 < rmax2)
339                     {
340                         bin = gmx::roundToInt(std::sqrt(d2) * invbin);
341                         if (bin < nbin)
342                         {
343                             mat[mbin][bin] += 1;
344                         }
345                     }
346                     if (fbin <= ftmax && d2 <= rint2)
347                     {
348                         pt[fbin]++;
349                     }
350                 }
351                 if (matfile)
352                 {
353                     mcount[mbin]++;
354                 }
355                 if (otfile)
356                 {
357                     tcount[fbin]++;
358                 }
359             }
360         }
361         if (orfile)
362         {
363             for (fbin = 0; fbin < nr; fbin++)
364             {
365                 ff = f - (fbin + 1) * fshift;
366                 if (ff >= 0)
367                 {
368                     for (i = 0; i < isize; i++)
369                     {
370                         d2  = distance2(sx[f][i], sx[ff][i]);
371                         bin = gmx::roundToInt(std::sqrt(d2) * invbin);
372                         if (bin >= nalloc)
373                         {
374                             nallocn = 10 * (bin / 10) + 11;
375                             for (m = 0; m < nr; m++)
376                             {
377                                 srenew(pr[m], nallocn);
378                                 for (i = nalloc; i < nallocn; i++)
379                                 {
380                                     pr[m][i] = 0;
381                                 }
382                             }
383                             nalloc = nallocn;
384                         }
385                         pr[fbin][bin]++;
386                     }
387                     rcount[fbin]++;
388                 }
389             }
390         }
391     }
392     fprintf(stderr, "\n");
393
394     if (matfile)
395     {
396         matmax = 0;
397         for (f = 0; f < mat_nx; f++)
398         {
399             normfac = 1.0 / (mcount[f] * isize * rbin);
400             for (i = 0; i < nbin; i++)
401             {
402                 mat[f][i] *= normfac;
403                 if (mat[f][i] > matmax && (f != 0 || i != 0))
404                 {
405                     matmax = mat[f][i];
406                 }
407             }
408         }
409         fprintf(stdout, "Value at (0,0): %.3f, maximum of the rest %.3f\n", mat[0][0], matmax);
410         if (mmax > 0)
411         {
412             matmax = mmax;
413         }
414         snew(tickx, mat_nx);
415         for (f = 0; f < mat_nx; f++)
416         {
417             if (sbin == 0)
418             {
419                 tickx[f] = f * dt;
420             }
421             else
422             {
423                 tickx[f] = f * sbin;
424             }
425         }
426         snew(ticky, nbin + 1);
427         for (i = 0; i <= nbin; i++)
428         {
429             ticky[i] = i * rbin;
430         }
431         fp = gmx_ffopen(matfile, "w");
432         write_xpm(fp, MAT_SPATIAL_Y, "Van Hove function", "G (1/nm)",
433                   sbin == 0 ? "time (ps)" : "sqrt(time) (ps^1/2)", "r (nm)", mat_nx, nbin, tickx,
434                   ticky, mat, 0, matmax, rlo, rhi, &nlev);
435         gmx_ffclose(fp);
436     }
437
438     if (orfile)
439     {
440         fp = xvgropen(orfile, "Van Hove function", "r (nm)", "G (nm\\S-1\\N)", oenv);
441         if (output_env_get_print_xvgr_codes(oenv))
442         {
443             fprintf(fp, "@ subtitle \"for particles in group %s\"\n", grpname);
444         }
445         snew(legend, nr);
446         for (fbin = 0; fbin < nr; fbin++)
447         {
448             sprintf(buf, "%g ps", (fbin + 1) * fshift * dt);
449             legend[fbin] = gmx_strdup(buf);
450         }
451         xvgr_legend(fp, nr, legend, oenv);
452         for (i = 0; i < nalloc; i++)
453         {
454             fprintf(fp, "%g", i * rbin);
455             for (fbin = 0; fbin < nr; fbin++)
456             {
457                 fprintf(fp, " %g",
458                         static_cast<real>(pr[fbin][i]
459                                           / (rcount[fbin] * isize * rbin * (i == 0 ? 0.5 : 1.0))));
460             }
461             fprintf(fp, "\n");
462         }
463         xvgrclose(fp);
464     }
465
466     if (otfile)
467     {
468         sprintf(buf, "Probability of moving less than %g nm", rint);
469         fp = xvgropen(otfile, buf, "t (ps)", "", oenv);
470         if (output_env_get_print_xvgr_codes(oenv))
471         {
472             fprintf(fp, "@ subtitle \"for particles in group %s\"\n", grpname);
473         }
474         for (f = 0; f <= ftmax; f++)
475         {
476             fprintf(fp, "%g %g\n", f * dt, static_cast<real>(pt[f]) / (tcount[f] * isize));
477         }
478         xvgrclose(fp);
479     }
480
481     do_view(oenv, matfile, nullptr);
482     do_view(oenv, orfile, nullptr);
483     do_view(oenv, otfile, nullptr);
484
485     return 0;
486 }