Redefine the default boolean type to gmx_bool.
[alexxy/gromacs.git] / src / tools / gmx_rotacf.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 <math.h>
40 #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 "index.h"
49 #include "macros.h"
50 #include "gmx_fatal.h"
51 #include "xvgr.h"
52 #include "gstat.h"
53 #include "vec.h"
54 #include "gmx_ana.h"
55
56
57 int gmx_rotacf(int argc,char *argv[])
58 {
59   const char *desc[] = {
60     "g_rotacf calculates the rotational correlation function",
61     "for molecules. Three atoms (i,j,k) must be given in the index",
62     "file, defining two vectors ij and jk. The rotational acf",
63     "is calculated as the autocorrelation function of the vector",
64     "n = ij x jk, i.e. the cross product of the two vectors.",
65     "Since three atoms span a plane, the order of the three atoms",
66     "does not matter. Optionally, controlled by the -d switch, you can",
67     "calculate the rotational correlation function for linear molecules",
68     "by specifying two atoms (i,j) in the index file.",
69     "[PAR]",
70     "EXAMPLES[PAR]",
71     "g_rotacf -P 1 -nparm 2 -fft -n index -o rotacf-x-P1",
72     "-fa expfit-x-P1 -beginfit 2.5 -endfit 20.0[PAR]",
73     "This will calculate the rotational correlation function using a first",
74     "order Legendre polynomial of the angle of a vector defined by the index",
75     "file. The correlation function will be fitted from 2.5 ps till 20.0 ps",
76     "to a two parameter exponential.",
77
78
79     ""
80   };
81   static gmx_bool bVec    = FALSE,bAver=TRUE;
82
83   t_pargs pa[] = {
84     { "-d",   FALSE, etBOOL, {&bVec},
85       "Use index doublets (vectors) for correlation function instead of triplets (planes)" },
86     { "-aver",FALSE, etBOOL, {&bAver},
87       "Average over molecules" }
88   };
89
90   t_trxstatus *status;
91   int        isize;
92   atom_id    *index;
93   char       *grpname;
94   rvec       *x,*x_s;
95   matrix     box;
96   real       **c1;
97   rvec       xij,xjk,n;
98   int        i,m,teller,n_alloc,natoms,nvec,ai,aj,ak;
99   unsigned long mode;
100   real       t,t0,t1,dt;
101   gmx_rmpbc_t  gpbc=NULL;
102   t_topology *top;
103   int        ePBC;
104   t_filenm   fnm[] = {
105     { efTRX, "-f", NULL,  ffREAD  },
106     { efTPX, NULL, NULL,  ffREAD },
107     { efNDX, NULL, NULL,  ffREAD  },
108     { efXVG, "-o", "rotacf",  ffWRITE }
109   };
110 #define NFILE asize(fnm)
111   int     npargs;
112   t_pargs *ppa;
113
114   output_env_t oenv;
115   
116   CopyRight(stderr,argv[0]);
117   npargs = asize(pa);
118   ppa    = add_acf_pargs(&npargs,pa);
119   
120   parse_common_args(&argc,argv,PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
121                     NFILE,fnm,npargs,ppa,asize(desc),desc,0,NULL,&oenv);
122   
123   rd_index(ftp2fn(efNDX,NFILE,fnm),1,&isize,&index,&grpname);
124   
125   if (bVec) 
126     nvec = isize/2;
127   else
128     nvec = isize/3;
129   
130   if (((isize % 3) != 0) && !bVec)
131     gmx_fatal(FARGS,"number of index elements not multiple of 3, "
132                 "these can not be atom triplets\n");
133   if (((isize % 2) != 0) && bVec)
134     gmx_fatal(FARGS,"number of index elements not multiple of 2, "
135                 "these can not be atom doublets\n");
136   
137   top=read_top(ftp2fn(efTPX,NFILE,fnm),&ePBC);
138   
139   snew(c1,nvec);
140   for (i=0; (i<nvec); i++)
141     c1[i]=NULL;
142   n_alloc=0;
143
144   natoms=read_first_x(oenv,&status,ftp2fn(efTRX,NFILE,fnm),&t,&x,box);
145   snew(x_s,natoms);
146   
147   gpbc = gmx_rmpbc_init(&(top->idef),ePBC,natoms,box);
148
149   /* Start the loop over frames */
150   t1 = t0 = t;
151   teller  = 0;
152   do {
153     if (teller >= n_alloc) {
154       n_alloc+=100;
155       for (i=0; (i<nvec); i++)
156         srenew(c1[i],DIM*n_alloc);
157     }
158     t1 = t;
159     
160     /* Remove periodicity */
161     gmx_rmpbc_copy(gpbc,natoms,box,x,x_s);
162   
163     /* Compute crossproducts for all vectors, if triplets.
164      * else, just get the vectors in case of doublets.
165      */
166     if (bVec == FALSE) {
167       for (i=0; (i<nvec); i++) {
168         ai=index[3*i];
169         aj=index[3*i+1];
170         ak=index[3*i+2];
171         rvec_sub(x_s[ai],x_s[aj],xij);
172         rvec_sub(x_s[aj],x_s[ak],xjk);
173         cprod(xij,xjk,n);
174         for(m=0; (m<DIM); m++)
175           c1[i][DIM*teller+m]=n[m];
176       }
177     }
178     else {
179       for (i=0; (i<nvec); i++) {
180         ai=index[2*i];
181         aj=index[2*i+1];
182         rvec_sub(x_s[ai],x_s[aj],n);
183         for(m=0; (m<DIM); m++)
184           c1[i][DIM*teller+m]=n[m];
185       }
186     }
187     /* Increment loop counter */
188     teller++;
189   } while (read_next_x(oenv,status,&t,natoms,x,box));  
190   close_trj(status); 
191   fprintf(stderr,"\nDone with trajectory\n");
192   
193   gmx_rmpbc_done(gpbc);
194
195
196   /* Autocorrelation function */
197   if (teller < 2)
198     fprintf(stderr,"Not enough frames for correlation function\n");
199   else {
200     dt=(t1 - t0)/(teller-1);
201     
202     mode = eacVector;
203     
204     do_autocorr(ftp2fn(efXVG,NFILE,fnm),oenv,"Rotational Correlation Function",
205                 teller,nvec,c1,dt,mode,bAver);
206   }
207
208   do_view(oenv,ftp2fn(efXVG,NFILE,fnm),NULL);
209     
210   thanx(stderr);
211     
212   return 0;
213 }