Cleaned up memory usage in gmx traj and trjconv
[gromacs.git] / src / gromacs / gmxana / gmx_velacc.cpp
blob902b9617cf1dfa504045e67f8f08bcb8b3c8b550
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 <cstdio>
41 #include <cstring>
43 #include "gromacs/commandline/pargs.h"
44 #include "gromacs/commandline/viewit.h"
45 #include "gromacs/correlationfunctions/autocorr.h"
46 #include "gromacs/fft/fft.h"
47 #include "gromacs/fileio/confio.h"
48 #include "gromacs/fileio/trxio.h"
49 #include "gromacs/fileio/xvgr.h"
50 #include "gromacs/gmxana/gmx_ana.h"
51 #include "gromacs/gmxana/gstat.h"
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/units.h"
54 #include "gromacs/math/utilities.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/topology/index.h"
57 #include "gromacs/topology/topology.h"
58 #include "gromacs/trajectory/trajectoryframe.h"
59 #include "gromacs/utility/arraysize.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/futil.h"
62 #include "gromacs/utility/smalloc.h"
64 static void index_atom2mol(int *n, int *index, const t_block *mols)
66 int nat, i, nmol, mol, j;
68 nat = *n;
69 i = 0;
70 nmol = 0;
71 mol = 0;
72 while (i < nat)
74 while (index[i] > mols->index[mol])
76 mol++;
77 if (mol >= mols->nr)
79 gmx_fatal(FARGS, "Atom index out of range: %d", index[i]+1);
82 for (j = mols->index[mol]; j < mols->index[mol+1]; j++)
84 if (i >= nat || index[i] != j)
86 gmx_fatal(FARGS, "The index group does not consist of whole molecules");
88 i++;
90 index[nmol++] = mol;
93 fprintf(stderr, "\nSplit group of %d atoms into %d molecules\n", nat, nmol);
95 *n = nmol;
98 static void precalc(const t_topology &top, real normm[])
101 real mtot;
102 int i, j, k, l;
104 for (i = 0; i < top.mols.nr; i++)
106 k = top.mols.index[i];
107 l = top.mols.index[i+1];
108 mtot = 0.0;
110 for (j = k; j < l; j++)
112 mtot += top.atoms.atom[j].m;
115 for (j = k; j < l; j++)
117 normm[j] = top.atoms.atom[j].m/mtot;
124 static void calc_spectrum(int n, real c[], real dt, const char *fn,
125 gmx_output_env_t *oenv, gmx_bool bRecip)
127 FILE *fp;
128 gmx_fft_t fft;
129 int i, status;
130 real *data;
131 real nu, omega, recip_fac;
133 snew(data, n*2);
134 for (i = 0; (i < n); i++)
136 data[i] = c[i];
139 if ((status = gmx_fft_init_1d_real(&fft, n, GMX_FFT_FLAG_NONE)) != 0)
141 gmx_fatal(FARGS, "Invalid fft return status %d", status);
143 if ((status = gmx_fft_1d_real(fft, GMX_FFT_REAL_TO_COMPLEX, data, data)) != 0)
145 gmx_fatal(FARGS, "Invalid fft return status %d", status);
147 fp = xvgropen(fn, "Vibrational Power Spectrum",
148 bRecip ? "\\f{12}w\\f{4} (cm\\S-1\\N)" :
149 "\\f{12}n\\f{4} (ps\\S-1\\N)",
150 "a.u.", oenv);
151 /* This is difficult.
152 * The length of the ACF is dt (as passed to this routine).
153 * We pass the vacf with N time steps from 0 to dt.
154 * That means that after FFT we have lowest frequency = 1/dt
155 * then 1/(2 dt) etc. (this is the X-axis of the data after FFT).
156 * To convert to 1/cm we need to have to realize that
157 * E = hbar w = h nu = h c/lambda. We want to have reciprokal cm
158 * on the x-axis, that is 1/lambda, so we then have
159 * 1/lambda = nu/c. Since nu has units of 1/ps and c has gromacs units
160 * of nm/ps, we need to multiply by 1e7.
161 * The timestep between saving the trajectory is
162 * 1e7 is to convert nanometer to cm
164 recip_fac = bRecip ? (1e7/SPEED_OF_LIGHT) : 1.0;
165 for (i = 0; (i < n); i += 2)
167 nu = i/(2*dt);
168 omega = nu*recip_fac;
169 /* Computing the square magnitude of a complex number, since this is a power
170 * spectrum.
172 fprintf(fp, "%10g %10g\n", omega, gmx::square(data[i])+gmx::square(data[i+1]));
174 xvgrclose(fp);
175 gmx_fft_destroy(fft);
176 sfree(data);
179 int gmx_velacc(int argc, char *argv[])
181 const char *desc[] = {
182 "[THISMODULE] computes the velocity autocorrelation function.",
183 "When the [TT]-m[tt] option is used, the momentum autocorrelation",
184 "function is calculated.[PAR]",
185 "With option [TT]-mol[tt] the velocity autocorrelation function of",
186 "molecules is calculated. In this case the index group should consist",
187 "of molecule numbers instead of atom numbers.[PAR]",
188 "Be sure that your trajectory contains frames with velocity information",
189 "(i.e. [TT]nstvout[tt] was set in your original [REF].mdp[ref] file),",
190 "and that the time interval between data collection points is",
191 "much shorter than the time scale of the autocorrelation."
194 static gmx_bool bMass = FALSE, bMol = FALSE, bRecip = TRUE;
195 t_pargs pa[] = {
196 { "-m", FALSE, etBOOL, {&bMass},
197 "Calculate the momentum autocorrelation function" },
198 { "-recip", FALSE, etBOOL, {&bRecip},
199 "Use cm^-1 on X-axis instead of 1/ps for spectra." },
200 { "-mol", FALSE, etBOOL, {&bMol},
201 "Calculate the velocity acf of molecules" }
204 t_topology top;
205 int ePBC = -1;
206 t_trxframe fr;
207 matrix box;
208 gmx_bool bTPS = FALSE, bTop = FALSE;
209 int gnx;
210 int *index;
211 char *grpname;
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;
216 t_trxstatus *status;
217 int counter, n_alloc, i, j, counter_dim, k, l;
218 rvec mv_mol;
219 /* Array for the correlation function */
220 real **c1;
221 real *normm = nullptr;
222 gmx_output_env_t *oenv;
224 #define NHISTO 360
226 t_filenm fnm[] = {
227 { efTRN, "-f", nullptr, ffREAD },
228 { efTPS, nullptr, nullptr, ffOPTRD },
229 { efNDX, nullptr, nullptr, ffOPTRD },
230 { efXVG, "-o", "vac", ffWRITE },
231 { efXVG, "-os", "spectrum", ffOPTWR }
233 #define NFILE asize(fnm)
234 int npargs;
235 t_pargs *ppa;
237 npargs = asize(pa);
238 ppa = add_acf_pargs(&npargs, pa);
239 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME,
240 NFILE, fnm, npargs, ppa, asize(desc), desc, 0, nullptr, &oenv))
242 sfree(ppa);
243 return 0;
246 if (bMol || bMass)
248 bTPS = ftp2bSet(efTPS, NFILE, fnm) || !ftp2bSet(efNDX, NFILE, fnm);
251 if (bTPS)
253 bTop = read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &ePBC, nullptr, nullptr, box,
254 TRUE);
255 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
257 else
259 rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
262 if (bMol)
264 if (!bTop)
266 gmx_fatal(FARGS, "Need a topology to determine the molecules");
268 snew(normm, top.atoms.nr);
269 precalc(top, normm);
270 index_atom2mol(&gnx, index, &top.mols);
273 /* Correlation stuff */
274 snew(c1, gnx);
275 for (i = 0; (i < gnx); i++)
277 c1[i] = nullptr;
280 read_first_frame(oenv, &status, ftp2fn(efTRN, NFILE, fnm), &fr, TRX_NEED_V);
281 t0 = fr.time;
283 n_alloc = 0;
284 counter = 0;
287 if (counter >= n_alloc)
289 n_alloc += 100;
290 for (i = 0; i < gnx; i++)
292 srenew(c1[i], DIM*n_alloc);
295 counter_dim = DIM*counter;
296 if (bMol)
298 for (i = 0; i < gnx; i++)
300 clear_rvec(mv_mol);
301 k = top.mols.index[index[i]];
302 l = top.mols.index[index[i]+1];
303 for (j = k; j < l; j++)
305 if (bMass)
307 mass = top.atoms.atom[j].m;
309 else
311 mass = normm[j];
313 mv_mol[XX] += mass*fr.v[j][XX];
314 mv_mol[YY] += mass*fr.v[j][YY];
315 mv_mol[ZZ] += mass*fr.v[j][ZZ];
317 c1[i][counter_dim+XX] = mv_mol[XX];
318 c1[i][counter_dim+YY] = mv_mol[YY];
319 c1[i][counter_dim+ZZ] = mv_mol[ZZ];
322 else
324 for (i = 0; i < gnx; i++)
326 if (bMass)
328 mass = top.atoms.atom[index[i]].m;
330 else
332 mass = 1;
334 c1[i][counter_dim+XX] = mass*fr.v[index[i]][XX];
335 c1[i][counter_dim+YY] = mass*fr.v[index[i]][YY];
336 c1[i][counter_dim+ZZ] = mass*fr.v[index[i]][ZZ];
340 t1 = fr.time;
342 counter++;
344 while (read_next_frame(oenv, status, &fr));
346 close_trx(status);
348 if (counter >= 4)
350 /* Compute time step between frames */
351 dt = (t1-t0)/(counter-1);
352 do_autocorr(opt2fn("-o", NFILE, fnm), oenv,
353 bMass ?
354 "Momentum Autocorrelation Function" :
355 "Velocity Autocorrelation Function",
356 counter, gnx, c1, dt, eacVector, TRUE);
358 do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy");
360 if (opt2bSet("-os", NFILE, fnm))
362 calc_spectrum(counter/2, (real *) (c1[0]), (t1-t0)/2, opt2fn("-os", NFILE, fnm),
363 oenv, bRecip);
364 do_view(oenv, opt2fn("-os", NFILE, fnm), "-nxy");
367 else
369 fprintf(stderr, "Not enough frames in trajectory - no output generated.\n");
372 return 0;