Move physics.* to math/units.*
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_rotacf.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <math.h>
42 #include <string.h>
43
44 #include "typedefs.h"
45 #include "gromacs/utility/smalloc.h"
46 #include "gromacs/utility/futil.h"
47 #include "gromacs/commandline/pargs.h"
48 #include "index.h"
49 #include "macros.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gstat.h"
52 #include "gromacs/math/vec.h"
53 #include "viewit.h"
54 #include "gmx_ana.h"
55 #include "gromacs/fileio/trxio.h"
56
57
58 int gmx_rotacf(int argc, char *argv[])
59 {
60     const char     *desc[] = {
61         "[THISMODULE] calculates the rotational correlation function",
62         "for molecules. Atom triplets (i,j,k) must be given in the index",
63         "file, defining two vectors ij and jk. The rotational ACF",
64         "is calculated as the autocorrelation function of the vector",
65         "n = ij x jk, i.e. the cross product of the two vectors.",
66         "Since three atoms span a plane, the order of the three atoms",
67         "does not matter. Optionally, by invoking the [TT]-d[tt] switch, you can",
68         "calculate the rotational correlation function for linear molecules",
69         "by specifying atom pairs (i,j) in the index file.",
70         "[PAR]",
71         "EXAMPLES[PAR]",
72         "[TT]gmx rotacf -P 1 -nparm 2 -fft -n index -o rotacf-x-P1",
73         "-fa expfit-x-P1 -beginfit 2.5 -endfit 20.0[tt][PAR]",
74         "This will calculate the rotational correlation function using a first",
75         "order Legendre polynomial of the angle of a vector defined by the index",
76         "file. The correlation function will be fitted from 2.5 ps until 20.0 ps",
77         "to a two-parameter exponential."
78     };
79     static gmx_bool bVec    = FALSE, bAver = TRUE;
80
81     t_pargs         pa[] = {
82         { "-d",   FALSE, etBOOL, {&bVec},
83           "Use index doublets (vectors) for correlation function instead of triplets (planes)" },
84         { "-aver", FALSE, etBOOL, {&bAver},
85           "Average over molecules" }
86     };
87
88     t_trxstatus    *status;
89     int             isize;
90     atom_id        *index;
91     char           *grpname;
92     rvec           *x, *x_s;
93     matrix          box;
94     real          **c1;
95     rvec            xij, xjk, n;
96     int             i, m, teller, n_alloc, natoms, nvec, ai, aj, ak;
97     unsigned long   mode;
98     real            t, t0, t1, dt;
99     gmx_rmpbc_t     gpbc = NULL;
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 }
107     };
108 #define NFILE asize(fnm)
109     int             npargs;
110     t_pargs        *ppa;
111
112     output_env_t    oenv;
113
114     npargs = asize(pa);
115     ppa    = add_acf_pargs(&npargs, pa);
116
117     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
118                            NFILE, fnm, npargs, ppa, asize(desc), desc, 0, NULL, &oenv))
119     {
120         return 0;
121     }
122
123     rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
124
125     if (bVec)
126     {
127         nvec = isize/2;
128     }
129     else
130     {
131         nvec = isize/3;
132     }
133
134     if (((isize % 3) != 0) && !bVec)
135     {
136         gmx_fatal(FARGS, "number of index elements not multiple of 3, "
137                   "these can not be atom triplets\n");
138     }
139     if (((isize % 2) != 0) && bVec)
140     {
141         gmx_fatal(FARGS, "number of index elements not multiple of 2, "
142                   "these can not be atom doublets\n");
143     }
144
145     top = read_top(ftp2fn(efTPX, NFILE, fnm), &ePBC);
146
147     snew(c1, nvec);
148     for (i = 0; (i < nvec); i++)
149     {
150         c1[i] = NULL;
151     }
152     n_alloc = 0;
153
154     natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
155     snew(x_s, natoms);
156
157     gpbc = gmx_rmpbc_init(&(top->idef), ePBC, natoms);
158
159     /* Start the loop over frames */
160     t1      = t0 = t;
161     teller  = 0;
162     do
163     {
164         if (teller >= n_alloc)
165         {
166             n_alloc += 100;
167             for (i = 0; (i < nvec); i++)
168             {
169                 srenew(c1[i], DIM*n_alloc);
170             }
171         }
172         t1 = t;
173
174         /* Remove periodicity */
175         gmx_rmpbc_copy(gpbc, natoms, box, x, x_s);
176
177         /* Compute crossproducts for all vectors, if triplets.
178          * else, just get the vectors in case of doublets.
179          */
180         if (bVec == FALSE)
181         {
182             for (i = 0; (i < nvec); i++)
183             {
184                 ai = index[3*i];
185                 aj = index[3*i+1];
186                 ak = index[3*i+2];
187                 rvec_sub(x_s[ai], x_s[aj], xij);
188                 rvec_sub(x_s[aj], x_s[ak], xjk);
189                 cprod(xij, xjk, n);
190                 for (m = 0; (m < DIM); m++)
191                 {
192                     c1[i][DIM*teller+m] = n[m];
193                 }
194             }
195         }
196         else
197         {
198             for (i = 0; (i < nvec); i++)
199             {
200                 ai = index[2*i];
201                 aj = index[2*i+1];
202                 rvec_sub(x_s[ai], x_s[aj], n);
203                 for (m = 0; (m < DIM); m++)
204                 {
205                     c1[i][DIM*teller+m] = n[m];
206                 }
207             }
208         }
209         /* Increment loop counter */
210         teller++;
211     }
212     while (read_next_x(oenv, status, &t, x, box));
213     close_trj(status);
214     fprintf(stderr, "\nDone with trajectory\n");
215
216     gmx_rmpbc_done(gpbc);
217
218
219     /* Autocorrelation function */
220     if (teller < 2)
221     {
222         fprintf(stderr, "Not enough frames for correlation function\n");
223     }
224     else
225     {
226         dt = (t1 - t0)/(teller-1);
227
228         mode = eacVector;
229
230         do_autocorr(ftp2fn(efXVG, NFILE, fnm), oenv, "Rotational Correlation Function",
231                     teller, nvec, c1, dt, mode, bAver);
232     }
233
234     do_view(oenv, ftp2fn(efXVG, NFILE, fnm), NULL);
235
236     return 0;
237 }