Merge "Improve master-specific CMake behavior."
[alexxy/gromacs.git] / src / gromacs / gmxana / 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         "[TT]g_rotacf[tt] calculates the rotational correlation function",
61         "for molecules. Atom triplets (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, by invoking the [TT]-d[tt] switch, you can",
67         "calculate the rotational correlation function for linear molecules",
68         "by specifying atom pairs (i,j) in the index file.",
69         "[PAR]",
70         "EXAMPLES[PAR]",
71         "[TT]g_rotacf -P 1 -nparm 2 -fft -n index -o rotacf-x-P1",
72         "-fa expfit-x-P1 -beginfit 2.5 -endfit 20.0[tt][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 until 20.0 ps",
76         "to a two-parameter exponential."
77     };
78     static gmx_bool bVec    = FALSE, bAver = TRUE;
79
80     t_pargs         pa[] = {
81         { "-d",   FALSE, etBOOL, {&bVec},
82           "Use index doublets (vectors) for correlation function instead of triplets (planes)" },
83         { "-aver", FALSE, etBOOL, {&bAver},
84           "Average over molecules" }
85     };
86
87     t_trxstatus    *status;
88     int             isize;
89     atom_id        *index;
90     char           *grpname;
91     rvec           *x, *x_s;
92     matrix          box;
93     real          **c1;
94     rvec            xij, xjk, n;
95     int             i, m, teller, n_alloc, natoms, nvec, ai, aj, ak;
96     unsigned long   mode;
97     real            t, t0, t1, dt;
98     gmx_rmpbc_t     gpbc = NULL;
99     t_topology     *top;
100     int             ePBC;
101     t_filenm        fnm[] = {
102         { efTRX, "-f", NULL,  ffREAD  },
103         { efTPX, NULL, NULL,  ffREAD },
104         { efNDX, NULL, NULL,  ffREAD  },
105         { efXVG, "-o", "rotacf",  ffWRITE }
106     };
107 #define NFILE asize(fnm)
108     int             npargs;
109     t_pargs        *ppa;
110
111     output_env_t    oenv;
112
113     npargs = asize(pa);
114     ppa    = add_acf_pargs(&npargs, pa);
115
116     parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
117                       NFILE, fnm, npargs, ppa, asize(desc), desc, 0, NULL, &oenv);
118
119     rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
120
121     if (bVec)
122     {
123         nvec = isize/2;
124     }
125     else
126     {
127         nvec = isize/3;
128     }
129
130     if (((isize % 3) != 0) && !bVec)
131     {
132         gmx_fatal(FARGS, "number of index elements not multiple of 3, "
133                   "these can not be atom triplets\n");
134     }
135     if (((isize % 2) != 0) && bVec)
136     {
137         gmx_fatal(FARGS, "number of index elements not multiple of 2, "
138                   "these can not be atom doublets\n");
139     }
140
141     top = read_top(ftp2fn(efTPX, NFILE, fnm), &ePBC);
142
143     snew(c1, nvec);
144     for (i = 0; (i < nvec); i++)
145     {
146         c1[i] = NULL;
147     }
148     n_alloc = 0;
149
150     natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
151     snew(x_s, natoms);
152
153     gpbc = gmx_rmpbc_init(&(top->idef), ePBC, natoms, box);
154
155     /* Start the loop over frames */
156     t1      = t0 = t;
157     teller  = 0;
158     do
159     {
160         if (teller >= n_alloc)
161         {
162             n_alloc += 100;
163             for (i = 0; (i < nvec); i++)
164             {
165                 srenew(c1[i], DIM*n_alloc);
166             }
167         }
168         t1 = t;
169
170         /* Remove periodicity */
171         gmx_rmpbc_copy(gpbc, natoms, box, x, x_s);
172
173         /* Compute crossproducts for all vectors, if triplets.
174          * else, just get the vectors in case of doublets.
175          */
176         if (bVec == FALSE)
177         {
178             for (i = 0; (i < nvec); i++)
179             {
180                 ai = index[3*i];
181                 aj = index[3*i+1];
182                 ak = index[3*i+2];
183                 rvec_sub(x_s[ai], x_s[aj], xij);
184                 rvec_sub(x_s[aj], x_s[ak], xjk);
185                 cprod(xij, xjk, n);
186                 for (m = 0; (m < DIM); m++)
187                 {
188                     c1[i][DIM*teller+m] = n[m];
189                 }
190             }
191         }
192         else
193         {
194             for (i = 0; (i < nvec); i++)
195             {
196                 ai = index[2*i];
197                 aj = index[2*i+1];
198                 rvec_sub(x_s[ai], x_s[aj], n);
199                 for (m = 0; (m < DIM); m++)
200                 {
201                     c1[i][DIM*teller+m] = n[m];
202                 }
203             }
204         }
205         /* Increment loop counter */
206         teller++;
207     }
208     while (read_next_x(oenv, status, &t, natoms, x, box));
209     close_trj(status);
210     fprintf(stderr, "\nDone with trajectory\n");
211
212     gmx_rmpbc_done(gpbc);
213
214
215     /* Autocorrelation function */
216     if (teller < 2)
217     {
218         fprintf(stderr, "Not enough frames for correlation function\n");
219     }
220     else
221     {
222         dt = (t1 - t0)/(teller-1);
223
224         mode = eacVector;
225
226         do_autocorr(ftp2fn(efXVG, NFILE, fnm), oenv, "Rotational Correlation Function",
227                     teller, nvec, c1, dt, mode, bAver);
228     }
229
230     do_view(oenv, ftp2fn(efXVG, NFILE, fnm), NULL);
231
232     thanx(stderr);
233
234     return 0;
235 }