3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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
33 * Green Red Orange Magenta Azure Cyan Skyblue
51 #include "gmx_fatal.h"
58 static void dump_dih_trn(int nframes,int nangles,real **dih,const char *fn,
64 matrix box = {{2,0,0},{0,2,0},{0,0,2}};
71 printf("There are %d dihedrals. Will fill %d atom positions with cos/sin\n",
74 trn = open_trn(fn,"w");
75 for(i=0; (i<nframes); i++) {
77 for(j=0; (j<nangles); j++) {
78 for(m=0; (m<2); m++) {
79 x[k][l] = (m == 0) ? cos(dih[j][i]) : sin(dih[j][i]);
87 fwrite_trn(trn,i,time[i],0,box,na,x,NULL,NULL);
93 int gmx_g_angle(int argc,char *argv[])
95 static const char *desc[] = {
96 "[TT]g_angle[tt] computes the angle distribution for a number of angles",
97 "or dihedrals. This way you can check whether your simulation",
98 "is correct. With option [TT]-ov[tt] you can plot the average angle of",
99 "a group of angles as a function of time. With the [TT]-all[tt] option",
100 "the first graph is the average, the rest are the individual angles.[PAR]",
101 "With the [TT]-of[tt] option, [TT]g_angle[tt] also calculates the fraction of trans",
102 "dihedrals (only for dihedrals) as function of time, but this is",
103 "probably only fun for a selected few.[PAR]",
104 "With option [TT]-oc[tt] a dihedral correlation function is calculated.[PAR]",
105 "It should be noted that the index file should contain",
106 "atom-triples for angles or atom-quadruplets for dihedrals.",
107 "If this is not the case, the program will crash.[PAR]",
108 "With option [TT]-or[tt] a trajectory file is dumped containing cos and",
109 "sin of selected dihedral angles which subsequently can be used as",
110 "input for a PCA analysis using [TT]g_covar[tt].[PAR]",
111 "Option [TT]-ot[tt] plots when transitions occur between",
112 "dihedral rotamers of multiplicity 3 and [TT]-oh[tt]",
113 "records a histogram of the times between such transitions,",
114 "assuming the input trajectory frames are equally spaced in time."
116 static const char *opt[] = { NULL, "angle", "dihedral", "improper", "ryckaert-bellemans", NULL };
117 static gmx_bool bALL=FALSE,bChandler=FALSE,bAverCorr=FALSE,bPBC=TRUE;
118 static real binwidth=1;
120 { "-type", FALSE, etENUM, {opt},
121 "Type of angle to analyse" },
122 { "-all", FALSE, etBOOL, {&bALL},
123 "Plot all angles separately in the averages file, in the order of appearance in the index file." },
124 { "-binwidth", FALSE, etREAL, {&binwidth},
125 "binwidth (degrees) for calculating the distribution" },
126 { "-periodic", FALSE, etBOOL, {&bPBC},
127 "Print dihedral angles modulo 360 degrees" },
128 { "-chandler", FALSE, etBOOL, {&bChandler},
129 "Use Chandler correlation function (N[trans] = 1, N[gauche] = 0) rather than cosine correlation function. Trans is defined as phi < -60 or phi > 60." },
130 { "-avercorr", FALSE, etBOOL, {&bAverCorr},
131 "Average the correlation functions for the individual angles/dihedrals" }
133 static const char *bugs[] = {
134 "Counting transitions only works for dihedrals with multiplicity 3"
142 real maxang,Jc,S2,norm_fac,maxstat;
144 int nframes,maxangstat,mult,*angstat;
145 int i,j,total,nangles,natoms,nat2,first,last,angind;
146 gmx_bool bAver,bRb,bPeriodic,
147 bFrac, /* calculate fraction too? */
148 bTrans, /* worry about transtions too? */
149 bCorr; /* correlation function ? */
150 real t,aa,aver,aver2,aversig,fraction; /* fraction trans dihedrals */
153 real **dih=NULL; /* mega array with all dih. angles at all times*/
155 real *time,*trans_frac,*aver_angle;
157 { efTRX, "-f", NULL, ffREAD },
158 { efNDX, NULL, "angle", ffREAD },
159 { efXVG, "-od", "angdist", ffWRITE },
160 { efXVG, "-ov", "angaver", ffOPTWR },
161 { efXVG, "-of", "dihfrac", ffOPTWR },
162 { efXVG, "-ot", "dihtrans", ffOPTWR },
163 { efXVG, "-oh", "trhisto", ffOPTWR },
164 { efXVG, "-oc", "dihcorr", ffOPTWR },
165 { efTRR, "-or", NULL, ffOPTWR }
167 #define NFILE asize(fnm)
172 CopyRight(stderr,argv[0]);
174 ppa = add_acf_pargs(&npargs,pa);
175 parse_common_args(&argc,argv,PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
176 NFILE,fnm,npargs,ppa,asize(desc),desc,asize(bugs),bugs,
196 if (opt2bSet("-or",NFILE,fnm)) {
198 gmx_fatal(FARGS,"Can not combine angles with trn dump");
200 please_cite(stdout,"Mu2005a");
203 /* Calculate bin size */
204 maxangstat=(int)(maxang/binwidth+0.5);
205 binwidth=maxang/maxangstat;
207 rd_index(ftp2fn(efNDX,NFILE,fnm),1,&isize,&index,&grpname);
209 if ((isize % mult) != 0)
210 gmx_fatal(FARGS,"number of index elements not multiple of %d, "
211 "these can not be %s\n",
212 mult,(mult==3) ? "angle triplets" : "dihedral quadruplets");
215 /* Check whether specific analysis has to be performed */
216 bCorr=opt2bSet("-oc",NFILE,fnm);
217 bAver=opt2bSet("-ov",NFILE,fnm);
218 bTrans=opt2bSet("-ot",NFILE,fnm);
219 bFrac=opt2bSet("-of",NFILE,fnm);
220 if (bTrans && opt[0][0] != 'd') {
221 fprintf(stderr, "Option -ot should only accompany -type dihedral. Disabling -ot.\n");
225 if (bChandler && !bCorr)
229 fprintf(stderr,"Warning:"
230 " calculating fractions as defined in this program\n"
231 "makes sense for Ryckaert Bellemans dihs. only. Ignoring -of\n\n");
235 if ( (bTrans || bFrac || bCorr) && mult==3)
236 gmx_fatal(FARGS,"Can only do transition, fraction or correlation\n"
237 "on dihedrals. Select -d\n");
240 * We need to know the nr of frames so we can allocate memory for an array
241 * with all dihedral angles at all timesteps. Works for me.
243 if (bTrans || bCorr || bALL || opt2bSet("-or",NFILE,fnm))
246 snew(angstat,maxangstat);
248 read_ang_dih(ftp2fn(efTRX,NFILE,fnm),(mult == 3),
249 bALL || bCorr || bTrans || opt2bSet("-or",NFILE,fnm),
250 bRb,bPBC,maxangstat,angstat,
251 &nframes,&time,isize,index,&trans_frac,&aver_angle,dih,
254 dt=(time[nframes-1]-time[0])/(nframes-1);
257 sprintf(title,"Average Angle: %s",grpname);
258 out=xvgropen(opt2fn("-ov",NFILE,fnm),
259 title,"Time (ps)","Angle (degrees)",oenv);
260 for(i=0; (i<nframes); i++) {
261 fprintf(out,"%10.5f %8.3f",time[i],aver_angle[i]*RAD2DEG);
263 for(j=0; (j<nangles); j++)
266 fprintf(out," %8.3f",atan2(sin(dd),cos(dd))*RAD2DEG);
269 fprintf(out," %8.3f",dih[j][i]*RAD2DEG);
275 if (opt2bSet("-or",NFILE,fnm))
276 dump_dih_trn(nframes,nangles,dih,opt2fn("-or",NFILE,fnm),time);
279 sprintf(title,"Trans fraction: %s",grpname);
280 out=xvgropen(opt2fn("-of",NFILE,fnm),
281 title,"Time (ps)","Fraction",oenv);
283 for(i=0; (i<nframes); i++) {
284 fprintf(out,"%10.5f %10.3f\n",time[i],trans_frac[i]);
285 tfrac += trans_frac[i];
290 fprintf(stderr,"Average trans fraction: %g\n",tfrac);
295 ana_dih_trans(opt2fn("-ot",NFILE,fnm),opt2fn("-oh",NFILE,fnm),
296 dih,nframes,nangles,grpname,time,bRb,oenv);
299 /* Autocorrelation function */
301 fprintf(stderr,"Not enough frames for correlation function\n");
305 real dval,sixty=DEG2RAD*60;
308 for(i=0; (i<nangles); i++)
309 for(j=0; (j<nframes); j++) {
312 bTest=(dval > -sixty) && (dval < sixty);
314 bTest=(dval < -sixty) || (dval > sixty);
316 dih[i][j] = dval-tfrac;
325 do_autocorr(opt2fn("-oc",NFILE,fnm), oenv,
326 "Dihedral Autocorrelation Function",
327 nframes,nangles,dih,dt,mode,bAverCorr);
332 /* Determine the non-zero part of the distribution */
333 for(first=0; (first < maxangstat-1) && (angstat[first+1] == 0); first++)
335 for(last=maxangstat-1; (last > 0) && (angstat[last-1] == 0) ; last--)
339 for(i=0; (i<nframes); i++) {
340 aver += RAD2DEG*aver_angle[i];
341 aver2 += sqr(RAD2DEG*aver_angle[i]);
343 aver /= (real) nframes;
344 aver2 /= (real) nframes;
345 aversig = sqrt(aver2-sqr(aver));
346 printf("Found points in the range from %d to %d (max %d)\n",
347 first,last,maxangstat);
348 printf(" < angle > = %g\n",aver);
349 printf("< angle^2 > = %g\n",aver2);
350 printf("Std. Dev. = %g\n",aversig);
353 sprintf(title,"Angle Distribution: %s",grpname);
355 sprintf(title,"Dihedral Distribution: %s",grpname);
357 calc_distribution_props(maxangstat,angstat,-180.0,0,NULL,&S2);
358 fprintf(stderr,"Order parameter S^2 = %g\n",S2);
361 bPeriodic=(mult==4) && (first==0) && (last==maxangstat-1);
363 out=xvgropen(opt2fn("-od",NFILE,fnm),title,"Degrees","",oenv);
364 if (output_env_get_print_xvgr_codes(oenv))
365 fprintf(out,"@ subtitle \"average angle: %g\\So\\N\"\n",aver);
366 norm_fac=1.0/(nangles*nframes*binwidth);
369 for(i=first; (i<=last); i++)
370 maxstat=max(maxstat,angstat[i]*norm_fac);
371 fprintf(out,"@with g0\n");
372 fprintf(out,"@ world xmin -180\n");
373 fprintf(out,"@ world xmax 180\n");
374 fprintf(out,"@ world ymin 0\n");
375 fprintf(out,"@ world ymax %g\n",maxstat*1.05);
376 fprintf(out,"@ xaxis tick major 60\n");
377 fprintf(out,"@ xaxis tick minor 30\n");
378 fprintf(out,"@ yaxis tick major 0.005\n");
379 fprintf(out,"@ yaxis tick minor 0.0025\n");
381 for(i=first; (i<=last); i++)
382 fprintf(out,"%10g %10f\n",i*binwidth+180.0-maxang,angstat[i]*norm_fac);
384 /* print first bin again as last one */
385 fprintf(out,"%10g %10f\n",180.0,angstat[0]*norm_fac);
389 do_view(oenv,opt2fn("-od",NFILE,fnm),"-nxy");
391 do_view(oenv,opt2fn("-ov",NFILE,fnm),"-nxy");