Merge branch release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_velacc.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
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.
14
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.
19  *
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.
26  *
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.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Green Red Orange Magenta Azure Cyan Skyblue
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38 #include <stdio.h>
39 #include <math.h>
40
41 #include "confio.h"
42 #include "gmx_fatal.h"
43 #include "futil.h"
44 #include "gstat.h"
45 #include "macros.h"
46 #include "maths.h"
47 #include "physics.h"
48 #include "index.h"
49 #include "smalloc.h"
50 #include "statutil.h"
51 #include <string.h>
52 #include "sysstuff.h"
53 #include "txtdump.h"
54 #include "typedefs.h"
55 #include "vec.h"
56 #include "strdb.h"
57 #include "xvgr.h"
58 #include "gmx_ana.h"
59 #include "gromacs/fft/fft.h"
60
61 static void index_atom2mol(int *n, atom_id *index, t_block *mols)
62 {
63     int nat, i, nmol, mol, j;
64
65     nat  = *n;
66     i    = 0;
67     nmol = 0;
68     mol  = 0;
69     while (i < nat)
70     {
71         while (index[i] > mols->index[mol])
72         {
73             mol++;
74             if (mol >= mols->nr)
75             {
76                 gmx_fatal(FARGS, "Atom index out of range: %d", index[i]+1);
77             }
78         }
79         for (j = mols->index[mol]; j < mols->index[mol+1]; j++)
80         {
81             if (i >= nat || index[i] != j)
82             {
83                 gmx_fatal(FARGS, "The index group does not consist of whole molecules");
84             }
85             i++;
86         }
87         index[nmol++] = mol;
88     }
89
90     fprintf(stderr, "\nSplit group of %d atoms into %d molecules\n", nat, nmol);
91
92     *n = nmol;
93 }
94
95 static void precalc(t_topology top, real normm[])
96 {
97
98     real mtot;
99     int  i, j, k, l;
100
101     for (i = 0; i < top.mols.nr; i++)
102     {
103         k    = top.mols.index[i];
104         l    = top.mols.index[i+1];
105         mtot = 0.0;
106
107         for (j = k; j < l; j++)
108         {
109             mtot += top.atoms.atom[j].m;
110         }
111
112         for (j = k; j < l; j++)
113         {
114             normm[j] = top.atoms.atom[j].m/mtot;
115         }
116
117     }
118
119 }
120
121 static void calc_spectrum(int n, real c[], real dt, const char *fn,
122                           output_env_t oenv, gmx_bool bRecip)
123 {
124     FILE     *fp;
125     gmx_fft_t fft;
126     int       i, status;
127     real     *data;
128     real      nu, omega, recip_fac;
129
130     snew(data, n*2);
131     for (i = 0; (i < n); i++)
132     {
133         data[i] = c[i];
134     }
135
136     if ((status = gmx_fft_init_1d_real(&fft, n, GMX_FFT_FLAG_NONE)) != 0)
137     {
138         gmx_fatal(FARGS, "Invalid fft return status %d", status);
139     }
140     if ((status = gmx_fft_1d_real(fft, GMX_FFT_REAL_TO_COMPLEX, data, data)) != 0)
141     {
142         gmx_fatal(FARGS, "Invalid fft return status %d", status);
143     }
144     fp = xvgropen(fn, "Vibrational Power Spectrum",
145                   bRecip ? "\\f{12}w\\f{4} (cm\\S-1\\N)" :
146                   "\\f{12}n\\f{4} (ps\\S-1\\N)",
147                   "a.u.", oenv);
148     /* This is difficult.
149      * The length of the ACF is dt (as passed to this routine).
150      * We pass the vacf with N time steps from 0 to dt.
151      * That means that after FFT we have lowest frequency = 1/dt
152      * then 1/(2 dt) etc. (this is the X-axis of the data after FFT).
153      * To convert to 1/cm we need to have to realize that
154      * E = hbar w = h nu = h c/lambda. We want to have reciprokal cm
155      * on the x-axis, that is 1/lambda, so we then have
156      * 1/lambda = nu/c. Since nu has units of 1/ps and c has gromacs units
157      * of nm/ps, we need to multiply by 1e7.
158      * The timestep between saving the trajectory is
159      * 1e7 is to convert nanometer to cm
160      */
161     recip_fac = bRecip ? (1e7/SPEED_OF_LIGHT) : 1.0;
162     for (i = 0; (i < n); i += 2)
163     {
164         nu    = i/(2*dt);
165         omega = nu*recip_fac;
166         /* Computing the square magnitude of a complex number, since this is a power
167          * spectrum.
168          */
169         fprintf(fp, "%10g  %10g\n", omega, sqr(data[i])+sqr(data[i+1]));
170     }
171     xvgrclose(fp);
172     gmx_fft_destroy(fft);
173     sfree(data);
174 }
175
176 int gmx_velacc(int argc, char *argv[])
177 {
178     const char     *desc[] = {
179         "[TT]g_velacc[tt] computes the velocity autocorrelation function.",
180         "When the [TT]-m[tt] option is used, the momentum autocorrelation",
181         "function is calculated.[PAR]",
182         "With option [TT]-mol[tt] the velocity autocorrelation function of",
183         "molecules is calculated. In this case the index group should consist",
184         "of molecule numbers instead of atom numbers.[PAR]",
185         "Be sure that your trajectory contains frames with velocity information",
186         "(i.e. [TT]nstvout[tt] was set in your original [TT].mdp[tt] file),",
187         "and that the time interval between data collection points is",
188         "much shorter than the time scale of the autocorrelation."
189     };
190
191     static gmx_bool bMass = FALSE, bMol = FALSE, bRecip = TRUE;
192     t_pargs         pa[]  = {
193         { "-m", FALSE, etBOOL, {&bMass},
194           "Calculate the momentum autocorrelation function" },
195         { "-recip", FALSE, etBOOL, {&bRecip},
196           "Use cm^-1 on X-axis instead of 1/ps for spectra." },
197         { "-mol", FALSE, etBOOL, {&bMol},
198           "Calculate the velocity acf of molecules" }
199     };
200
201     t_topology      top;
202     int             ePBC = -1;
203     t_trxframe      fr;
204     matrix          box;
205     gmx_bool        bTPS = FALSE, bTop = FALSE;
206     int             gnx;
207     atom_id        *index;
208     char           *grpname;
209     char            title[256];
210     /* t0, t1 are the beginning and end time respectively.
211      * dt is the time step, mass is temp variable for atomic mass.
212      */
213     real              t0, t1, dt, mass;
214     t_trxstatus      *status;
215     int               counter, n_alloc, i, j, counter_dim, k, l;
216     rvec              mv_mol;
217     /* Array for the correlation function */
218     real            **c1;
219     real             *normm = NULL;
220     output_env_t      oenv;
221
222 #define NHISTO 360
223
224     t_filenm  fnm[] = {
225         { efTRN, "-f",    NULL,   ffREAD  },
226         { efTPS, NULL,    NULL,   ffOPTRD },
227         { efNDX, NULL,    NULL,   ffOPTRD },
228         { efXVG, "-o",    "vac",  ffWRITE },
229         { efXVG, "-os",   "spectrum", ffOPTWR }
230     };
231 #define NFILE asize(fnm)
232     int       npargs;
233     t_pargs  *ppa;
234
235     npargs = asize(pa);
236     ppa    = add_acf_pargs(&npargs, pa);
237     parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
238                       NFILE, fnm, npargs, ppa, asize(desc), desc, 0, NULL, &oenv);
239
240     if (bMol || bMass)
241     {
242         bTPS = ftp2bSet(efTPS, NFILE, fnm) || !ftp2bSet(efNDX, NFILE, fnm);
243     }
244
245     if (bTPS)
246     {
247         bTop = read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, NULL, NULL, box,
248                              TRUE);
249         get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
250     }
251     else
252     {
253         rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
254     }
255
256     if (bMol)
257     {
258         if (!bTop)
259         {
260             gmx_fatal(FARGS, "Need a topology to determine the molecules");
261         }
262         snew(normm, top.atoms.nr);
263         precalc(top, normm);
264         index_atom2mol(&gnx, index, &top.mols);
265     }
266
267     /* Correlation stuff */
268     snew(c1, gnx);
269     for (i = 0; (i < gnx); i++)
270     {
271         c1[i] = NULL;
272     }
273
274     read_first_frame(oenv, &status, ftp2fn(efTRN, NFILE, fnm), &fr, TRX_NEED_V);
275     t0 = fr.time;
276
277     n_alloc = 0;
278     counter = 0;
279     do
280     {
281         if (counter >= n_alloc)
282         {
283             n_alloc += 100;
284             for (i = 0; i < gnx; i++)
285             {
286                 srenew(c1[i], DIM*n_alloc);
287             }
288         }
289         counter_dim = DIM*counter;
290         if (bMol)
291         {
292             for (i = 0; i < gnx; i++)
293             {
294                 clear_rvec(mv_mol);
295                 k = top.mols.index[index[i]];
296                 l = top.mols.index[index[i]+1];
297                 for (j = k; j < l; j++)
298                 {
299                     if (bMass)
300                     {
301                         mass = top.atoms.atom[j].m;
302                     }
303                     else
304                     {
305                         mass = normm[j];
306                     }
307                     mv_mol[XX] += mass*fr.v[j][XX];
308                     mv_mol[YY] += mass*fr.v[j][YY];
309                     mv_mol[ZZ] += mass*fr.v[j][ZZ];
310                 }
311                 c1[i][counter_dim+XX] = mv_mol[XX];
312                 c1[i][counter_dim+YY] = mv_mol[YY];
313                 c1[i][counter_dim+ZZ] = mv_mol[ZZ];
314             }
315         }
316         else
317         {
318             for (i = 0; i < gnx; i++)
319             {
320                 if (bMass)
321                 {
322                     mass = top.atoms.atom[index[i]].m;
323                 }
324                 else
325                 {
326                     mass = 1;
327                 }
328                 c1[i][counter_dim+XX] = mass*fr.v[index[i]][XX];
329                 c1[i][counter_dim+YY] = mass*fr.v[index[i]][YY];
330                 c1[i][counter_dim+ZZ] = mass*fr.v[index[i]][ZZ];
331             }
332         }
333
334         t1 = fr.time;
335
336         counter++;
337     }
338     while (read_next_frame(oenv, status, &fr));
339
340     close_trj(status);
341
342     if (counter >= 4)
343     {
344         /* Compute time step between frames */
345         dt = (t1-t0)/(counter-1);
346         do_autocorr(opt2fn("-o", NFILE, fnm), oenv,
347                     bMass ?
348                     "Momentum Autocorrelation Function" :
349                     "Velocity Autocorrelation Function",
350                     counter, gnx, c1, dt, eacVector, TRUE);
351
352         do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy");
353
354         if (opt2bSet("-os", NFILE, fnm))
355         {
356             calc_spectrum(counter/2, (real *) (c1[0]), (t1-t0)/2, opt2fn("-os", NFILE, fnm),
357                           oenv, bRecip);
358             do_view(oenv, opt2fn("-os", NFILE, fnm), "-nxy");
359         }
360     }
361     else
362     {
363         fprintf(stderr, "Not enough frames in trajectory - no output generated.\n");
364     }
365
366     return 0;
367 }