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, 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.
43 #include "gromacs/fileio/confio.h"
44 #include "gmx_fatal.h"
45 #include "gromacs/fileio/futil.h"
48 #include "gromacs/math/utilities.h"
51 #include "gromacs/utility/smalloc.h"
52 #include "gromacs/commandline/pargs.h"
60 #include "gromacs/fft/fft.h"
61 #include "gromacs/fileio/trxio.h"
63 static void index_atom2mol(int *n, atom_id *index, t_block *mols)
65 int nat, i, nmol, mol, j;
73 while (index[i] > mols->index[mol])
78 gmx_fatal(FARGS, "Atom index out of range: %d", index[i]+1);
81 for (j = mols->index[mol]; j < mols->index[mol+1]; j++)
83 if (i >= nat || index[i] != j)
85 gmx_fatal(FARGS, "The index group does not consist of whole molecules");
92 fprintf(stderr, "\nSplit group of %d atoms into %d molecules\n", nat, nmol);
97 static void precalc(t_topology top, real normm[])
103 for (i = 0; i < top.mols.nr; i++)
105 k = top.mols.index[i];
106 l = top.mols.index[i+1];
109 for (j = k; j < l; j++)
111 mtot += top.atoms.atom[j].m;
114 for (j = k; j < l; j++)
116 normm[j] = top.atoms.atom[j].m/mtot;
123 static void calc_spectrum(int n, real c[], real dt, const char *fn,
124 output_env_t oenv, gmx_bool bRecip)
130 real nu, omega, recip_fac;
133 for (i = 0; (i < n); i++)
138 if ((status = gmx_fft_init_1d_real(&fft, n, GMX_FFT_FLAG_NONE)) != 0)
140 gmx_fatal(FARGS, "Invalid fft return status %d", status);
142 if ((status = gmx_fft_1d_real(fft, GMX_FFT_REAL_TO_COMPLEX, data, data)) != 0)
144 gmx_fatal(FARGS, "Invalid fft return status %d", status);
146 fp = xvgropen(fn, "Vibrational Power Spectrum",
147 bRecip ? "\\f{12}w\\f{4} (cm\\S-1\\N)" :
148 "\\f{12}n\\f{4} (ps\\S-1\\N)",
150 /* This is difficult.
151 * The length of the ACF is dt (as passed to this routine).
152 * We pass the vacf with N time steps from 0 to dt.
153 * That means that after FFT we have lowest frequency = 1/dt
154 * then 1/(2 dt) etc. (this is the X-axis of the data after FFT).
155 * To convert to 1/cm we need to have to realize that
156 * E = hbar w = h nu = h c/lambda. We want to have reciprokal cm
157 * on the x-axis, that is 1/lambda, so we then have
158 * 1/lambda = nu/c. Since nu has units of 1/ps and c has gromacs units
159 * of nm/ps, we need to multiply by 1e7.
160 * The timestep between saving the trajectory is
161 * 1e7 is to convert nanometer to cm
163 recip_fac = bRecip ? (1e7/SPEED_OF_LIGHT) : 1.0;
164 for (i = 0; (i < n); i += 2)
167 omega = nu*recip_fac;
168 /* Computing the square magnitude of a complex number, since this is a power
171 fprintf(fp, "%10g %10g\n", omega, sqr(data[i])+sqr(data[i+1]));
174 gmx_fft_destroy(fft);
178 int gmx_velacc(int argc, char *argv[])
180 const char *desc[] = {
181 "[THISMODULE] computes the velocity autocorrelation function.",
182 "When the [TT]-m[tt] option is used, the momentum autocorrelation",
183 "function is calculated.[PAR]",
184 "With option [TT]-mol[tt] the velocity autocorrelation function of",
185 "molecules is calculated. In this case the index group should consist",
186 "of molecule numbers instead of atom numbers.[PAR]",
187 "Be sure that your trajectory contains frames with velocity information",
188 "(i.e. [TT]nstvout[tt] was set in your original [TT].mdp[tt] file),",
189 "and that the time interval between data collection points is",
190 "much shorter than the time scale of the autocorrelation."
193 static gmx_bool bMass = FALSE, bMol = FALSE, bRecip = TRUE;
195 { "-m", FALSE, etBOOL, {&bMass},
196 "Calculate the momentum autocorrelation function" },
197 { "-recip", FALSE, etBOOL, {&bRecip},
198 "Use cm^-1 on X-axis instead of 1/ps for spectra." },
199 { "-mol", FALSE, etBOOL, {&bMol},
200 "Calculate the velocity acf of molecules" }
207 gmx_bool bTPS = FALSE, bTop = FALSE;
212 /* t0, t1 are the beginning and end time respectively.
213 * dt is the time step, mass is temp variable for atomic mass.
215 real t0, t1, dt, mass;
217 int counter, n_alloc, i, j, counter_dim, k, l;
219 /* Array for the correlation function */
227 { efTRN, "-f", NULL, ffREAD },
228 { efTPS, NULL, NULL, ffOPTRD },
229 { efNDX, NULL, NULL, ffOPTRD },
230 { efXVG, "-o", "vac", ffWRITE },
231 { efXVG, "-os", "spectrum", ffOPTWR }
233 #define NFILE asize(fnm)
238 ppa = add_acf_pargs(&npargs, pa);
239 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
240 NFILE, fnm, npargs, ppa, asize(desc), desc, 0, NULL, &oenv))
247 bTPS = ftp2bSet(efTPS, NFILE, fnm) || !ftp2bSet(efNDX, NFILE, fnm);
252 bTop = read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, NULL, NULL, box,
254 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
258 rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
265 gmx_fatal(FARGS, "Need a topology to determine the molecules");
267 snew(normm, top.atoms.nr);
269 index_atom2mol(&gnx, index, &top.mols);
272 /* Correlation stuff */
274 for (i = 0; (i < gnx); i++)
279 read_first_frame(oenv, &status, ftp2fn(efTRN, NFILE, fnm), &fr, TRX_NEED_V);
286 if (counter >= n_alloc)
289 for (i = 0; i < gnx; i++)
291 srenew(c1[i], DIM*n_alloc);
294 counter_dim = DIM*counter;
297 for (i = 0; i < gnx; i++)
300 k = top.mols.index[index[i]];
301 l = top.mols.index[index[i]+1];
302 for (j = k; j < l; j++)
306 mass = top.atoms.atom[j].m;
312 mv_mol[XX] += mass*fr.v[j][XX];
313 mv_mol[YY] += mass*fr.v[j][YY];
314 mv_mol[ZZ] += mass*fr.v[j][ZZ];
316 c1[i][counter_dim+XX] = mv_mol[XX];
317 c1[i][counter_dim+YY] = mv_mol[YY];
318 c1[i][counter_dim+ZZ] = mv_mol[ZZ];
323 for (i = 0; i < gnx; i++)
327 mass = top.atoms.atom[index[i]].m;
333 c1[i][counter_dim+XX] = mass*fr.v[index[i]][XX];
334 c1[i][counter_dim+YY] = mass*fr.v[index[i]][YY];
335 c1[i][counter_dim+ZZ] = mass*fr.v[index[i]][ZZ];
343 while (read_next_frame(oenv, status, &fr));
349 /* Compute time step between frames */
350 dt = (t1-t0)/(counter-1);
351 do_autocorr(opt2fn("-o", NFILE, fnm), oenv,
353 "Momentum Autocorrelation Function" :
354 "Velocity Autocorrelation Function",
355 counter, gnx, c1, dt, eacVector, TRUE);
357 do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy");
359 if (opt2bSet("-os", NFILE, fnm))
361 calc_spectrum(counter/2, (real *) (c1[0]), (t1-t0)/2, opt2fn("-os", NFILE, fnm),
363 do_view(oenv, opt2fn("-os", NFILE, fnm), "-nxy");
368 fprintf(stderr, "Not enough frames in trajectory - no output generated.\n");