Merge branch 'release-4-6'
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_vanhove.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
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.
14
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.
19  *
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.
26  *
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.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Green Red Orange Magenta Azure Cyan Skyblue
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <string.h>
40 #include <ctype.h>
41 #include <math.h>
42
43 #include "sysstuff.h"
44 #include "smalloc.h"
45 #include "macros.h"
46 #include "statutil.h"
47 #include "maths.h"
48 #include "futil.h"
49 #include "index.h"
50 #include "typedefs.h"
51 #include "xvgr.h"
52 #include "gstat.h"
53 #include "tpxio.h"
54 #include "vec.h"
55 #include "matio.h"
56 #include "gmx_ana.h"
57
58
59 int gmx_vanhove(int argc, char *argv[])
60 {
61     const char *desc[] = {
62         "[TT]g_vanhove[tt] computes the Van Hove correlation function.",
63         "The Van Hove G(r,t) is the probability that a particle that is at r[SUB]0[sub]",
64         "at time zero can be found at position r[SUB]0[sub]+r at time t.",
65         "[TT]g_vanhove[tt] determines G not for a vector r, but for the length of r.",
66         "Thus it gives the probability that a particle moves a distance of r",
67         "in time t.",
68         "Jumps across the periodic boundaries are removed.",
69         "Corrections are made for scaling due to isotropic",
70         "or anisotropic pressure coupling.",
71         "[PAR]",
72         "With option [TT]-om[tt] the whole matrix can be written as a function",
73         "of t and r or as a function of [SQRT]t[sqrt] and r (option [TT]-sqrt[tt]).",
74         "[PAR]",
75         "With option [TT]-or[tt] the Van Hove function is plotted for one",
76         "or more values of t. Option [TT]-nr[tt] sets the number of times,",
77         "option [TT]-fr[tt] the number spacing between the times.",
78         "The binwidth is set with option [TT]-rbin[tt]. The number of bins",
79         "is determined automatically.",
80         "[PAR]",
81         "With option [TT]-ot[tt] the integral up to a certain distance",
82         "(option [TT]-rt[tt]) is plotted as a function of time.",
83         "[PAR]",
84         "For all frames that are read the coordinates of the selected particles",
85         "are stored in memory. Therefore the program may use a lot of memory.",
86         "For options [TT]-om[tt] and [TT]-ot[tt] the program may be slow.",
87         "This is because the calculation scales as the number of frames times",
88         "[TT]-fm[tt] or [TT]-ft[tt].",
89         "Note that with the [TT]-dt[tt] option the memory usage and calculation",
90         "time can be reduced."
91     };
92     static int  fmmax = 0, ftmax = 0, nlev = 81, nr = 1, fshift = 0;
93     static real sbin  = 0, rmax = 2, rbin = 0.01, mmax = 0, rint = 0;
94     t_pargs     pa[]  = {
95         { "-sqrt",    FALSE, etREAL, {&sbin},
96           "Use [SQRT]t[sqrt] on the matrix axis which binspacing # in [SQRT]ps[sqrt]" },
97         { "-fm",      FALSE, etINT, {&fmmax},
98           "Number of frames in the matrix, 0 is plot all" },
99         { "-rmax",    FALSE, etREAL, {&rmax},
100           "Maximum r in the matrix (nm)" },
101         { "-rbin",    FALSE, etREAL, {&rbin},
102           "Binwidth in the matrix and for [TT]-or[tt] (nm)" },
103         { "-mmax",    FALSE, etREAL, {&mmax},
104           "Maximum density in the matrix, 0 is calculate (1/nm)" },
105         { "-nlevels", FALSE, etINT,  {&nlev},
106           "Number of levels in the matrix" },
107         { "-nr",      FALSE, etINT, {&nr},
108           "Number of curves for the [TT]-or[tt] output" },
109         { "-fr",      FALSE, etINT, {&fshift},
110           "Frame spacing for the [TT]-or[tt] output" },
111         { "-rt",      FALSE, etREAL, {&rint},
112           "Integration limit for the [TT]-ot[tt] output (nm)" },
113         { "-ft",      FALSE, etINT, {&ftmax},
114           "Number of frames in the [TT]-ot[tt] output, 0 is plot all" }
115     };
116 #define NPA asize(pa)
117
118     t_filenm fnm[] = {
119         { efTRX, NULL, NULL,  ffREAD },
120         { efTPS, NULL, NULL,  ffREAD },
121         { efNDX, NULL, NULL,  ffOPTRD },
122         { efXPM, "-om", "vanhove", ffOPTWR },
123         { efXVG, "-or", "vanhove_r", ffOPTWR },
124         { efXVG, "-ot", "vanhove_t", ffOPTWR }
125     };
126 #define NFILE asize(fnm)
127
128     output_env_t oenv;
129     const char  *matfile, *otfile, *orfile;
130     char         title[256];
131     t_topology   top;
132     int          ePBC;
133     matrix       boxtop, box, *sbox, avbox, corr;
134     rvec        *xtop, *x, **sx;
135     int          isize, nalloc, nallocn, natom;
136     t_trxstatus *status;
137     atom_id     *index;
138     char        *grpname;
139     int          nfr, f, ff, i, m, mat_nx = 0, nbin = 0, bin, mbin, fbin;
140     real        *time, t, invbin = 0, rmax2 = 0, rint2 = 0, d2;
141     real         invsbin = 0, matmax, normfac, dt, *tickx, *ticky;
142     char         buf[STRLEN], **legend;
143     real       **mat = NULL;
144     int         *pt  = NULL, **pr = NULL, *mcount = NULL, *tcount = NULL, *rcount = NULL;
145     FILE        *fp;
146     t_rgb        rlo = {1, 1, 1}, rhi = {0, 0, 0};
147
148     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
149                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
150     {
151         return 0;
152     }
153
154     matfile = opt2fn_null("-om", NFILE, fnm);
155     if (opt2parg_bSet("-fr", NPA, pa))
156     {
157         orfile  = opt2fn("-or", NFILE, fnm);
158     }
159     else
160     {
161         orfile  = opt2fn_null("-or", NFILE, fnm);
162     }
163     if (opt2parg_bSet("-rt", NPA, pa))
164     {
165         otfile  = opt2fn("-ot", NFILE, fnm);
166     }
167     else
168     {
169         otfile  = opt2fn_null("-ot", NFILE, fnm);
170     }
171
172     if (!matfile && !otfile && !orfile)
173     {
174         fprintf(stderr,
175                 "For output set one (or more) of the output file options\n");
176         exit(0);
177     }
178
179     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, boxtop,
180                   FALSE);
181     get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
182
183     nalloc = 0;
184     time   = NULL;
185     sbox   = NULL;
186     sx     = NULL;
187     clear_mat(avbox);
188
189     natom = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
190     nfr   = 0;
191     do
192     {
193         if (nfr >= nalloc)
194         {
195             nalloc += 100;
196             srenew(time, nalloc);
197             srenew(sbox, nalloc);
198             srenew(sx, nalloc);
199         }
200
201         time[nfr] = t;
202         copy_mat(box, sbox[nfr]);
203         /* This assumes that the off-diagonal box elements
204          * are not affected by jumps across the periodic boundaries.
205          */
206         m_add(avbox, box, avbox);
207         snew(sx[nfr], isize);
208         for (i = 0; i < isize; i++)
209         {
210             copy_rvec(x[index[i]], sx[nfr][i]);
211         }
212
213         nfr++;
214     }
215     while (read_next_x(oenv, status, &t, x, box));
216
217     /* clean up */
218     sfree(x);
219     close_trj(status);
220
221     fprintf(stderr, "Read %d frames\n", nfr);
222
223     dt = (time[nfr-1] - time[0])/(nfr - 1);
224     /* Some ugly rounding to get nice nice times in the output */
225     dt = (int)(10000.0*dt + 0.5)/10000.0;
226
227     invbin = 1.0/rbin;
228
229     if (matfile)
230     {
231         if (fmmax <= 0 || fmmax >= nfr)
232         {
233             fmmax = nfr - 1;
234         }
235         snew(mcount, fmmax);
236         nbin = (int)(rmax*invbin + 0.5);
237         if (sbin == 0)
238         {
239             mat_nx = fmmax + 1;
240         }
241         else
242         {
243             invsbin = 1.0/sbin;
244             mat_nx  = sqrt(fmmax*dt)*invsbin + 1;
245         }
246         snew(mat, mat_nx);
247         for (f = 0; f < mat_nx; f++)
248         {
249             snew(mat[f], nbin);
250         }
251         rmax2 = sqr(nbin*rbin);
252         /* Initialize time zero */
253         mat[0][0]  = nfr*isize;
254         mcount[0] += nfr;
255     }
256     else
257     {
258         fmmax = 0;
259     }
260
261     if (orfile)
262     {
263         snew(pr, nr);
264         nalloc = 0;
265         snew(rcount, nr);
266     }
267
268     if (otfile)
269     {
270         if (ftmax <= 0)
271         {
272             ftmax = nfr - 1;
273         }
274         snew(tcount, ftmax);
275         snew(pt, nfr);
276         rint2 = rint*rint;
277         /* Initialize time zero */
278         pt[0]      = nfr*isize;
279         tcount[0] += nfr;
280     }
281     else
282     {
283         ftmax = 0;
284     }
285
286     msmul(avbox, 1.0/nfr, avbox);
287     for (f = 0; f < nfr; f++)
288     {
289         if (f % 100 == 0)
290         {
291             fprintf(stderr, "\rProcessing frame %d", f);
292         }
293         /* Scale all the configuration to the average box */
294         m_inv_ur0(sbox[f], corr);
295         mmul_ur0(avbox, corr, corr);
296         for (i = 0; i < isize; i++)
297         {
298             mvmul_ur0(corr, sx[f][i], sx[f][i]);
299             if (f > 0)
300             {
301                 /* Correct for periodic jumps */
302                 for (m = DIM-1; m >= 0; m--)
303                 {
304                     while (sx[f][i][m] - sx[f-1][i][m] > 0.5*avbox[m][m])
305                     {
306                         rvec_dec(sx[f][i], avbox[m]);
307                     }
308                     while (sx[f][i][m] - sx[f-1][i][m] <= -0.5*avbox[m][m])
309                     {
310                         rvec_inc(sx[f][i], avbox[m]);
311                     }
312                 }
313             }
314         }
315         for (ff = 0; ff < f; ff++)
316         {
317             fbin = f - ff;
318             if (fbin <= fmmax || fbin <= ftmax)
319             {
320                 if (sbin == 0)
321                 {
322                     mbin = fbin;
323                 }
324                 else
325                 {
326                     mbin = (int)(sqrt(fbin*dt)*invsbin + 0.5);
327                 }
328                 for (i = 0; i < isize; i++)
329                 {
330                     d2 = distance2(sx[f][i], sx[ff][i]);
331                     if (mbin < mat_nx && d2 < rmax2)
332                     {
333                         bin = (int)(sqrt(d2)*invbin + 0.5);
334                         if (bin < nbin)
335                         {
336                             mat[mbin][bin] += 1;
337                         }
338                     }
339                     if (fbin <= ftmax && d2 <= rint2)
340                     {
341                         pt[fbin]++;
342                     }
343                 }
344                 if (matfile)
345                 {
346                     mcount[mbin]++;
347                 }
348                 if (otfile)
349                 {
350                     tcount[fbin]++;
351                 }
352             }
353         }
354         if (orfile)
355         {
356             for (fbin = 0; fbin < nr; fbin++)
357             {
358                 ff = f - (fbin + 1)*fshift;
359                 if (ff >= 0)
360                 {
361                     for (i = 0; i < isize; i++)
362                     {
363                         d2  = distance2(sx[f][i], sx[ff][i]);
364                         bin = (int)(sqrt(d2)*invbin);
365                         if (bin >= nalloc)
366                         {
367                             nallocn = 10*(bin/10) + 11;
368                             for (m = 0; m < nr; m++)
369                             {
370                                 srenew(pr[m], nallocn);
371                                 for (i = nalloc; i < nallocn; i++)
372                                 {
373                                     pr[m][i] = 0;
374                                 }
375                             }
376                             nalloc = nallocn;
377                         }
378                         pr[fbin][bin]++;
379                     }
380                     rcount[fbin]++;
381                 }
382             }
383         }
384     }
385     fprintf(stderr, "\n");
386
387     if (matfile)
388     {
389         matmax = 0;
390         for (f = 0; f < mat_nx; f++)
391         {
392             normfac = 1.0/(mcount[f]*isize*rbin);
393             for (i = 0; i < nbin; i++)
394             {
395                 mat[f][i] *= normfac;
396                 if (mat[f][i] > matmax && (f != 0 || i != 0))
397                 {
398                     matmax = mat[f][i];
399                 }
400             }
401         }
402         fprintf(stdout, "Value at (0,0): %.3f, maximum of the rest %.3f\n",
403                 mat[0][0], matmax);
404         if (mmax > 0)
405         {
406             matmax = mmax;
407         }
408         snew(tickx, mat_nx);
409         for (f = 0; f < mat_nx; f++)
410         {
411             if (sbin == 0)
412             {
413                 tickx[f] = f*dt;
414             }
415             else
416             {
417                 tickx[f] = f*sbin;
418             }
419         }
420         snew(ticky, nbin+1);
421         for (i = 0; i <= nbin; i++)
422         {
423             ticky[i] = i*rbin;
424         }
425         fp = ffopen(matfile, "w");
426         write_xpm(fp, MAT_SPATIAL_Y, "Van Hove function", "G (1/nm)",
427                   sbin == 0 ? "time (ps)" : "sqrt(time) (ps^1/2)", "r (nm)",
428                   mat_nx, nbin, tickx, ticky, mat, 0, matmax, rlo, rhi, &nlev);
429         ffclose(fp);
430     }
431
432     if (orfile)
433     {
434         fp = xvgropen(orfile, "Van Hove function", "r (nm)", "G (nm\\S-1\\N)", oenv);
435         fprintf(fp, "@ subtitle \"for particles in group %s\"\n", grpname);
436         snew(legend, nr);
437         for (fbin = 0; fbin < nr; fbin++)
438         {
439             sprintf(buf, "%g ps", (fbin + 1)*fshift*dt);
440             legend[fbin] = strdup(buf);
441         }
442         xvgr_legend(fp, nr, (const char**)legend, oenv);
443         for (i = 0; i < nalloc; i++)
444         {
445             fprintf(fp, "%g", i*rbin);
446             for (fbin = 0; fbin < nr; fbin++)
447             {
448                 fprintf(fp, " %g",
449                         (real)pr[fbin][i]/(rcount[fbin]*isize*rbin*(i == 0 ? 0.5 : 1)));
450             }
451             fprintf(fp, "\n");
452         }
453         ffclose(fp);
454     }
455
456     if (otfile)
457     {
458         sprintf(buf, "Probability of moving less than %g nm", rint);
459         fp = xvgropen(otfile, buf, "t (ps)", "", oenv);
460         fprintf(fp, "@ subtitle \"for particles in group %s\"\n", grpname);
461         for (f = 0; f <= ftmax; f++)
462         {
463             fprintf(fp, "%g %g\n", f*dt, (real)pt[f]/(tcount[f]*isize));
464         }
465         ffclose(fp);
466     }
467
468     do_view(oenv, matfile, NULL);
469     do_view(oenv, orfile, NULL);
470     do_view(oenv, otfile, NULL);
471
472     return 0;
473 }