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
41 #include "gromacs/fileio/confio.h"
42 #include "gmx_fatal.h"
43 #include "gromacs/fileio/futil.h"
50 #include "gromacs/commandline/pargs.h"
59 #include "gromacs/fft/fft.h"
60 #include "gromacs/fileio/trxio.h"
62 static void index_atom2mol(int *n, atom_id *index, t_block *mols)
64 int nat, i, nmol, mol, j;
72 while (index[i] > mols->index[mol])
77 gmx_fatal(FARGS, "Atom index out of range: %d", index[i]+1);
80 for (j = mols->index[mol]; j < mols->index[mol+1]; j++)
82 if (i >= nat || index[i] != j)
84 gmx_fatal(FARGS, "The index group does not consist of whole molecules");
91 fprintf(stderr, "\nSplit group of %d atoms into %d molecules\n", nat, nmol);
96 static void precalc(t_topology top, real normm[])
102 for (i = 0; i < top.mols.nr; i++)
104 k = top.mols.index[i];
105 l = top.mols.index[i+1];
108 for (j = k; j < l; j++)
110 mtot += top.atoms.atom[j].m;
113 for (j = k; j < l; j++)
115 normm[j] = top.atoms.atom[j].m/mtot;
122 static void calc_spectrum(int n, real c[], real dt, const char *fn,
123 output_env_t oenv, gmx_bool bRecip)
129 real nu, omega, recip_fac;
132 for (i = 0; (i < n); i++)
137 if ((status = gmx_fft_init_1d_real(&fft, n, GMX_FFT_FLAG_NONE)) != 0)
139 gmx_fatal(FARGS, "Invalid fft return status %d", status);
141 if ((status = gmx_fft_1d_real(fft, GMX_FFT_REAL_TO_COMPLEX, data, data)) != 0)
143 gmx_fatal(FARGS, "Invalid fft return status %d", status);
145 fp = xvgropen(fn, "Vibrational Power Spectrum",
146 bRecip ? "\\f{12}w\\f{4} (cm\\S-1\\N)" :
147 "\\f{12}n\\f{4} (ps\\S-1\\N)",
149 /* This is difficult.
150 * The length of the ACF is dt (as passed to this routine).
151 * We pass the vacf with N time steps from 0 to dt.
152 * That means that after FFT we have lowest frequency = 1/dt
153 * then 1/(2 dt) etc. (this is the X-axis of the data after FFT).
154 * To convert to 1/cm we need to have to realize that
155 * E = hbar w = h nu = h c/lambda. We want to have reciprokal cm
156 * on the x-axis, that is 1/lambda, so we then have
157 * 1/lambda = nu/c. Since nu has units of 1/ps and c has gromacs units
158 * of nm/ps, we need to multiply by 1e7.
159 * The timestep between saving the trajectory is
160 * 1e7 is to convert nanometer to cm
162 recip_fac = bRecip ? (1e7/SPEED_OF_LIGHT) : 1.0;
163 for (i = 0; (i < n); i += 2)
166 omega = nu*recip_fac;
167 /* Computing the square magnitude of a complex number, since this is a power
170 fprintf(fp, "%10g %10g\n", omega, sqr(data[i])+sqr(data[i+1]));
173 gmx_fft_destroy(fft);
177 int gmx_velacc(int argc, char *argv[])
179 const char *desc[] = {
180 "[THISMODULE] computes the velocity autocorrelation function.",
181 "When the [TT]-m[tt] option is used, the momentum autocorrelation",
182 "function is calculated.[PAR]",
183 "With option [TT]-mol[tt] the velocity autocorrelation function of",
184 "molecules is calculated. In this case the index group should consist",
185 "of molecule numbers instead of atom numbers.[PAR]",
186 "Be sure that your trajectory contains frames with velocity information",
187 "(i.e. [TT]nstvout[tt] was set in your original [TT].mdp[tt] file),",
188 "and that the time interval between data collection points is",
189 "much shorter than the time scale of the autocorrelation."
192 static gmx_bool bMass = FALSE, bMol = FALSE, bRecip = TRUE;
194 { "-m", FALSE, etBOOL, {&bMass},
195 "Calculate the momentum autocorrelation function" },
196 { "-recip", FALSE, etBOOL, {&bRecip},
197 "Use cm^-1 on X-axis instead of 1/ps for spectra." },
198 { "-mol", FALSE, etBOOL, {&bMol},
199 "Calculate the velocity acf of molecules" }
206 gmx_bool bTPS = FALSE, bTop = FALSE;
211 /* t0, t1 are the beginning and end time respectively.
212 * dt is the time step, mass is temp variable for atomic mass.
214 real t0, t1, dt, mass;
216 int counter, n_alloc, i, j, counter_dim, k, l;
218 /* Array for the correlation function */
226 { efTRN, "-f", NULL, ffREAD },
227 { efTPS, NULL, NULL, ffOPTRD },
228 { efNDX, NULL, NULL, ffOPTRD },
229 { efXVG, "-o", "vac", ffWRITE },
230 { efXVG, "-os", "spectrum", ffOPTWR }
232 #define NFILE asize(fnm)
237 ppa = add_acf_pargs(&npargs, pa);
238 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
239 NFILE, fnm, npargs, ppa, asize(desc), desc, 0, NULL, &oenv))
246 bTPS = ftp2bSet(efTPS, NFILE, fnm) || !ftp2bSet(efNDX, NFILE, fnm);
251 bTop = read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, NULL, NULL, box,
253 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
257 rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
264 gmx_fatal(FARGS, "Need a topology to determine the molecules");
266 snew(normm, top.atoms.nr);
268 index_atom2mol(&gnx, index, &top.mols);
271 /* Correlation stuff */
273 for (i = 0; (i < gnx); i++)
278 read_first_frame(oenv, &status, ftp2fn(efTRN, NFILE, fnm), &fr, TRX_NEED_V);
285 if (counter >= n_alloc)
288 for (i = 0; i < gnx; i++)
290 srenew(c1[i], DIM*n_alloc);
293 counter_dim = DIM*counter;
296 for (i = 0; i < gnx; i++)
299 k = top.mols.index[index[i]];
300 l = top.mols.index[index[i]+1];
301 for (j = k; j < l; j++)
305 mass = top.atoms.atom[j].m;
311 mv_mol[XX] += mass*fr.v[j][XX];
312 mv_mol[YY] += mass*fr.v[j][YY];
313 mv_mol[ZZ] += mass*fr.v[j][ZZ];
315 c1[i][counter_dim+XX] = mv_mol[XX];
316 c1[i][counter_dim+YY] = mv_mol[YY];
317 c1[i][counter_dim+ZZ] = mv_mol[ZZ];
322 for (i = 0; i < gnx; i++)
326 mass = top.atoms.atom[index[i]].m;
332 c1[i][counter_dim+XX] = mass*fr.v[index[i]][XX];
333 c1[i][counter_dim+YY] = mass*fr.v[index[i]][YY];
334 c1[i][counter_dim+ZZ] = mass*fr.v[index[i]][ZZ];
342 while (read_next_frame(oenv, status, &fr));
348 /* Compute time step between frames */
349 dt = (t1-t0)/(counter-1);
350 do_autocorr(opt2fn("-o", NFILE, fnm), oenv,
352 "Momentum Autocorrelation Function" :
353 "Velocity Autocorrelation Function",
354 counter, gnx, c1, dt, eacVector, TRUE);
356 do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy");
358 if (opt2bSet("-os", NFILE, fnm))
360 calc_spectrum(counter/2, (real *) (c1[0]), (t1-t0)/2, opt2fn("-os", NFILE, fnm),
362 do_view(oenv, opt2fn("-os", NFILE, fnm), "-nxy");
367 fprintf(stderr, "Not enough frames in trajectory - no output generated.\n");