minor fixes in ditribution files
[gromacs/qmmm-gamess-us.git] / src / tools / gmx_angle.c
blob0c2d69ed9162699ad00818b88b1a9f764f4b9686
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
38 #include <math.h>
39 #include <string.h>
41 #include "sysstuff.h"
42 #include "physics.h"
43 #include "typedefs.h"
44 #include "smalloc.h"
45 #include "futil.h"
46 #include "statutil.h"
47 #include "copyrite.h"
48 #include "vec.h"
49 #include "index.h"
50 #include "macros.h"
51 #include "gmx_fatal.h"
52 #include "xvgr.h"
53 #include "gstat.h"
54 #include "trnio.h"
55 #include "gmx_ana.h"
58 static void dump_dih_trn(int nframes,int nangles,real **dih,const char *fn,
59 real dt)
61 int i,j,k,l,m,na,trn;
62 rvec *x;
63 matrix box = {{2,0,0},{0,2,0},{0,0,2}};
65 na = (nangles*2);
66 if ((na % 3) != 0)
67 na = 1+na/3;
68 else
69 na = na/3;
70 printf("There are %d dihedrals. Will fill %d atom positions with cos/sin\n",
71 nangles,na);
72 snew(x,na);
73 trn = open_trn(fn,"w");
74 for(i=0; (i<nframes); i++) {
75 k = l = 0;
76 for(j=0; (j<nangles); j++) {
77 for(m=0; (m<2); m++) {
78 x[k][l] = (m == 0) ? cos(dih[j][i]) : sin(dih[j][i]);
79 l++;
80 if (l == DIM) {
81 l = 0;
82 k++;
86 fwrite_trn(trn,i,(real)i*dt,0,box,na,x,NULL,NULL);
88 close_trn(trn);
89 sfree(x);
92 int gmx_g_angle(int argc,char *argv[])
94 static const char *desc[] = {
95 "g_angle computes the angle distribution for a number of angles",
96 "or dihedrals. This way you can check whether your simulation",
97 "is correct. With option -ov you can plot the average angle of",
98 "a group of angles as a function of time. With the -all option",
99 "the first graph is the average, the rest are the individual angles.[PAR]",
100 "With the -of option g_angle also calculates the fraction of trans",
101 "dihedrals (only for dihedrals) as function of time, but this is",
102 "probably only fun for a selected few.[PAR]",
103 "With option -oc a dihedral correlation function is calculated.[PAR]",
104 "It should be noted that the indexfile should contain",
105 "atom-triples for angles or atom-quadruplets for dihedrals.",
106 "If this is not the case, the program will crash.[PAR]",
107 "With option [TT]-or[tt] a trajectory file is dumped containing cos and"
108 "sin of selected dihedral angles which subsequently can be used as",
109 "input for a PCA analysis using [TT]g_covar[tt]."
111 static const char *opt[] = { NULL, "angle", "dihedral", "improper", "ryckaert-bellemans", NULL };
112 static bool bALL=FALSE,bChandler=FALSE,bAverCorr=FALSE,bPBC=TRUE;
113 static real binwidth=1;
114 t_pargs pa[] = {
115 { "-type", FALSE, etENUM, {opt},
116 "Type of angle to analyse" },
117 { "-all", FALSE, etBOOL, {&bALL},
118 "Plot all angles separately in the averages file, in the order of appearance in the index file." },
119 { "-binwidth", FALSE, etREAL, {&binwidth},
120 "binwidth (degrees) for calculating the distribution" },
121 { "-periodic", FALSE, etBOOL, {&bPBC},
122 "Print dihedral angles modulo 360 degrees" },
123 { "-chandler", FALSE, etBOOL, {&bChandler},
124 "Use Chandler correlation function (N[trans] = 1, N[gauche] = 0) rather than cosine correlation function. Trans is defined as phi < -60 || phi > 60." },
125 { "-avercorr", FALSE, etBOOL, {&bAverCorr},
126 "Average the correlation functions for the individual angles/dihedrals" }
128 static const char *bugs[] = {
129 "Counting transitions only works for dihedrals with multiplicity 3"
132 FILE *out;
133 real tmp,dt;
134 int status,isize;
135 atom_id *index;
136 char *grpname;
137 real maxang,Jc,S2,norm_fac,maxstat;
138 unsigned long mode;
139 int nframes,maxangstat,mult,*angstat;
140 int i,j,total,nangles,natoms,nat2,first,last,angind;
141 bool bAver,bRb,bPeriodic,
142 bFrac, /* calculate fraction too? */
143 bTrans, /* worry about transtions too? */
144 bCorr; /* correlation function ? */
145 real t,aa,aver,aver2,aversig,fraction; /* fraction trans dihedrals */
146 double tfrac=0;
147 char title[256];
148 real **dih=NULL; /* mega array with all dih. angles at all times*/
149 char buf[80];
150 real *time,*trans_frac,*aver_angle;
151 t_filenm fnm[] = {
152 { efTRX, "-f", NULL, ffREAD },
153 { efNDX, NULL, "angle", ffREAD },
154 { efXVG, "-od", "angdist", ffWRITE },
155 { efXVG, "-ov", "angaver", ffOPTWR },
156 { efXVG, "-of", "dihfrac", ffOPTWR },
157 { efXVG, "-ot", "dihtrans", ffOPTWR },
158 { efXVG, "-oh", "trhisto", ffOPTWR },
159 { efXVG, "-oc", "dihcorr", ffOPTWR },
160 { efTRR, "-or", NULL, ffOPTWR }
162 #define NFILE asize(fnm)
163 int npargs;
164 t_pargs *ppa;
165 output_env_t oenv;
167 CopyRight(stderr,argv[0]);
168 npargs = asize(pa);
169 ppa = add_acf_pargs(&npargs,pa);
170 parse_common_args(&argc,argv,PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
171 NFILE,fnm,npargs,ppa,asize(desc),desc,asize(bugs),bugs,
172 &oenv);
174 mult = 4;
175 maxang = 360.0;
176 bRb = FALSE;
177 switch(opt[0][0]) {
178 case 'a':
179 mult = 3;
180 maxang = 180.0;
181 break;
182 case 'd':
183 break;
184 case 'i':
185 break;
186 case 'r':
187 bRb = TRUE;
188 break;
191 if (opt2bSet("-or",NFILE,fnm)) {
192 if (mult != 4)
193 gmx_fatal(FARGS,"Can not combine angles with trn dump");
194 else
195 please_cite(stdout,"Mu2005a");
198 /* Calculate bin size */
199 maxangstat=(int)(maxang/binwidth+0.5);
200 binwidth=maxang/maxangstat;
202 rd_index(ftp2fn(efNDX,NFILE,fnm),1,&isize,&index,&grpname);
203 nangles=isize/mult;
204 if ((isize % mult) != 0)
205 gmx_fatal(FARGS,"number of index elements not multiple of %d, "
206 "these can not be %s\n",
207 mult,(mult==3) ? "angle triplets" : "dihedral quadruplets");
210 /* Check whether specific analysis has to be performed */
211 bCorr=opt2bSet("-oc",NFILE,fnm);
212 bAver=opt2bSet("-ov",NFILE,fnm);
213 bTrans=opt2bSet("-ot",NFILE,fnm);
214 bFrac=opt2bSet("-of",NFILE,fnm);
216 if (bChandler && !bCorr)
217 bCorr=TRUE;
219 if (bFrac && !bRb) {
220 fprintf(stderr,"Warning:"
221 " calculating fractions as defined in this program\n"
222 "makes sense for Ryckaert Bellemans dihs. only. Ignoring -of\n\n");
223 bFrac = FALSE;
226 if ( (bTrans || bFrac || bCorr) && mult==3)
227 gmx_fatal(FARGS,"Can only do transition, fraction or correlation\n"
228 "on dihedrals. Select -d\n");
231 * We need to know the nr of frames so we can allocate memory for an array
232 * with all dihedral angles at all timesteps. Works for me.
234 if (bTrans || bCorr || bALL || opt2bSet("-or",NFILE,fnm))
235 snew(dih,nangles);
237 snew(angstat,maxangstat);
239 read_ang_dih(ftp2fn(efTRX,NFILE,fnm),(mult == 3),
240 bALL || bCorr || bTrans || opt2bSet("-or",NFILE,fnm),
241 bRb,bPBC,maxangstat,angstat,
242 &nframes,&time,isize,index,&trans_frac,&aver_angle,dih,
243 oenv);
245 dt=(time[nframes-1]-time[0])/(nframes-1);
247 if (bAver) {
248 sprintf(title,"Average Angle: %s",grpname);
249 out=xvgropen(opt2fn("-ov",NFILE,fnm),
250 title,"Time (ps)","Angle (degrees)",oenv);
251 for(i=0; (i<nframes); i++) {
252 fprintf(out,"%10.5f %8.3f",time[i],aver_angle[i]*RAD2DEG);
253 if (bALL) {
254 for(j=0; (j<nangles); j++)
255 if (bPBC) {
256 real dd = dih[j][i];
257 fprintf(out," %8.3f",atan2(sin(dd),cos(dd))*RAD2DEG);
259 else
260 fprintf(out," %8.3f",dih[j][i]*RAD2DEG);
262 fprintf(out,"\n");
264 ffclose(out);
266 if (opt2bSet("-or",NFILE,fnm))
267 dump_dih_trn(nframes,nangles,dih,opt2fn("-or",NFILE,fnm),dt);
269 if (bFrac) {
270 sprintf(title,"Trans fraction: %s",grpname);
271 out=xvgropen(opt2fn("-of",NFILE,fnm),
272 title,"Time (ps)","Fraction",oenv);
273 tfrac = 0.0;
274 for(i=0; (i<nframes); i++) {
275 fprintf(out,"%10.5f %10.3f\n",time[i],trans_frac[i]);
276 tfrac += trans_frac[i];
278 ffclose(out);
280 tfrac/=nframes;
281 fprintf(stderr,"Average trans fraction: %g\n",tfrac);
283 sfree(trans_frac);
285 if (bTrans)
286 ana_dih_trans(opt2fn("-ot",NFILE,fnm),opt2fn("-oh",NFILE,fnm),
287 dih,nframes,nangles,grpname,time[0],dt,bRb,oenv);
289 if (bCorr) {
290 /* Autocorrelation function */
291 if (nframes < 2)
292 fprintf(stderr,"Not enough frames for correlation function\n");
293 else {
295 if (bChandler) {
296 real dval,sixty=DEG2RAD*60;
297 bool bTest;
299 for(i=0; (i<nangles); i++)
300 for(j=0; (j<nframes); j++) {
301 dval = dih[i][j];
302 if (bRb)
303 bTest=(dval > -sixty) && (dval < sixty);
304 else
305 bTest=(dval < -sixty) || (dval > sixty);
306 if (bTest)
307 dih[i][j] = dval-tfrac;
308 else
309 dih[i][j] = -tfrac;
312 if (bChandler)
313 mode = eacNormal;
314 else
315 mode = eacCos;
316 do_autocorr(opt2fn("-oc",NFILE,fnm), oenv,
317 "Dihedral Autocorrelation Function",
318 nframes,nangles,dih,dt,mode,bAverCorr);
323 /* Determine the non-zero part of the distribution */
324 for(first=0; (first < maxangstat-1) && (angstat[first+1] == 0); first++)
326 for(last=maxangstat-1; (last > 0) && (angstat[last-1] == 0) ; last--)
329 aver=aver2=0;
330 for(i=0; (i<nframes); i++) {
331 aver += RAD2DEG*aver_angle[i];
332 aver2 += sqr(RAD2DEG*aver_angle[i]);
334 aver /= (real) nframes;
335 aver2 /= (real) nframes;
336 aversig = sqrt(aver2-sqr(aver));
337 printf("Found points in the range from %d to %d (max %d)\n",
338 first,last,maxangstat);
339 printf(" < angle > = %g\n",aver);
340 printf("< angle^2 > = %g\n",aver2);
341 printf("Std. Dev. = %g\n",aversig);
343 if (mult == 3)
344 sprintf(title,"Angle Distribution: %s",grpname);
345 else {
346 sprintf(title,"Dihedral Distribution: %s",grpname);
348 calc_distribution_props(maxangstat,angstat,-180.0,0,NULL,&S2);
349 fprintf(stderr,"Order parameter S^2 = %g\n",S2);
352 bPeriodic=(mult==4) && (first==0) && (last==maxangstat-1);
354 out=xvgropen(opt2fn("-od",NFILE,fnm),title,"Degrees","",oenv);
355 if (output_env_get_print_xvgr_codes(oenv))
356 fprintf(out,"@ subtitle \"average angle: %g\\So\\N\"\n",aver);
357 norm_fac=1.0/(nangles*nframes*binwidth);
358 if (bPeriodic) {
359 maxstat=0;
360 for(i=first; (i<=last); i++)
361 maxstat=max(maxstat,angstat[i]*norm_fac);
362 fprintf(out,"@with g0\n");
363 fprintf(out,"@ world xmin -180\n");
364 fprintf(out,"@ world xmax 180\n");
365 fprintf(out,"@ world ymin 0\n");
366 fprintf(out,"@ world ymax %g\n",maxstat*1.05);
367 fprintf(out,"@ xaxis tick major 60\n");
368 fprintf(out,"@ xaxis tick minor 30\n");
369 fprintf(out,"@ yaxis tick major 0.005\n");
370 fprintf(out,"@ yaxis tick minor 0.0025\n");
372 for(i=first; (i<=last); i++)
373 fprintf(out,"%10g %10f\n",i*binwidth+180.0-maxang,angstat[i]*norm_fac);
374 if ( bPeriodic )
375 /* print first bin again as last one */
376 fprintf(out,"%10g %10f\n",180.0,angstat[0]*norm_fac);
378 ffclose(out);
380 do_view(oenv,opt2fn("-od",NFILE,fnm),"-nxy");
381 if (bAver)
382 do_view(oenv,opt2fn("-ov",NFILE,fnm),"-nxy");
384 thanx(stderr);
386 return 0;