af00d9175755ad7b0e2051ab769aa95cc57864bb
[alexxy/gromacs.git] / src / tools / gmx_angle.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 #include <math.h>
39 #include <string.h>
40
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"
56
57
58 static void dump_dih_trn(int nframes,int nangles,real **dih,const char *fn,
59                          real *time)
60 {
61   int    i,j,k,l,m,na;
62   t_fileio *trn;
63   rvec   *x;
64   matrix box = {{2,0,0},{0,2,0},{0,0,2}};  
65   
66   na = (nangles*2);
67   if ((na % 3) != 0)
68     na = 1+na/3;
69   else
70     na = na/3;
71   printf("There are %d dihedrals. Will fill %d atom positions with cos/sin\n",
72          nangles,na);
73   snew(x,na);
74   trn = open_trn(fn,"w");
75   for(i=0; (i<nframes); i++) {
76     k = l = 0;
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]);
80         l++;
81         if (l == DIM) {
82           l = 0;
83           k++;
84         }
85       }
86     }
87     fwrite_trn(trn,i,time[i],0,box,na,x,NULL,NULL);
88   }
89   close_trn(trn);
90   sfree(x);
91 }
92
93 int gmx_g_angle(int argc,char *argv[])
94 {
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."
115   };
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;
119   t_pargs pa[] = {
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" }
132   };
133   static const char *bugs[] = {
134     "Counting transitions only works for dihedrals with multiplicity 3"
135   };
136   
137   FILE       *out;
138   real       tmp,dt;
139   int        status,isize;
140   atom_id    *index;
141   char       *grpname;
142   real       maxang,Jc,S2,norm_fac,maxstat;
143   unsigned long mode;
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 */
151   double     tfrac=0;
152   char       title[256];
153   real       **dih=NULL;          /* mega array with all dih. angles at all times*/
154   char       buf[80];       
155   real       *time,*trans_frac,*aver_angle;
156   t_filenm   fnm[] = {
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 }
166   };
167 #define NFILE asize(fnm)
168   int     npargs;
169   t_pargs *ppa;
170   output_env_t oenv;
171   
172   CopyRight(stderr,argv[0]);
173   npargs = asize(pa);
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,
177                     &oenv);
178                     
179   mult   = 4;
180   maxang = 360.0;
181   bRb    = FALSE;
182   switch(opt[0][0]) {
183   case 'a':
184     mult   = 3;
185     maxang = 180.0;
186     break;
187   case 'd':
188     break;
189   case 'i':
190     break;
191   case 'r':
192     bRb = TRUE;
193     break;
194   }
195
196   if (opt2bSet("-or",NFILE,fnm)) {
197     if (mult != 4)
198       gmx_fatal(FARGS,"Can not combine angles with trn dump");
199     else
200       please_cite(stdout,"Mu2005a");
201   }
202     
203   /* Calculate bin size */
204   maxangstat=(int)(maxang/binwidth+0.5);
205   binwidth=maxang/maxangstat;
206     
207   rd_index(ftp2fn(efNDX,NFILE,fnm),1,&isize,&index,&grpname);
208   nangles=isize/mult;
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");
213   
214
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");
222     bTrans = FALSE;
223   }
224   
225   if (bChandler && !bCorr)
226     bCorr=TRUE;
227     
228   if (bFrac && !bRb) {
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"); 
232     bFrac = FALSE;
233   }
234   
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");
238   
239   /* 
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.
242    */
243   if (bTrans || bCorr  || bALL || opt2bSet("-or",NFILE,fnm))
244     snew(dih,nangles);
245   
246   snew(angstat,maxangstat);
247
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,
252                oenv);
253   
254   dt=(time[nframes-1]-time[0])/(nframes-1);
255   
256   if (bAver) {
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);
262       if (bALL) {
263         for(j=0; (j<nangles); j++)
264           if (bPBC) {
265             real dd = dih[j][i];
266             fprintf(out,"  %8.3f",atan2(sin(dd),cos(dd))*RAD2DEG);
267           }
268           else
269             fprintf(out,"  %8.3f",dih[j][i]*RAD2DEG);
270       }
271       fprintf(out,"\n");
272     }   
273     ffclose(out);
274   }
275   if (opt2bSet("-or",NFILE,fnm)) 
276     dump_dih_trn(nframes,nangles,dih,opt2fn("-or",NFILE,fnm),time);
277   
278   if (bFrac) {
279     sprintf(title,"Trans fraction: %s",grpname);
280     out=xvgropen(opt2fn("-of",NFILE,fnm),
281                   title,"Time (ps)","Fraction",oenv);
282     tfrac = 0.0;
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];
286     }
287     ffclose(out);
288     
289     tfrac/=nframes;
290     fprintf(stderr,"Average trans fraction: %g\n",tfrac);
291   }
292   sfree(trans_frac);
293   
294   if (bTrans) 
295     ana_dih_trans(opt2fn("-ot",NFILE,fnm),opt2fn("-oh",NFILE,fnm),
296                   dih,nframes,nangles,grpname,time,bRb,oenv);
297                   
298   if (bCorr) {
299     /* Autocorrelation function */
300     if (nframes < 2)
301       fprintf(stderr,"Not enough frames for correlation function\n");
302     else {
303       
304       if (bChandler) {
305         real dval,sixty=DEG2RAD*60;
306         gmx_bool bTest;
307
308         for(i=0; (i<nangles); i++)
309           for(j=0; (j<nframes); j++) {
310             dval = dih[i][j];
311             if (bRb)
312               bTest=(dval > -sixty) && (dval < sixty);
313             else
314               bTest=(dval < -sixty) || (dval > sixty);
315             if (bTest)
316               dih[i][j] = dval-tfrac;
317             else
318               dih[i][j] = -tfrac;
319           }
320       }
321       if (bChandler)
322         mode = eacNormal;
323       else
324         mode = eacCos;
325       do_autocorr(opt2fn("-oc",NFILE,fnm), oenv,
326                   "Dihedral Autocorrelation Function",
327                   nframes,nangles,dih,dt,mode,bAverCorr);
328     }
329   }
330
331   
332   /* Determine the non-zero part of the distribution */
333   for(first=0; (first < maxangstat-1) && (angstat[first+1] == 0); first++)
334     ;
335   for(last=maxangstat-1; (last > 0) && (angstat[last-1] == 0) ; last--)
336     ;
337
338   aver=aver2=0;
339   for(i=0; (i<nframes); i++) {
340     aver  += RAD2DEG*aver_angle[i];
341     aver2 += sqr(RAD2DEG*aver_angle[i]);
342   }
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);
351     
352   if (mult == 3)
353     sprintf(title,"Angle Distribution: %s",grpname);
354   else {
355     sprintf(title,"Dihedral Distribution: %s",grpname);
356     
357     calc_distribution_props(maxangstat,angstat,-180.0,0,NULL,&S2);
358     fprintf(stderr,"Order parameter S^2 = %g\n",S2);
359   }
360   
361   bPeriodic=(mult==4) && (first==0) && (last==maxangstat-1);
362   
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);
367   if (bPeriodic) {
368     maxstat=0;
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");
380   }
381   for(i=first; (i<=last); i++) 
382     fprintf(out,"%10g  %10f\n",i*binwidth+180.0-maxang,angstat[i]*norm_fac);
383   if ( bPeriodic )
384     /* print first bin again as last one */
385     fprintf(out,"%10g  %10f\n",180.0,angstat[0]*norm_fac);
386   
387   ffclose(out);
388
389   do_view(oenv,opt2fn("-od",NFILE,fnm),"-nxy");
390   if (bAver)
391     do_view(oenv,opt2fn("-ov",NFILE,fnm),"-nxy");
392     
393   thanx(stderr);
394     
395   return 0;
396 }