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
49 #include "gmx_fatal.h"
56 int gmx_rotacf(int argc, char *argv[])
58 const char *desc[] = {
59 "[TT]g_rotacf[tt] 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.",
70 "[TT]g_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;
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" }
94 int i, m, teller, n_alloc, natoms, nvec, ai, aj, ak;
97 gmx_rmpbc_t gpbc = NULL;
101 { efTRX, "-f", NULL, ffREAD },
102 { efTPX, NULL, NULL, ffREAD },
103 { efNDX, NULL, NULL, ffREAD },
104 { efXVG, "-o", "rotacf", ffWRITE }
106 #define NFILE asize(fnm)
113 ppa = add_acf_pargs(&npargs, pa);
115 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
116 NFILE, fnm, npargs, ppa, asize(desc), desc, 0, NULL, &oenv))
121 rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
132 if (((isize % 3) != 0) && !bVec)
134 gmx_fatal(FARGS, "number of index elements not multiple of 3, "
135 "these can not be atom triplets\n");
137 if (((isize % 2) != 0) && bVec)
139 gmx_fatal(FARGS, "number of index elements not multiple of 2, "
140 "these can not be atom doublets\n");
143 top = read_top(ftp2fn(efTPX, NFILE, fnm), &ePBC);
146 for (i = 0; (i < nvec); i++)
152 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
155 gpbc = gmx_rmpbc_init(&(top->idef), ePBC, natoms);
157 /* Start the loop over frames */
162 if (teller >= n_alloc)
165 for (i = 0; (i < nvec); i++)
167 srenew(c1[i], DIM*n_alloc);
172 /* Remove periodicity */
173 gmx_rmpbc_copy(gpbc, natoms, box, x, x_s);
175 /* Compute crossproducts for all vectors, if triplets.
176 * else, just get the vectors in case of doublets.
180 for (i = 0; (i < nvec); i++)
185 rvec_sub(x_s[ai], x_s[aj], xij);
186 rvec_sub(x_s[aj], x_s[ak], xjk);
188 for (m = 0; (m < DIM); m++)
190 c1[i][DIM*teller+m] = n[m];
196 for (i = 0; (i < nvec); i++)
200 rvec_sub(x_s[ai], x_s[aj], n);
201 for (m = 0; (m < DIM); m++)
203 c1[i][DIM*teller+m] = n[m];
207 /* Increment loop counter */
210 while (read_next_x(oenv, status, &t, x, box));
212 fprintf(stderr, "\nDone with trajectory\n");
214 gmx_rmpbc_done(gpbc);
217 /* Autocorrelation function */
220 fprintf(stderr, "Not enough frames for correlation function\n");
224 dt = (t1 - t0)/(teller-1);
228 do_autocorr(ftp2fn(efXVG, NFILE, fnm), oenv, "Rotational Correlation Function",
229 teller, nvec, c1, dt, mode, bAver);
232 do_view(oenv, ftp2fn(efXVG, NFILE, fnm), NULL);