3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
50 #include "gmx_fatal.h"
57 int gmx_rotacf(int argc, char *argv[])
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.",
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."
78 static gmx_bool bVec = FALSE, bAver = TRUE;
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" }
95 int i, m, teller, n_alloc, natoms, nvec, ai, aj, ak;
98 gmx_rmpbc_t gpbc = NULL;
102 { efTRX, "-f", NULL, ffREAD },
103 { efTPX, NULL, NULL, ffREAD },
104 { efNDX, NULL, NULL, ffREAD },
105 { efXVG, "-o", "rotacf", ffWRITE }
107 #define NFILE asize(fnm)
114 ppa = add_acf_pargs(&npargs, pa);
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);
119 rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
130 if (((isize % 3) != 0) && !bVec)
132 gmx_fatal(FARGS, "number of index elements not multiple of 3, "
133 "these can not be atom triplets\n");
135 if (((isize % 2) != 0) && bVec)
137 gmx_fatal(FARGS, "number of index elements not multiple of 2, "
138 "these can not be atom doublets\n");
141 top = read_top(ftp2fn(efTPX, NFILE, fnm), &ePBC);
144 for (i = 0; (i < nvec); i++)
150 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
153 gpbc = gmx_rmpbc_init(&(top->idef), ePBC, natoms, box);
155 /* Start the loop over frames */
160 if (teller >= n_alloc)
163 for (i = 0; (i < nvec); i++)
165 srenew(c1[i], DIM*n_alloc);
170 /* Remove periodicity */
171 gmx_rmpbc_copy(gpbc, natoms, box, x, x_s);
173 /* Compute crossproducts for all vectors, if triplets.
174 * else, just get the vectors in case of doublets.
178 for (i = 0; (i < nvec); i++)
183 rvec_sub(x_s[ai], x_s[aj], xij);
184 rvec_sub(x_s[aj], x_s[ak], xjk);
186 for (m = 0; (m < DIM); m++)
188 c1[i][DIM*teller+m] = n[m];
194 for (i = 0; (i < nvec); i++)
198 rvec_sub(x_s[ai], x_s[aj], n);
199 for (m = 0; (m < DIM); m++)
201 c1[i][DIM*teller+m] = n[m];
205 /* Increment loop counter */
208 while (read_next_x(oenv, status, &t, natoms, x, box));
210 fprintf(stderr, "\nDone with trajectory\n");
212 gmx_rmpbc_done(gpbc);
215 /* Autocorrelation function */
218 fprintf(stderr, "Not enough frames for correlation function\n");
222 dt = (t1 - t0)/(teller-1);
226 do_autocorr(ftp2fn(efXVG, NFILE, fnm), oenv, "Rotational Correlation Function",
227 teller, nvec, c1, dt, mode, bAver);
230 do_view(oenv, ftp2fn(efXVG, NFILE, fnm), NULL);