Upped the version to 3.2.0
[gromacs.git] / src / tools / g_rotacf.c
blobc46f94d99dd803009b84c58f6eec5895b15ca049
1 /*
2 * $Id$
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.2.0
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * Green Red Orange Magenta Azure Cyan Skyblue
36 #include <math.h>
37 #include <string.h>
38 #include "sysstuff.h"
39 #include "physics.h"
40 #include "typedefs.h"
41 #include "smalloc.h"
42 #include "futil.h"
43 #include "statutil.h"
44 #include "copyrite.h"
45 #include "index.h"
46 #include "macros.h"
47 #include "fatal.h"
48 #include "xvgr.h"
49 #include "gstat.h"
50 #include "vec.h"
52 int gmx_rotacf(int argc,char *argv[])
54 static char *desc[] = {
55 "g_rotacf calculates the rotational correlation function",
56 "for molecules. Three atoms (i,j,k) must be given in the index",
57 "file, defining two vectors ij and jk. The rotational acf",
58 "is calculated as the autocorrelation function of the vector",
59 "n = ij x jk, i.e. the cross product of the two vectors.",
60 "Since three atoms span a plane, the order of the three atoms",
61 "does not matter. Optionally, controlled by the -d switch, you can",
62 "calculate the rotational correlation function for linear molecules",
63 "by specifying two atoms (i,j) in the index file.",
64 "[PAR]",
65 "EXAMPLES[PAR]",
66 "g_rotacf -P 1 -nparm 2 -fft -n index -o rotacf-x-P1",
67 "-fa expfit-x-P1 -beginfit 2.5 -endfit 20.0[PAR]",
68 "This will calculate the rotational correlation function using a first",
69 "order Legendre polynomial of the angle of a vector defined by the index",
70 "file. The correlation function will be fitted from 2.5 ps till 20.0 ps",
71 "to a two parameter exponential",
76 static bool bVec = FALSE,bAver=TRUE;
78 t_pargs pa[] = {
79 { "-d", FALSE, etBOOL, {&bVec},
80 "Use index doublets (vectors) for correlation function instead of triplets (planes)" },
81 { "-aver",FALSE, etBOOL, {&bAver},
82 "Average over molecules" }
85 int status,isize;
86 atom_id *index;
87 char *grpname;
88 rvec *x,*x_s;
89 matrix box;
90 real **c1;
91 rvec xij,xjk,n;
92 int i,m,teller,n_alloc,natoms,nvec,ai,aj,ak;
93 unsigned long mode;
94 real t,t0,t1,dt;
95 t_topology *top;
96 t_filenm fnm[] = {
97 { efTRX, "-f", NULL, ffREAD },
98 { efTPX, NULL, NULL, ffREAD },
99 { efNDX, NULL, NULL, ffREAD },
100 { efXVG, "-o", "rotacf", ffWRITE }
102 #define NFILE asize(fnm)
103 int npargs;
104 t_pargs *ppa;
106 CopyRight(stderr,argv[0]);
107 npargs = asize(pa);
108 ppa = add_acf_pargs(&npargs,pa);
110 parse_common_args(&argc,argv,PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
111 NFILE,fnm,npargs,ppa,asize(desc),desc,0,NULL);
113 rd_index(ftp2fn(efNDX,NFILE,fnm),1,&isize,&index,&grpname);
115 if (bVec)
116 nvec = isize/2;
117 else
118 nvec = isize/3;
120 if (((isize % 3) != 0) && !bVec)
121 fatal_error(0,"number of index elements not multiple of 3, "
122 "these can not be atom triplets\n");
123 if (((isize % 2) != 0) && bVec)
124 fatal_error(0,"number of index elements not multiple of 2, "
125 "these can not be atom doublets\n");
127 top=read_top(ftp2fn(efTPX,NFILE,fnm));
129 snew(c1,nvec);
130 for (i=0; (i<nvec); i++)
131 c1[i]=NULL;
132 n_alloc=0;
134 natoms=read_first_x(&status,ftp2fn(efTRX,NFILE,fnm),&t,&x,box);
135 snew(x_s,natoms);
137 /* Start the loop over frames */
138 t1 = t0 = t;
139 teller = 0;
140 do {
141 if (teller >= n_alloc) {
142 n_alloc+=100;
143 for (i=0; (i<nvec); i++)
144 srenew(c1[i],DIM*n_alloc);
146 t1 = t;
148 /* Remove periodicity */
149 rm_pbc(&(top->idef),natoms,box,x,x_s);
151 /* Compute crossproducts for all vectors, if triplets.
152 * else, just get the vectors in case of doublets.
154 if (bVec == FALSE) {
155 for (i=0; (i<nvec); i++) {
156 ai=index[3*i];
157 aj=index[3*i+1];
158 ak=index[3*i+2];
159 rvec_sub(x_s[ai],x_s[aj],xij);
160 rvec_sub(x_s[aj],x_s[ak],xjk);
161 oprod(xij,xjk,n);
162 for(m=0; (m<DIM); m++)
163 c1[i][DIM*teller+m]=n[m];
166 else {
167 for (i=0; (i<nvec); i++) {
168 ai=index[2*i];
169 aj=index[2*i+1];
170 rvec_sub(x_s[ai],x_s[aj],n);
171 for(m=0; (m<DIM); m++)
172 c1[i][DIM*teller+m]=n[m];
175 /* Increment loop counter */
176 teller++;
177 } while (read_next_x(status,&t,natoms,x,box));
178 close_trj(status);
179 fprintf(stderr,"\nDone with trajectory\n");
181 /* Autocorrelation function */
182 if (teller < 2)
183 fprintf(stderr,"Not enough frames for correlation function\n");
184 else {
185 dt=(t1 - t0)/(teller-1);
187 mode = eacVector;
189 do_autocorr(ftp2fn(efXVG,NFILE,fnm),"Rotational Correlation Function",
190 teller,nvec,c1,dt,mode,bAver);
193 do_view(ftp2fn(efXVG,NFILE,fnm),NULL);
195 thanx(stderr);
197 return 0;