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