Fixed make_edi.c
[gromacs/rigid-bodies.git] / src / tools / gmx_rotacf.c
blobbeedaabfc2b6a56e89214a132842d2251c24c5c0
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 <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"
57 int gmx_rotacf(int argc,char *argv[])
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.",
81 static bool bVec = FALSE,bAver=TRUE;
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" }
90 int status,isize;
91 atom_id *index;
92 char *grpname;
93 rvec *x,*x_s;
94 matrix box;
95 real **c1;
96 rvec xij,xjk,n;
97 int i,m,teller,n_alloc,natoms,nvec,ai,aj,ak;
98 unsigned long mode;
99 real t,t0,t1,dt;
100 t_topology *top;
101 int ePBC;
102 t_filenm fnm[] = {
103 { efTRX, "-f", NULL, ffREAD },
104 { efTPX, NULL, NULL, ffREAD },
105 { efNDX, NULL, NULL, ffREAD },
106 { efXVG, "-o", "rotacf", ffWRITE }
108 #define NFILE asize(fnm)
109 int npargs;
110 t_pargs *ppa;
112 output_env_t oenv;
114 CopyRight(stderr,argv[0]);
115 npargs = asize(pa);
116 ppa = add_acf_pargs(&npargs,pa);
118 parse_common_args(&argc,argv,PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
119 NFILE,fnm,npargs,ppa,asize(desc),desc,0,NULL,&oenv);
121 rd_index(ftp2fn(efNDX,NFILE,fnm),1,&isize,&index,&grpname);
123 if (bVec)
124 nvec = isize/2;
125 else
126 nvec = isize/3;
128 if (((isize % 3) != 0) && !bVec)
129 gmx_fatal(FARGS,"number of index elements not multiple of 3, "
130 "these can not be atom triplets\n");
131 if (((isize % 2) != 0) && bVec)
132 gmx_fatal(FARGS,"number of index elements not multiple of 2, "
133 "these can not be atom doublets\n");
135 top=read_top(ftp2fn(efTPX,NFILE,fnm),&ePBC);
137 snew(c1,nvec);
138 for (i=0; (i<nvec); i++)
139 c1[i]=NULL;
140 n_alloc=0;
142 natoms=read_first_x(oenv,&status,ftp2fn(efTRX,NFILE,fnm),&t,&x,box);
143 snew(x_s,natoms);
145 /* Start the loop over frames */
146 t1 = t0 = t;
147 teller = 0;
148 do {
149 if (teller >= n_alloc) {
150 n_alloc+=100;
151 for (i=0; (i<nvec); i++)
152 srenew(c1[i],DIM*n_alloc);
154 t1 = t;
156 /* Remove periodicity */
157 rm_pbc(&(top->idef),ePBC,natoms,box,x,x_s);
159 /* Compute crossproducts for all vectors, if triplets.
160 * else, just get the vectors in case of doublets.
162 if (bVec == FALSE) {
163 for (i=0; (i<nvec); i++) {
164 ai=index[3*i];
165 aj=index[3*i+1];
166 ak=index[3*i+2];
167 rvec_sub(x_s[ai],x_s[aj],xij);
168 rvec_sub(x_s[aj],x_s[ak],xjk);
169 cprod(xij,xjk,n);
170 for(m=0; (m<DIM); m++)
171 c1[i][DIM*teller+m]=n[m];
174 else {
175 for (i=0; (i<nvec); i++) {
176 ai=index[2*i];
177 aj=index[2*i+1];
178 rvec_sub(x_s[ai],x_s[aj],n);
179 for(m=0; (m<DIM); m++)
180 c1[i][DIM*teller+m]=n[m];
183 /* Increment loop counter */
184 teller++;
185 } while (read_next_x(oenv,status,&t,natoms,x,box));
186 close_trj(status);
187 fprintf(stderr,"\nDone with trajectory\n");
189 /* Autocorrelation function */
190 if (teller < 2)
191 fprintf(stderr,"Not enough frames for correlation function\n");
192 else {
193 dt=(t1 - t0)/(teller-1);
195 mode = eacVector;
197 do_autocorr(ftp2fn(efXVG,NFILE,fnm),oenv,"Rotational Correlation Function",
198 teller,nvec,c1,dt,mode,bAver);
201 do_view(oenv,ftp2fn(efXVG,NFILE,fnm),NULL);
203 thanx(stderr);
205 return 0;