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