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