63d21d22c70c4ee29f099068143d8e36ccb3cdd6
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_rotacf.cpp
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,2015,2016,2017,2018,2019, 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 #include "gmxpre.h"
38
39 #include <cmath>
40 #include <cstring>
41
42 #include "gromacs/commandline/pargs.h"
43 #include "gromacs/commandline/viewit.h"
44 #include "gromacs/correlationfunctions/autocorr.h"
45 #include "gromacs/fileio/trxio.h"
46 #include "gromacs/gmxana/gmx_ana.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/pbcutil/rmpbc.h"
49 #include "gromacs/topology/index.h"
50 #include "gromacs/topology/topology.h"
51 #include "gromacs/utility/arraysize.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/futil.h"
54 #include "gromacs/utility/smalloc.h"
55
56 int gmx_rotacf(int argc, char* argv[])
57 {
58     const char* desc[] = {
59         "[THISMODULE] calculates the rotational correlation function",
60         "for molecules. Atom triplets (i,j,k) must be given in the index",
61         "file, defining two vectors ij and jk. The rotational ACF",
62         "is calculated as the autocorrelation function of the vector",
63         "n = ij x jk, i.e. the cross product of the two vectors.",
64         "Since three atoms span a plane, the order of the three atoms",
65         "does not matter. Optionally, by invoking the [TT]-d[tt] switch, you can",
66         "calculate the rotational correlation function for linear molecules",
67         "by specifying atom pairs (i,j) in the index file.",
68         "[PAR]",
69         "EXAMPLES[PAR]",
70         "[TT]gmx rotacf -P 1 -nparm 2 -fft -n index -o rotacf-x-P1",
71         "-fa expfit-x-P1 -beginfit 2.5 -endfit 20.0[tt][PAR]",
72         "This will calculate the rotational correlation function using a first",
73         "order Legendre polynomial of the angle of a vector defined by the index",
74         "file. The correlation function will be fitted from 2.5 ps until 20.0 ps",
75         "to a two-parameter exponential."
76     };
77     static gmx_bool bVec = FALSE, bAver = TRUE;
78
79     t_pargs pa[] = {
80         { "-d",
81           FALSE,
82           etBOOL,
83           { &bVec },
84           "Use index doublets (vectors) for correlation function instead of triplets (planes)" },
85         { "-aver", FALSE, etBOOL, { &bAver }, "Average over molecules" }
86     };
87
88     t_trxstatus*  status;
89     int           isize;
90     int*          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 = nullptr;
100     t_topology*   top;
101     int           ePBC;
102     t_filenm      fnm[] = { { efTRX, "-f", nullptr, ffREAD },
103                        { efTPR, nullptr, nullptr, ffREAD },
104                        { efNDX, nullptr, nullptr, ffREAD },
105                        { efXVG, "-o", "rotacf", ffWRITE } };
106 #define NFILE asize(fnm)
107     int      npargs;
108     t_pargs* ppa;
109
110     gmx_output_env_t* oenv;
111
112     npargs = asize(pa);
113     ppa    = add_acf_pargs(&npargs, pa);
114
115     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME, NFILE, fnm, npargs, ppa,
116                            asize(desc), desc, 0, nullptr, &oenv))
117     {
118         sfree(ppa);
119         return 0;
120     }
121
122     rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
123
124     if (bVec)
125     {
126         nvec = isize / 2;
127     }
128     else
129     {
130         nvec = isize / 3;
131     }
132
133     if (((isize % 3) != 0) && !bVec)
134     {
135         gmx_fatal(FARGS,
136                   "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,
142                   "number of index elements not multiple of 2, "
143                   "these can not be atom doublets\n");
144     }
145
146     top = read_top(ftp2fn(efTPR, NFILE, fnm), &ePBC);
147
148     snew(c1, nvec);
149     for (i = 0; (i < nvec); i++)
150     {
151         c1[i] = nullptr;
152     }
153     n_alloc = 0;
154
155     natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
156     snew(x_s, natoms);
157
158     gpbc = gmx_rmpbc_init(&(top->idef), ePBC, natoms);
159
160     /* Start the loop over frames */
161     t0     = t;
162     teller = 0;
163     do
164     {
165         if (teller >= n_alloc)
166         {
167             n_alloc += 100;
168             for (i = 0; (i < nvec); i++)
169             {
170                 srenew(c1[i], DIM * n_alloc);
171             }
172         }
173         t1 = t;
174
175         /* Remove periodicity */
176         gmx_rmpbc_copy(gpbc, natoms, box, x, x_s);
177
178         /* Compute crossproducts for all vectors, if triplets.
179          * else, just get the vectors in case of doublets.
180          */
181         if (!bVec)
182         {
183             for (i = 0; (i < nvec); i++)
184             {
185                 ai = index[3 * i];
186                 aj = index[3 * i + 1];
187                 ak = index[3 * i + 2];
188                 rvec_sub(x_s[ai], x_s[aj], xij);
189                 rvec_sub(x_s[aj], x_s[ak], xjk);
190                 cprod(xij, xjk, n);
191                 for (m = 0; (m < DIM); m++)
192                 {
193                     c1[i][DIM * teller + m] = n[m];
194                 }
195             }
196         }
197         else
198         {
199             for (i = 0; (i < nvec); i++)
200             {
201                 ai = index[2 * i];
202                 aj = index[2 * i + 1];
203                 rvec_sub(x_s[ai], x_s[aj], n);
204                 for (m = 0; (m < DIM); m++)
205                 {
206                     c1[i][DIM * teller + m] = n[m];
207                 }
208             }
209         }
210         /* Increment loop counter */
211         teller++;
212     } while (read_next_x(oenv, status, &t, x, box));
213     close_trx(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) / (static_cast<real>(teller - 1));
227
228         mode = eacVector;
229
230         do_autocorr(ftp2fn(efXVG, NFILE, fnm), oenv, "Rotational Correlation Function", teller,
231                     nvec, c1, dt, mode, bAver);
232     }
233
234     do_view(oenv, ftp2fn(efXVG, NFILE, fnm), nullptr);
235
236     return 0;
237 }