Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_velacc.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include <math.h>
40 #include <stdio.h>
41 #include <string.h>
42
43 #include "gromacs/commandline/pargs.h"
44 #include "gromacs/fft/fft.h"
45 #include "gromacs/fileio/confio.h"
46 #include "gromacs/fileio/trxio.h"
47 #include "gromacs/fileio/xvgr.h"
48 #include "gromacs/gmxana/gmx_ana.h"
49 #include "gromacs/gmxana/gstat.h"
50 #include "gromacs/legacyheaders/macros.h"
51 #include "gromacs/legacyheaders/txtdump.h"
52 #include "gromacs/legacyheaders/typedefs.h"
53 #include "gromacs/legacyheaders/viewit.h"
54 #include "gromacs/math/units.h"
55 #include "gromacs/math/utilities.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/topology/index.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/futil.h"
60 #include "gromacs/utility/smalloc.h"
61
62 static void index_atom2mol(int *n, atom_id *index, t_block *mols)
63 {
64     int nat, i, nmol, mol, j;
65
66     nat  = *n;
67     i    = 0;
68     nmol = 0;
69     mol  = 0;
70     while (i < nat)
71     {
72         while (index[i] > mols->index[mol])
73         {
74             mol++;
75             if (mol >= mols->nr)
76             {
77                 gmx_fatal(FARGS, "Atom index out of range: %d", index[i]+1);
78             }
79         }
80         for (j = mols->index[mol]; j < mols->index[mol+1]; j++)
81         {
82             if (i >= nat || index[i] != j)
83             {
84                 gmx_fatal(FARGS, "The index group does not consist of whole molecules");
85             }
86             i++;
87         }
88         index[nmol++] = mol;
89     }
90
91     fprintf(stderr, "\nSplit group of %d atoms into %d molecules\n", nat, nmol);
92
93     *n = nmol;
94 }
95
96 static void precalc(t_topology top, real normm[])
97 {
98
99     real mtot;
100     int  i, j, k, l;
101
102     for (i = 0; i < top.mols.nr; i++)
103     {
104         k    = top.mols.index[i];
105         l    = top.mols.index[i+1];
106         mtot = 0.0;
107
108         for (j = k; j < l; j++)
109         {
110             mtot += top.atoms.atom[j].m;
111         }
112
113         for (j = k; j < l; j++)
114         {
115             normm[j] = top.atoms.atom[j].m/mtot;
116         }
117
118     }
119
120 }
121
122 static void calc_spectrum(int n, real c[], real dt, const char *fn,
123                           output_env_t oenv, gmx_bool bRecip)
124 {
125     FILE     *fp;
126     gmx_fft_t fft;
127     int       i, status;
128     real     *data;
129     real      nu, omega, recip_fac;
130
131     snew(data, n*2);
132     for (i = 0; (i < n); i++)
133     {
134         data[i] = c[i];
135     }
136
137     if ((status = gmx_fft_init_1d_real(&fft, n, GMX_FFT_FLAG_NONE)) != 0)
138     {
139         gmx_fatal(FARGS, "Invalid fft return status %d", status);
140     }
141     if ((status = gmx_fft_1d_real(fft, GMX_FFT_REAL_TO_COMPLEX, data, data)) != 0)
142     {
143         gmx_fatal(FARGS, "Invalid fft return status %d", status);
144     }
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)",
148                   "a.u.", oenv);
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
161      */
162     recip_fac = bRecip ? (1e7/SPEED_OF_LIGHT) : 1.0;
163     for (i = 0; (i < n); i += 2)
164     {
165         nu    = i/(2*dt);
166         omega = nu*recip_fac;
167         /* Computing the square magnitude of a complex number, since this is a power
168          * spectrum.
169          */
170         fprintf(fp, "%10g  %10g\n", omega, sqr(data[i])+sqr(data[i+1]));
171     }
172     xvgrclose(fp);
173     gmx_fft_destroy(fft);
174     sfree(data);
175 }
176
177 int gmx_velacc(int argc, char *argv[])
178 {
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."
190     };
191
192     static gmx_bool bMass = FALSE, bMol = FALSE, bRecip = TRUE;
193     t_pargs         pa[]  = {
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" }
200     };
201
202     t_topology      top;
203     int             ePBC = -1;
204     t_trxframe      fr;
205     matrix          box;
206     gmx_bool        bTPS = FALSE, bTop = FALSE;
207     int             gnx;
208     atom_id        *index;
209     char           *grpname;
210     char            title[256];
211     /* t0, t1 are the beginning and end time respectively.
212      * dt is the time step, mass is temp variable for atomic mass.
213      */
214     real              t0, t1, dt, mass;
215     t_trxstatus      *status;
216     int               counter, n_alloc, i, j, counter_dim, k, l;
217     rvec              mv_mol;
218     /* Array for the correlation function */
219     real            **c1;
220     real             *normm = NULL;
221     output_env_t      oenv;
222
223 #define NHISTO 360
224
225     t_filenm  fnm[] = {
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 }
231     };
232 #define NFILE asize(fnm)
233     int       npargs;
234     t_pargs  *ppa;
235
236     npargs = asize(pa);
237     ppa    = add_acf_pargs(&npargs, pa);
238     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME,
239                            NFILE, fnm, npargs, ppa, asize(desc), desc, 0, NULL, &oenv))
240     {
241         return 0;
242     }
243
244     if (bMol || bMass)
245     {
246         bTPS = ftp2bSet(efTPS, NFILE, fnm) || !ftp2bSet(efNDX, NFILE, fnm);
247     }
248
249     if (bTPS)
250     {
251         bTop = read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, NULL, NULL, box,
252                              TRUE);
253         get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
254     }
255     else
256     {
257         rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
258     }
259
260     if (bMol)
261     {
262         if (!bTop)
263         {
264             gmx_fatal(FARGS, "Need a topology to determine the molecules");
265         }
266         snew(normm, top.atoms.nr);
267         precalc(top, normm);
268         index_atom2mol(&gnx, index, &top.mols);
269     }
270
271     /* Correlation stuff */
272     snew(c1, gnx);
273     for (i = 0; (i < gnx); i++)
274     {
275         c1[i] = NULL;
276     }
277
278     read_first_frame(oenv, &status, ftp2fn(efTRN, NFILE, fnm), &fr, TRX_NEED_V);
279     t0 = fr.time;
280
281     n_alloc = 0;
282     counter = 0;
283     do
284     {
285         if (counter >= n_alloc)
286         {
287             n_alloc += 100;
288             for (i = 0; i < gnx; i++)
289             {
290                 srenew(c1[i], DIM*n_alloc);
291             }
292         }
293         counter_dim = DIM*counter;
294         if (bMol)
295         {
296             for (i = 0; i < gnx; i++)
297             {
298                 clear_rvec(mv_mol);
299                 k = top.mols.index[index[i]];
300                 l = top.mols.index[index[i]+1];
301                 for (j = k; j < l; j++)
302                 {
303                     if (bMass)
304                     {
305                         mass = top.atoms.atom[j].m;
306                     }
307                     else
308                     {
309                         mass = normm[j];
310                     }
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];
314                 }
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];
318             }
319         }
320         else
321         {
322             for (i = 0; i < gnx; i++)
323             {
324                 if (bMass)
325                 {
326                     mass = top.atoms.atom[index[i]].m;
327                 }
328                 else
329                 {
330                     mass = 1;
331                 }
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];
335             }
336         }
337
338         t1 = fr.time;
339
340         counter++;
341     }
342     while (read_next_frame(oenv, status, &fr));
343
344     close_trj(status);
345
346     if (counter >= 4)
347     {
348         /* Compute time step between frames */
349         dt = (t1-t0)/(counter-1);
350         do_autocorr(opt2fn("-o", NFILE, fnm), oenv,
351                     bMass ?
352                     "Momentum Autocorrelation Function" :
353                     "Velocity Autocorrelation Function",
354                     counter, gnx, c1, dt, eacVector, TRUE);
355
356         do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy");
357
358         if (opt2bSet("-os", NFILE, fnm))
359         {
360             calc_spectrum(counter/2, (real *) (c1[0]), (t1-t0)/2, opt2fn("-os", NFILE, fnm),
361                           oenv, bRecip);
362             do_view(oenv, opt2fn("-os", NFILE, fnm), "-nxy");
363         }
364     }
365     else
366     {
367         fprintf(stderr, "Not enough frames in trajectory - no output generated.\n");
368     }
369
370     return 0;
371 }