Remove gmx custom fixed int (e.g. gmx_int64_t) types
[gromacs.git] / src / gromacs / gmxana / gmx_rotacf.cpp
blob771efcaf656333dc1ce68581ce12d9bc8d1473ff
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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, 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.
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.
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.
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.
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.
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.
37 #include "gmxpre.h"
39 #include <cmath>
40 #include <cstring>
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"
56 int gmx_rotacf(int argc, char *argv[])
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."
77 static gmx_bool bVec = FALSE, bAver = TRUE;
79 t_pargs pa[] = {
80 { "-d", FALSE, etBOOL, {&bVec},
81 "Use index doublets (vectors) for correlation function instead of triplets (planes)" },
82 { "-aver", FALSE, etBOOL, {&bAver},
83 "Average over molecules" }
86 t_trxstatus *status;
87 int isize;
88 int *index;
89 char *grpname;
90 rvec *x, *x_s;
91 matrix box;
92 real **c1;
93 rvec xij, xjk, n;
94 int i, m, teller, n_alloc, natoms, nvec, ai, aj, ak;
95 unsigned long mode;
96 real t, t0, t1, dt;
97 gmx_rmpbc_t gpbc = nullptr;
98 t_topology *top;
99 int ePBC;
100 t_filenm fnm[] = {
101 { efTRX, "-f", nullptr, ffREAD },
102 { efTPR, nullptr, nullptr, ffREAD },
103 { efNDX, nullptr, nullptr, ffREAD },
104 { efXVG, "-o", "rotacf", ffWRITE }
106 #define NFILE asize(fnm)
107 int npargs;
108 t_pargs *ppa;
110 gmx_output_env_t *oenv;
112 npargs = asize(pa);
113 ppa = add_acf_pargs(&npargs, pa);
115 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME,
116 NFILE, fnm, npargs, ppa, asize(desc), desc, 0, nullptr, &oenv))
118 sfree(ppa);
119 return 0;
122 rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
124 if (bVec)
126 nvec = isize/2;
128 else
130 nvec = isize/3;
133 if (((isize % 3) != 0) && !bVec)
135 gmx_fatal(FARGS, "number of index elements not multiple of 3, "
136 "these can not be atom triplets\n");
138 if (((isize % 2) != 0) && bVec)
140 gmx_fatal(FARGS, "number of index elements not multiple of 2, "
141 "these can not be atom doublets\n");
144 top = read_top(ftp2fn(efTPR, NFILE, fnm), &ePBC);
146 snew(c1, nvec);
147 for (i = 0; (i < nvec); i++)
149 c1[i] = nullptr;
151 n_alloc = 0;
153 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
154 snew(x_s, natoms);
156 gpbc = gmx_rmpbc_init(&(top->idef), ePBC, natoms);
158 /* Start the loop over frames */
159 t0 = t;
160 teller = 0;
163 if (teller >= n_alloc)
165 n_alloc += 100;
166 for (i = 0; (i < nvec); i++)
168 srenew(c1[i], DIM*n_alloc);
171 t1 = t;
173 /* Remove periodicity */
174 gmx_rmpbc_copy(gpbc, natoms, box, x, x_s);
176 /* Compute crossproducts for all vectors, if triplets.
177 * else, just get the vectors in case of doublets.
179 if (bVec == FALSE)
181 for (i = 0; (i < nvec); i++)
183 ai = index[3*i];
184 aj = index[3*i+1];
185 ak = index[3*i+2];
186 rvec_sub(x_s[ai], x_s[aj], xij);
187 rvec_sub(x_s[aj], x_s[ak], xjk);
188 cprod(xij, xjk, n);
189 for (m = 0; (m < DIM); m++)
191 c1[i][DIM*teller+m] = n[m];
195 else
197 for (i = 0; (i < nvec); i++)
199 ai = index[2*i];
200 aj = index[2*i+1];
201 rvec_sub(x_s[ai], x_s[aj], n);
202 for (m = 0; (m < DIM); m++)
204 c1[i][DIM*teller+m] = n[m];
208 /* Increment loop counter */
209 teller++;
211 while (read_next_x(oenv, status, &t, x, box));
212 close_trx(status);
213 fprintf(stderr, "\nDone with trajectory\n");
215 gmx_rmpbc_done(gpbc);
218 /* Autocorrelation function */
219 if (teller < 2)
221 fprintf(stderr, "Not enough frames for correlation function\n");
223 else
225 dt = (t1 - t0)/(teller-1);
227 mode = eacVector;
229 do_autocorr(ftp2fn(efXVG, NFILE, fnm), oenv, "Rotational Correlation Function",
230 teller, nvec, c1, dt, mode, bAver);
233 do_view(oenv, ftp2fn(efXVG, NFILE, fnm), nullptr);
235 return 0;