Sort all includes in src/gromacs
[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, 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/futil.h"
58 #include "gromacs/utility/smalloc.h"
59
60
61 int gmx_vanhove(int argc, char *argv[])
62 {
63     const char *desc[] = {
64         "[THISMODULE] computes the Van Hove correlation function.",
65         "The Van Hove G(r,t) is the probability that a particle that is at r[SUB]0[sub]",
66         "at time zero can be found at position r[SUB]0[sub]+r at time t.",
67         "[THISMODULE] determines G not for a vector r, but for the length of r.",
68         "Thus it gives the probability that a particle moves a distance of r",
69         "in time t.",
70         "Jumps across the periodic boundaries are removed.",
71         "Corrections are made for scaling due to isotropic",
72         "or anisotropic pressure coupling.",
73         "[PAR]",
74         "With option [TT]-om[tt] the whole matrix can be written as a function",
75         "of t and r or as a function of [SQRT]t[sqrt] and r (option [TT]-sqrt[tt]).",
76         "[PAR]",
77         "With option [TT]-or[tt] the Van Hove function is plotted for one",
78         "or more values of t. Option [TT]-nr[tt] sets the number of times,",
79         "option [TT]-fr[tt] the number spacing between the times.",
80         "The binwidth is set with option [TT]-rbin[tt]. The number of bins",
81         "is determined automatically.",
82         "[PAR]",
83         "With option [TT]-ot[tt] the integral up to a certain distance",
84         "(option [TT]-rt[tt]) is plotted as a function of time.",
85         "[PAR]",
86         "For all frames that are read the coordinates of the selected particles",
87         "are stored in memory. Therefore the program may use a lot of memory.",
88         "For options [TT]-om[tt] and [TT]-ot[tt] the program may be slow.",
89         "This is because the calculation scales as the number of frames times",
90         "[TT]-fm[tt] or [TT]-ft[tt].",
91         "Note that with the [TT]-dt[tt] option the memory usage and calculation",
92         "time can be reduced."
93     };
94     static int  fmmax = 0, ftmax = 0, nlev = 81, nr = 1, fshift = 0;
95     static real sbin  = 0, rmax = 2, rbin = 0.01, mmax = 0, rint = 0;
96     t_pargs     pa[]  = {
97         { "-sqrt",    FALSE, etREAL, {&sbin},
98           "Use [SQRT]t[sqrt] on the matrix axis which binspacing # in [SQRT]ps[sqrt]" },
99         { "-fm",      FALSE, etINT, {&fmmax},
100           "Number of frames in the matrix, 0 is plot all" },
101         { "-rmax",    FALSE, etREAL, {&rmax},
102           "Maximum r in the matrix (nm)" },
103         { "-rbin",    FALSE, etREAL, {&rbin},
104           "Binwidth in the matrix and for [TT]-or[tt] (nm)" },
105         { "-mmax",    FALSE, etREAL, {&mmax},
106           "Maximum density in the matrix, 0 is calculate (1/nm)" },
107         { "-nlevels", FALSE, etINT,  {&nlev},
108           "Number of levels in the matrix" },
109         { "-nr",      FALSE, etINT, {&nr},
110           "Number of curves for the [TT]-or[tt] output" },
111         { "-fr",      FALSE, etINT, {&fshift},
112           "Frame spacing for the [TT]-or[tt] output" },
113         { "-rt",      FALSE, etREAL, {&rint},
114           "Integration limit for the [TT]-ot[tt] output (nm)" },
115         { "-ft",      FALSE, etINT, {&ftmax},
116           "Number of frames in the [TT]-ot[tt] output, 0 is plot all" }
117     };
118 #define NPA asize(pa)
119
120     t_filenm fnm[] = {
121         { efTRX, NULL, NULL,  ffREAD },
122         { efTPS, NULL, NULL,  ffREAD },
123         { efNDX, NULL, NULL,  ffOPTRD },
124         { efXPM, "-om", "vanhove", ffOPTWR },
125         { efXVG, "-or", "vanhove_r", ffOPTWR },
126         { efXVG, "-ot", "vanhove_t", ffOPTWR }
127     };
128 #define NFILE asize(fnm)
129
130     output_env_t oenv;
131     const char  *matfile, *otfile, *orfile;
132     char         title[256];
133     t_topology   top;
134     int          ePBC;
135     matrix       boxtop, box, *sbox, avbox, corr;
136     rvec        *xtop, *x, **sx;
137     int          isize, nalloc, nallocn, natom;
138     t_trxstatus *status;
139     atom_id     *index;
140     char        *grpname;
141     int          nfr, f, ff, i, m, mat_nx = 0, nbin = 0, bin, mbin, fbin;
142     real        *time, t, invbin = 0, rmax2 = 0, rint2 = 0, d2;
143     real         invsbin = 0, matmax, normfac, dt, *tickx, *ticky;
144     char         buf[STRLEN], **legend;
145     real       **mat = NULL;
146     int         *pt  = NULL, **pr = NULL, *mcount = NULL, *tcount = NULL, *rcount = NULL;
147     FILE        *fp;
148     t_rgb        rlo = {1, 1, 1}, rhi = {0, 0, 0};
149
150     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME,
151                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
152     {
153         return 0;
154     }
155
156     matfile = opt2fn_null("-om", NFILE, fnm);
157     if (opt2parg_bSet("-fr", NPA, pa))
158     {
159         orfile  = opt2fn("-or", NFILE, fnm);
160     }
161     else
162     {
163         orfile  = opt2fn_null("-or", NFILE, fnm);
164     }
165     if (opt2parg_bSet("-rt", NPA, pa))
166     {
167         otfile  = opt2fn("-ot", NFILE, fnm);
168     }
169     else
170     {
171         otfile  = opt2fn_null("-ot", NFILE, fnm);
172     }
173
174     if (!matfile && !otfile && !orfile)
175     {
176         fprintf(stderr,
177                 "For output set one (or more) of the output file options\n");
178         exit(0);
179     }
180
181     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, boxtop,
182                   FALSE);
183     get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
184
185     nalloc = 0;
186     time   = NULL;
187     sbox   = NULL;
188     sx     = NULL;
189     clear_mat(avbox);
190
191     natom = 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         assert(time != NULL); assert(sbox != NULL);
203
204         time[nfr] = t;
205         copy_mat(box, sbox[nfr]);
206         /* This assumes that the off-diagonal box elements
207          * are not affected by jumps across the periodic boundaries.
208          */
209         m_add(avbox, box, avbox);
210         snew(sx[nfr], isize);
211         for (i = 0; i < isize; i++)
212         {
213             copy_rvec(x[index[i]], sx[nfr][i]);
214         }
215
216         nfr++;
217     }
218     while (read_next_x(oenv, status, &t, x, box));
219
220     /* clean up */
221     sfree(x);
222     close_trj(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 = (int)(10000.0*dt + 0.5)/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 = (int)(rmax*invbin + 0.5);
240         if (sbin == 0)
241         {
242             mat_nx = fmmax + 1;
243         }
244         else
245         {
246             invsbin = 1.0/sbin;
247             mat_nx  = 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 = sqr(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         }
296         /* Scale all the configuration to the average box */
297         m_inv_ur0(sbox[f], corr);
298         mmul_ur0(avbox, corr, corr);
299         for (i = 0; i < isize; i++)
300         {
301             mvmul_ur0(corr, sx[f][i], sx[f][i]);
302             if (f > 0)
303             {
304                 /* Correct for periodic jumps */
305                 for (m = DIM-1; m >= 0; m--)
306                 {
307                     while (sx[f][i][m] - sx[f-1][i][m] > 0.5*avbox[m][m])
308                     {
309                         rvec_dec(sx[f][i], avbox[m]);
310                     }
311                     while (sx[f][i][m] - sx[f-1][i][m] <= -0.5*avbox[m][m])
312                     {
313                         rvec_inc(sx[f][i], avbox[m]);
314                     }
315                 }
316             }
317         }
318         for (ff = 0; ff < f; ff++)
319         {
320             fbin = f - ff;
321             if (fbin <= fmmax || fbin <= ftmax)
322             {
323                 if (sbin == 0)
324                 {
325                     mbin = fbin;
326                 }
327                 else
328                 {
329                     mbin = (int)(sqrt(fbin*dt)*invsbin + 0.5);
330                 }
331                 for (i = 0; i < isize; i++)
332                 {
333                     d2 = distance2(sx[f][i], sx[ff][i]);
334                     if (mbin < mat_nx && d2 < rmax2)
335                     {
336                         bin = (int)(sqrt(d2)*invbin + 0.5);
337                         if (bin < nbin)
338                         {
339                             mat[mbin][bin] += 1;
340                         }
341                     }
342                     if (fbin <= ftmax && d2 <= rint2)
343                     {
344                         pt[fbin]++;
345                     }
346                 }
347                 if (matfile)
348                 {
349                     mcount[mbin]++;
350                 }
351                 if (otfile)
352                 {
353                     tcount[fbin]++;
354                 }
355             }
356         }
357         if (orfile)
358         {
359             for (fbin = 0; fbin < nr; fbin++)
360             {
361                 ff = f - (fbin + 1)*fshift;
362                 if (ff >= 0)
363                 {
364                     for (i = 0; i < isize; i++)
365                     {
366                         d2  = distance2(sx[f][i], sx[ff][i]);
367                         bin = (int)(sqrt(d2)*invbin + 0.5);
368                         if (bin >= nalloc)
369                         {
370                             nallocn = 10*(bin/10) + 11;
371                             for (m = 0; m < nr; m++)
372                             {
373                                 srenew(pr[m], nallocn);
374                                 for (i = nalloc; i < nallocn; i++)
375                                 {
376                                     pr[m][i] = 0;
377                                 }
378                             }
379                             nalloc = nallocn;
380                         }
381                         pr[fbin][bin]++;
382                     }
383                     rcount[fbin]++;
384                 }
385             }
386         }
387     }
388     fprintf(stderr, "\n");
389
390     if (matfile)
391     {
392         matmax = 0;
393         for (f = 0; f < mat_nx; f++)
394         {
395             normfac = 1.0/(mcount[f]*isize*rbin);
396             for (i = 0; i < nbin; i++)
397             {
398                 mat[f][i] *= normfac;
399                 if (mat[f][i] > matmax && (f != 0 || i != 0))
400                 {
401                     matmax = mat[f][i];
402                 }
403             }
404         }
405         fprintf(stdout, "Value at (0,0): %.3f, maximum of the rest %.3f\n",
406                 mat[0][0], matmax);
407         if (mmax > 0)
408         {
409             matmax = mmax;
410         }
411         snew(tickx, mat_nx);
412         for (f = 0; f < mat_nx; f++)
413         {
414             if (sbin == 0)
415             {
416                 tickx[f] = f*dt;
417             }
418             else
419             {
420                 tickx[f] = f*sbin;
421             }
422         }
423         snew(ticky, nbin+1);
424         for (i = 0; i <= nbin; i++)
425         {
426             ticky[i] = i*rbin;
427         }
428         fp = gmx_ffopen(matfile, "w");
429         write_xpm(fp, MAT_SPATIAL_Y, "Van Hove function", "G (1/nm)",
430                   sbin == 0 ? "time (ps)" : "sqrt(time) (ps^1/2)", "r (nm)",
431                   mat_nx, nbin, tickx, ticky, mat, 0, matmax, rlo, rhi, &nlev);
432         gmx_ffclose(fp);
433     }
434
435     if (orfile)
436     {
437         fp = xvgropen(orfile, "Van Hove function", "r (nm)", "G (nm\\S-1\\N)", oenv);
438         if (output_env_get_print_xvgr_codes(oenv))
439         {
440             fprintf(fp, "@ subtitle \"for particles in group %s\"\n", grpname);
441         }
442         snew(legend, nr);
443         for (fbin = 0; fbin < nr; fbin++)
444         {
445             sprintf(buf, "%g ps", (fbin + 1)*fshift*dt);
446             legend[fbin] = gmx_strdup(buf);
447         }
448         xvgr_legend(fp, nr, (const char**)legend, oenv);
449         for (i = 0; i < nalloc; i++)
450         {
451             fprintf(fp, "%g", i*rbin);
452             for (fbin = 0; fbin < nr; fbin++)
453             {
454                 fprintf(fp, " %g",
455                         (real)pr[fbin][i]/(rcount[fbin]*isize*rbin*(i == 0 ? 0.5 : 1)));
456             }
457             fprintf(fp, "\n");
458         }
459         gmx_ffclose(fp);
460     }
461
462     if (otfile)
463     {
464         sprintf(buf, "Probability of moving less than %g nm", rint);
465         fp = xvgropen(otfile, buf, "t (ps)", "", oenv);
466         if (output_env_get_print_xvgr_codes(oenv))
467         {
468             fprintf(fp, "@ subtitle \"for particles in group %s\"\n", grpname);
469         }
470         for (f = 0; f <= ftmax; f++)
471         {
472             fprintf(fp, "%g %g\n", f*dt, (real)pt[f]/(tcount[f]*isize));
473         }
474         gmx_ffclose(fp);
475     }
476
477     do_view(oenv, matfile, NULL);
478     do_view(oenv, orfile, NULL);
479     do_view(oenv, otfile, NULL);
480
481     return 0;
482 }