595ada61aeb15b7e23c21920fb6975c7aa411b0c
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_vanhove.c
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, 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 <assert.h>
40 #include <math.h>
41 #include <stdlib.h>
42 #include <string.h>
43
44 #include "gromacs/commandline/pargs.h"
45 #include "gromacs/fileio/matio.h"
46 #include "gromacs/fileio/tpxio.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/legacyheaders/macros.h"
52 #include "gromacs/legacyheaders/typedefs.h"
53 #include "gromacs/legacyheaders/viewit.h"
54 #include "gromacs/math/utilities.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/topology/index.h"
57 #include "gromacs/utility/cstringutil.h"
58 #include "gromacs/utility/futil.h"
59 #include "gromacs/utility/smalloc.h"
60
61
62 int gmx_vanhove(int argc, char *argv[])
63 {
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",
70         "in time t.",
71         "Jumps across the periodic boundaries are removed.",
72         "Corrections are made for scaling due to isotropic",
73         "or anisotropic pressure coupling.",
74         "[PAR]",
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]).",
77         "[PAR]",
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.",
83         "[PAR]",
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.",
86         "[PAR]",
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."
94     };
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;
97     t_pargs     pa[]  = {
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" }
118     };
119 #define NPA asize(pa)
120
121     t_filenm fnm[] = {
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 }
128     };
129 #define NFILE asize(fnm)
130
131     output_env_t oenv;
132     const char  *matfile, *otfile, *orfile;
133     char         title[256];
134     t_topology   top;
135     int          ePBC;
136     matrix       boxtop, box, *sbox, avbox, corr;
137     rvec        *xtop, *x, **sx;
138     int          isize, nalloc, nallocn, natom;
139     t_trxstatus *status;
140     atom_id     *index;
141     char        *grpname;
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;
146     real       **mat = NULL;
147     int         *pt  = NULL, **pr = NULL, *mcount = NULL, *tcount = NULL, *rcount = NULL;
148     FILE        *fp;
149     t_rgb        rlo = {1, 1, 1}, rhi = {0, 0, 0};
150
151     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME,
152                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
153     {
154         return 0;
155     }
156
157     matfile = opt2fn_null("-om", NFILE, fnm);
158     if (opt2parg_bSet("-fr", NPA, pa))
159     {
160         orfile  = opt2fn("-or", NFILE, fnm);
161     }
162     else
163     {
164         orfile  = opt2fn_null("-or", NFILE, fnm);
165     }
166     if (opt2parg_bSet("-rt", NPA, pa))
167     {
168         otfile  = opt2fn("-ot", NFILE, fnm);
169     }
170     else
171     {
172         otfile  = opt2fn_null("-ot", NFILE, fnm);
173     }
174
175     if (!matfile && !otfile && !orfile)
176     {
177         fprintf(stderr,
178                 "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), title, &top, &ePBC, &xtop, NULL, boxtop,
183                   FALSE);
184     get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
185
186     nalloc = 0;
187     time   = NULL;
188     sbox   = NULL;
189     sx     = NULL;
190     clear_mat(avbox);
191
192     natom = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
193     nfr   = 0;
194     do
195     {
196         if (nfr >= nalloc)
197         {
198             nalloc += 100;
199             srenew(time, nalloc);
200             srenew(sbox, nalloc);
201             srenew(sx, nalloc);
202         }
203         assert(time != NULL); assert(sbox != 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     }
219     while (read_next_x(oenv, status, &t, x, box));
220
221     /* clean up */
222     sfree(x);
223     close_trj(status);
224
225     fprintf(stderr, "Read %d frames\n", nfr);
226
227     dt = (time[nfr-1] - time[0])/(nfr - 1);
228     /* Some ugly rounding to get nice nice times in the output */
229     dt = (int)(10000.0*dt + 0.5)/10000.0;
230
231     invbin = 1.0/rbin;
232
233     if (matfile)
234     {
235         if (fmmax <= 0 || fmmax >= nfr)
236         {
237             fmmax = nfr - 1;
238         }
239         snew(mcount, fmmax);
240         nbin = (int)(rmax*invbin + 0.5);
241         if (sbin == 0)
242         {
243             mat_nx = fmmax + 1;
244         }
245         else
246         {
247             invsbin = 1.0/sbin;
248             mat_nx  = sqrt(fmmax*dt)*invsbin + 1;
249         }
250         snew(mat, mat_nx);
251         for (f = 0; f < mat_nx; f++)
252         {
253             snew(mat[f], nbin);
254         }
255         rmax2 = sqr(nbin*rbin);
256         /* Initialize time zero */
257         mat[0][0]  = nfr*isize;
258         mcount[0] += nfr;
259     }
260     else
261     {
262         fmmax = 0;
263     }
264
265     if (orfile)
266     {
267         snew(pr, nr);
268         nalloc = 0;
269         snew(rcount, nr);
270     }
271
272     if (otfile)
273     {
274         if (ftmax <= 0)
275         {
276             ftmax = nfr - 1;
277         }
278         snew(tcount, ftmax);
279         snew(pt, nfr);
280         rint2 = rint*rint;
281         /* Initialize time zero */
282         pt[0]      = nfr*isize;
283         tcount[0] += nfr;
284     }
285     else
286     {
287         ftmax = 0;
288     }
289
290     msmul(avbox, 1.0/nfr, avbox);
291     for (f = 0; f < nfr; f++)
292     {
293         if (f % 100 == 0)
294         {
295             fprintf(stderr, "\rProcessing frame %d", f);
296         }
297         /* Scale all the configuration to the average box */
298         m_inv_ur0(sbox[f], corr);
299         mmul_ur0(avbox, corr, corr);
300         for (i = 0; i < isize; i++)
301         {
302             mvmul_ur0(corr, sx[f][i], sx[f][i]);
303             if (f > 0)
304             {
305                 /* Correct for periodic jumps */
306                 for (m = DIM-1; m >= 0; m--)
307                 {
308                     while (sx[f][i][m] - sx[f-1][i][m] > 0.5*avbox[m][m])
309                     {
310                         rvec_dec(sx[f][i], avbox[m]);
311                     }
312                     while (sx[f][i][m] - sx[f-1][i][m] <= -0.5*avbox[m][m])
313                     {
314                         rvec_inc(sx[f][i], avbox[m]);
315                     }
316                 }
317             }
318         }
319         for (ff = 0; ff < f; ff++)
320         {
321             fbin = f - ff;
322             if (fbin <= fmmax || fbin <= ftmax)
323             {
324                 if (sbin == 0)
325                 {
326                     mbin = fbin;
327                 }
328                 else
329                 {
330                     mbin = (int)(sqrt(fbin*dt)*invsbin + 0.5);
331                 }
332                 for (i = 0; i < isize; i++)
333                 {
334                     d2 = distance2(sx[f][i], sx[ff][i]);
335                     if (mbin < mat_nx && d2 < rmax2)
336                     {
337                         bin = (int)(sqrt(d2)*invbin + 0.5);
338                         if (bin < nbin)
339                         {
340                             mat[mbin][bin] += 1;
341                         }
342                     }
343                     if (fbin <= ftmax && d2 <= rint2)
344                     {
345                         pt[fbin]++;
346                     }
347                 }
348                 if (matfile)
349                 {
350                     mcount[mbin]++;
351                 }
352                 if (otfile)
353                 {
354                     tcount[fbin]++;
355                 }
356             }
357         }
358         if (orfile)
359         {
360             for (fbin = 0; fbin < nr; fbin++)
361             {
362                 ff = f - (fbin + 1)*fshift;
363                 if (ff >= 0)
364                 {
365                     for (i = 0; i < isize; i++)
366                     {
367                         d2  = distance2(sx[f][i], sx[ff][i]);
368                         bin = (int)(sqrt(d2)*invbin + 0.5);
369                         if (bin >= nalloc)
370                         {
371                             nallocn = 10*(bin/10) + 11;
372                             for (m = 0; m < nr; m++)
373                             {
374                                 srenew(pr[m], nallocn);
375                                 for (i = nalloc; i < nallocn; i++)
376                                 {
377                                     pr[m][i] = 0;
378                                 }
379                             }
380                             nalloc = nallocn;
381                         }
382                         pr[fbin][bin]++;
383                     }
384                     rcount[fbin]++;
385                 }
386             }
387         }
388     }
389     fprintf(stderr, "\n");
390
391     if (matfile)
392     {
393         matmax = 0;
394         for (f = 0; f < mat_nx; f++)
395         {
396             normfac = 1.0/(mcount[f]*isize*rbin);
397             for (i = 0; i < nbin; i++)
398             {
399                 mat[f][i] *= normfac;
400                 if (mat[f][i] > matmax && (f != 0 || i != 0))
401                 {
402                     matmax = mat[f][i];
403                 }
404             }
405         }
406         fprintf(stdout, "Value at (0,0): %.3f, maximum of the rest %.3f\n",
407                 mat[0][0], matmax);
408         if (mmax > 0)
409         {
410             matmax = mmax;
411         }
412         snew(tickx, mat_nx);
413         for (f = 0; f < mat_nx; f++)
414         {
415             if (sbin == 0)
416             {
417                 tickx[f] = f*dt;
418             }
419             else
420             {
421                 tickx[f] = f*sbin;
422             }
423         }
424         snew(ticky, nbin+1);
425         for (i = 0; i <= nbin; i++)
426         {
427             ticky[i] = i*rbin;
428         }
429         fp = gmx_ffopen(matfile, "w");
430         write_xpm(fp, MAT_SPATIAL_Y, "Van Hove function", "G (1/nm)",
431                   sbin == 0 ? "time (ps)" : "sqrt(time) (ps^1/2)", "r (nm)",
432                   mat_nx, nbin, tickx, ticky, mat, 0, matmax, rlo, rhi, &nlev);
433         gmx_ffclose(fp);
434     }
435
436     if (orfile)
437     {
438         fp = xvgropen(orfile, "Van Hove function", "r (nm)", "G (nm\\S-1\\N)", oenv);
439         if (output_env_get_print_xvgr_codes(oenv))
440         {
441             fprintf(fp, "@ subtitle \"for particles in group %s\"\n", grpname);
442         }
443         snew(legend, nr);
444         for (fbin = 0; fbin < nr; fbin++)
445         {
446             sprintf(buf, "%g ps", (fbin + 1)*fshift*dt);
447             legend[fbin] = gmx_strdup(buf);
448         }
449         xvgr_legend(fp, nr, (const char**)legend, oenv);
450         for (i = 0; i < nalloc; i++)
451         {
452             fprintf(fp, "%g", i*rbin);
453             for (fbin = 0; fbin < nr; fbin++)
454             {
455                 fprintf(fp, " %g",
456                         (real)pr[fbin][i]/(rcount[fbin]*isize*rbin*(i == 0 ? 0.5 : 1)));
457             }
458             fprintf(fp, "\n");
459         }
460         xvgrclose(fp);
461     }
462
463     if (otfile)
464     {
465         sprintf(buf, "Probability of moving less than %g nm", rint);
466         fp = xvgropen(otfile, buf, "t (ps)", "", oenv);
467         if (output_env_get_print_xvgr_codes(oenv))
468         {
469             fprintf(fp, "@ subtitle \"for particles in group %s\"\n", grpname);
470         }
471         for (f = 0; f <= ftmax; f++)
472         {
473             fprintf(fp, "%g %g\n", f*dt, (real)pt[f]/(tcount[f]*isize));
474         }
475         xvgrclose(fp);
476     }
477
478     do_view(oenv, matfile, NULL);
479     do_view(oenv, orfile, NULL);
480     do_view(oenv, otfile, NULL);
481
482     return 0;
483 }