Fixed make_edi.c
[gromacs/rigid-bodies.git] / src / tools / gmx_vanhove.c
blob8e4c123edddd05cb4ed2b175afaddc574bab8464
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.
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.
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.
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.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * Green Red Orange Magenta Azure Cyan Skyblue
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
39 #include <string.h>
40 #include <ctype.h>
41 #include <math.h>
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"
60 int gmx_vanhove(int argc,char *argv[])
62 const char *desc[] = {
63 "g_vanhove computes the Van Hove correlation function.",
64 "The Van Hove G(r,t) is the probability that a particle that is at r0",
65 "at time zero can be found at position r0+r at time t.",
66 "g_vanhove 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) 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."
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) on the matrix axis which binspacing # in sqrt(ps)" },
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 -or (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 -or output" },
110 { "-fr", FALSE, etINT, {&fshift},
111 "Frame spacing for the -or output" },
112 { "-rt", FALSE, etREAL, {&rint},
113 "Integration limit for the -ot output (nm)" },
114 { "-ft", FALSE, etINT, {&ftmax},
115 "Number of frames in the -ot output, 0 is plot all" }
117 #define NPA asize(pa)
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 }
127 #define NFILE asize(fnm)
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 status,isize,nalloc,nallocn,natom;
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};
148 CopyRight(stderr,argv[0]);
150 parse_common_args(&argc,argv,PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
151 NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL,&oenv);
153 matfile = opt2fn_null("-om",NFILE,fnm);
154 if (opt2parg_bSet("-fr",NPA,pa))
155 orfile = opt2fn("-or",NFILE,fnm);
156 else
157 orfile = opt2fn_null("-or",NFILE,fnm);
158 if (opt2parg_bSet("-rt",NPA,pa))
159 otfile = opt2fn("-ot",NFILE,fnm);
160 else
161 otfile = opt2fn_null("-ot",NFILE,fnm);
163 if (!matfile && !otfile && !orfile) {
164 fprintf(stderr,
165 "For output set one (or more) of the output file options\n");
166 exit(0);
169 read_tps_conf(ftp2fn(efTPS,NFILE,fnm),title,&top,&ePBC,&xtop,NULL,boxtop,
170 FALSE);
171 get_index(&top.atoms,ftp2fn_null(efNDX,NFILE,fnm),1,&isize,&index,&grpname);
173 nalloc = 0;
174 time = NULL;
175 sbox = NULL;
176 sx = NULL;
177 clear_mat(avbox);
179 natom=read_first_x(oenv,&status,ftp2fn(efTRX,NFILE,fnm),&t,&x,box);
180 nfr = 0;
181 do {
182 if (nfr >= nalloc) {
183 nalloc += 100;
184 srenew(time,nalloc);
185 srenew(sbox,nalloc);
186 srenew(sx,nalloc);
189 time[nfr] = t;
190 copy_mat(box,sbox[nfr]);
191 /* This assumes that the off-diagonal box elements
192 * are not affected by jumps across the periodic boundaries.
194 m_add(avbox,box,avbox);
195 snew(sx[nfr],isize);
196 for(i=0; i<isize; i++)
197 copy_rvec(x[index[i]],sx[nfr][i]);
199 nfr++;
200 } while (read_next_x(oenv,status,&t,natom,x,box));
202 /* clean up */
203 sfree(x);
204 close_trj(status);
206 fprintf(stderr,"Read %d frames\n",nfr);
208 dt = (time[nfr-1] - time[0])/(nfr - 1);
209 /* Some ugly rounding to get nice nice times in the output */
210 dt = (int)(10000.0*dt + 0.5)/10000.0;
212 invbin = 1.0/rbin;
214 if (matfile) {
215 if (fmmax <= 0 || fmmax >= nfr)
216 fmmax = nfr - 1;
217 snew(mcount,fmmax);
218 nbin = (int)(rmax*invbin + 0.5);
219 if (sbin == 0) {
220 mat_nx = fmmax + 1;
221 } else {
222 invsbin = 1.0/sbin;
223 mat_nx = sqrt(fmmax*dt)*invsbin + 1;
225 snew(mat,mat_nx);
226 for(f=0; f<mat_nx; f++)
227 snew(mat[f],nbin);
228 rmax2 = sqr(nbin*rbin);
229 /* Initialize time zero */
230 mat[0][0] = nfr*isize;
231 mcount[0] += nfr;
232 } else {
233 fmmax = 0;
236 if (orfile) {
237 snew(pr,nr);
238 nalloc = 0;
239 snew(rcount,nr);
242 if (otfile) {
243 if (ftmax <= 0)
244 ftmax = nfr - 1;
245 snew(tcount,ftmax);
246 snew(pt,nfr);
247 rint2 = rint*rint;
248 /* Initialize time zero */
249 pt[0] = nfr*isize;
250 tcount[0] += nfr;
251 } else {
252 ftmax = 0;
255 msmul(avbox,1.0/nfr,avbox);
256 for(f=0; f<nfr; f++) {
257 if (f % 100 == 0)
258 fprintf(stderr,"\rProcessing frame %d",f);
259 /* Scale all the configuration to the average box */
260 m_inv_ur0(sbox[f],corr);
261 mmul_ur0(avbox,corr,corr);
262 for(i=0; i<isize; i++) {
263 mvmul_ur0(corr,sx[f][i],sx[f][i]);
264 if (f > 0) {
265 /* Correct for periodic jumps */
266 for(m=DIM-1; m>=0; m--) {
267 while(sx[f][i][m] - sx[f-1][i][m] > 0.5*avbox[m][m])
268 rvec_dec(sx[f][i],avbox[m]);
269 while(sx[f][i][m] - sx[f-1][i][m] <= -0.5*avbox[m][m])
270 rvec_inc(sx[f][i],avbox[m]);
274 for(ff=0; ff<f; ff++) {
275 fbin = f - ff;
276 if (fbin <= fmmax || fbin <= ftmax) {
277 if (sbin == 0)
278 mbin = fbin;
279 else
280 mbin = (int)(sqrt(fbin*dt)*invsbin + 0.5);
281 for(i=0; i<isize; i++) {
282 d2 = distance2(sx[f][i],sx[ff][i]);
283 if (mbin < mat_nx && d2 < rmax2) {
284 bin = (int)(sqrt(d2)*invbin + 0.5);
285 if (bin < nbin) {
286 mat[mbin][bin] += 1;
289 if (fbin <= ftmax && d2 <= rint2)
290 pt[fbin]++;
292 if (matfile)
293 mcount[mbin]++;
294 if (otfile)
295 tcount[fbin]++;
298 if (orfile) {
299 for(fbin=0; fbin<nr; fbin++) {
300 ff = f - (fbin + 1)*fshift;
301 if (ff >= 0) {
302 for(i=0; i<isize; i++) {
303 d2 = distance2(sx[f][i],sx[ff][i]);
304 bin = (int)(sqrt(d2)*invbin);
305 if (bin >= nalloc) {
306 nallocn = 10*(bin/10) + 11;
307 for(m=0; m<nr; m++) {
308 srenew(pr[m],nallocn);
309 for(i=nalloc; i<nallocn; i++)
310 pr[m][i] = 0;
312 nalloc = nallocn;
314 pr[fbin][bin]++;
316 rcount[fbin]++;
321 fprintf(stderr,"\n");
323 if (matfile) {
324 matmax = 0;
325 for(f=0; f<mat_nx; f++) {
326 normfac = 1.0/(mcount[f]*isize*rbin);
327 for(i=0; i<nbin; i++) {
328 mat[f][i] *= normfac;
329 if (mat[f][i] > matmax && (f!=0 || i!=0))
330 matmax = mat[f][i];
333 fprintf(stdout,"Value at (0,0): %.3f, maximum of the rest %.3f\n",
334 mat[0][0],matmax);
335 if (mmax > 0)
336 matmax = mmax;
337 snew(tickx,mat_nx);
338 for(f=0; f<mat_nx; f++) {
339 if (sbin == 0)
340 tickx[f] = f*dt;
341 else
342 tickx[f] = f*sbin;
344 snew(ticky,nbin+1);
345 for(i=0; i<=nbin; i++)
346 ticky[i] = i*rbin;
347 fp = ffopen(matfile,"w");
348 write_xpm(fp,MAT_SPATIAL_Y,"Van Hove function","G (1/nm)",
349 sbin==0 ? "time (ps)" : "sqrt(time) (ps^1/2)","r (nm)",
350 mat_nx,nbin,tickx,ticky,mat,0,matmax,rlo,rhi,&nlev);
351 ffclose(fp);
354 if (orfile) {
355 fp = xvgropen(orfile,"Van Hove function","r (nm)","G (nm\\S-1\\N)",oenv);
356 fprintf(fp,"@ subtitle \"for particles in group %s\"\n",grpname);
357 snew(legend,nr);
358 for(fbin=0; fbin<nr; fbin++) {
359 sprintf(buf,"%g ps",(fbin + 1)*fshift*dt);
360 legend[fbin] = strdup(buf);
362 xvgr_legend(fp,nr,legend,oenv);
363 for(i=0; i<nalloc; i++) {
364 fprintf(fp,"%g",i*rbin);
365 for(fbin=0; fbin<nr; fbin++)
366 fprintf(fp," %g",
367 (real)pr[fbin][i]/(rcount[fbin]*isize*rbin*(i==0 ? 0.5 : 1)));
368 fprintf(fp,"\n");
370 ffclose(fp);
373 if (otfile) {
374 sprintf(buf,"Probability of moving less than %g nm",rint);
375 fp = xvgropen(otfile,buf,"t (ps)","",oenv);
376 fprintf(fp,"@ subtitle \"for particles in group %s\"\n",grpname);
377 for(f=0; f<=ftmax; f++)
378 fprintf(fp,"%g %g\n",f*dt,(real)pt[f]/(tcount[f]*isize));
379 ffclose(fp);
382 do_view(oenv, matfile,NULL);
383 do_view(oenv, orfile,NULL);
384 do_view(oenv, otfile,NULL);
386 thanx(stderr);
388 return 0;