Fix remaining copyright headers
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_tcaf.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, 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 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40 #include <stdio.h>
41
42 #include <math.h>
43 #include "gromacs/fileio/confio.h"
44 #include "gmx_fatal.h"
45 #include "gromacs/fileio/futil.h"
46 #include "gstat.h"
47 #include "macros.h"
48 #include "maths.h"
49 #include "physics.h"
50 #include "index.h"
51 #include "smalloc.h"
52 #include "gromacs/commandline/pargs.h"
53 #include <string.h>
54 #include "sysstuff.h"
55 #include "txtdump.h"
56 #include "typedefs.h"
57 #include "vec.h"
58 #include "strdb.h"
59 #include "xvgr.h"
60 #include "pbc.h"
61 #include "gmx_ana.h"
62 #include "gromacs/fileio/trxio.h"
63
64
65 #define NK  24
66 #define NPK 4
67
68 #define NKC  6
69 #define NKC0 4
70 int  kset_c[NKC+1] = { 0, 3, 9, 13, 16, 19, NK };
71
72 rvec v0[NK] = {{1, 0, 0}, {0, 1, 0}, {0, 0, 1}, {1, 1, 0}, {1, -1, 0}, {1, 0, 1}, {1, 0, -1}, {0, 1, 1}, {0, 1, -1}, {1, 1, 1}, {1, 1, -1}, {1, -1, 1}, {-1, 1, 1}, {2, 0, 0}, {0, 2, 0}, {0, 0, 2}, {3, 0, 0}, {0, 3, 0}, {0, 0, 3}, {4, 0, 0}, {0, 4, 0}, {0, 0, 4}};
73 rvec v1[NK] = {{0, 1, 0}, {0, 0, 1}, {1, 0, 0}, {0, 0, 1}, {0, 0, 1}, {0, 1, 0}, {0, 1, 0}, {1, 0, 0}, {1, 0, 0}, {1, -1, 0}, {1, -1, 0}, {1, 0, -1}, { 0, 1, -1}, {0, 1, 0}, {0, 0, 1}, {1, 0, 0}, {0, 1, 0}, {0, 0, 1}, {1, 0, 0}, {0, 1, 0}, {0, 0, 1}, {1, 0, 0}};
74 rvec v2[NK] = {{0, 0, 1}, {1, 0, 0}, {0, 1, 0}, {1, -1, 0}, {1, 1, 0}, {1, 0, -1}, {1, 0, 1}, {0, 1, -1}, {0, 1, 1}, {1, 1, -2}, {1, 1, 2}, {1, 2, 1}, { 2, 1, 1}, {0, 0, 1}, {1, 0, 0}, {0, 1, 0}, {0, 0, 1}, {1, 0, 0}, {0, 1, 0}, {0, 0, 1}, {1, 0, 0}, {0, 1, 0}};
75
76 static void process_tcaf(int nframes, real dt, int nkc, real **tc, rvec *kfac,
77                          real rho, real wt, const char *fn_trans,
78                          const char *fn_tca, const char *fn_tc,
79                          const char *fn_tcf, const char *fn_cub,
80                          const char *fn_vk, const output_env_t oenv)
81 {
82     FILE  *fp, *fp_vk, *fp_cub = NULL;
83     int    nk, ntc;
84     real **tcaf, **tcafc = NULL, eta;
85     int    i, j, k, kc;
86     int    ncorr;
87     real   fitparms[3], *sig;
88
89     nk  = kset_c[nkc];
90     ntc = nk*NPK;
91
92     if (fn_trans)
93     {
94         fp = xvgropen(fn_trans, "Transverse Current", "Time (ps)", "TC (nm/ps)",
95                       oenv);
96         for (i = 0; i < nframes; i++)
97         {
98             fprintf(fp, "%g", i*dt);
99             for (j = 0; j < ntc; j++)
100             {
101                 fprintf(fp, " %g", tc[j][i]);
102             }
103             fprintf(fp, "\n");
104         }
105         ffclose(fp);
106         do_view(oenv, fn_trans, "-nxy");
107     }
108
109     ncorr = (nframes+1)/2;
110     if (ncorr > (int)(5*wt/dt+0.5))
111     {
112         ncorr = (int)(5*wt/dt+0.5)+1;
113     }
114     snew(tcaf, nk);
115     for (k = 0; k < nk; k++)
116     {
117         snew(tcaf[k], ncorr);
118     }
119     if (fn_cub)
120     {
121         snew(tcafc, nkc);
122         for (k = 0; k < nkc; k++)
123         {
124             snew(tcafc[k], ncorr);
125         }
126     }
127     snew(sig, ncorr);
128     for (i = 0; i < ncorr; i++)
129     {
130         sig[i] = exp(0.5*i*dt/wt);
131     }
132
133     low_do_autocorr(fn_tca, oenv, "Transverse Current Autocorrelation Functions",
134                     nframes, ntc, ncorr, tc, dt, eacNormal,
135                     1, FALSE, FALSE, FALSE, 0, 0, 0);
136     do_view(oenv, fn_tca, "-nxy");
137
138     fp = xvgropen(fn_tc, "Transverse Current Autocorrelation Functions",
139                   "Time (ps)", "TCAF", oenv);
140     for (i = 0; i < ncorr; i++)
141     {
142         kc = 0;
143         fprintf(fp, "%g", i*dt);
144         for (k = 0; k < nk; k++)
145         {
146             for (j = 0; j < NPK; j++)
147             {
148                 tcaf[k][i] += tc[NPK*k+j][i];
149             }
150             if (fn_cub)
151             {
152                 for (j = 0; j < NPK; j++)
153                 {
154                     tcafc[kc][i] += tc[NPK*k+j][i];
155                 }
156             }
157             if (i == 0)
158             {
159                 fprintf(fp, " %g", 1.0);
160             }
161             else
162             {
163                 tcaf[k][i] /= tcaf[k][0];
164                 fprintf(fp, " %g", tcaf[k][i]);
165             }
166             if (k+1 == kset_c[kc+1])
167             {
168                 kc++;
169             }
170         }
171         fprintf(fp, "\n");
172     }
173     ffclose(fp);
174     do_view(oenv, fn_tc, "-nxy");
175
176     if (fn_cub)
177     {
178         fp_cub = xvgropen(fn_cub, "TCAFs and fits", "Time (ps)", "TCAF", oenv);
179         for (kc = 0; kc < nkc; kc++)
180         {
181             fprintf(fp_cub, "%g %g\n", 0.0, 1.0);
182             for (i = 1; i < ncorr; i++)
183             {
184                 tcafc[kc][i] /= tcafc[kc][0];
185                 fprintf(fp_cub, "%g %g\n", i*dt, tcafc[kc][i]);
186             }
187             fprintf(fp_cub, "&\n");
188             tcafc[kc][0] = 1.0;
189         }
190     }
191
192     fp_vk = xvgropen(fn_vk, "Fits", "k (nm\\S-1\\N)",
193                      "\\8h\\4 (10\\S-3\\N kg m\\S-1\\N s\\S-1\\N)", oenv);
194     fprintf(fp_vk, "@    s0 symbol 2\n");
195     fprintf(fp_vk, "@    s0 symbol color 1\n");
196     fprintf(fp_vk, "@    s0 linestyle 0\n");
197     if (fn_cub)
198     {
199         fprintf(fp_vk, "@    s1 symbol 3\n");
200         fprintf(fp_vk, "@    s1 symbol color 2\n");
201     }
202     fp = xvgropen(fn_tcf, "TCAF Fits", "Time (ps)", "", oenv);
203     for (k = 0; k < nk; k++)
204     {
205         tcaf[k][0]   = 1.0;
206         fitparms[0]  = 1;
207         fitparms[1]  = 1;
208         do_lmfit(ncorr, tcaf[k], sig, dt, 0, 0, ncorr*dt,
209                  oenv, bDebugMode(), effnVAC, fitparms, 0);
210         eta = 1000*fitparms[1]*rho/
211             (4*fitparms[0]*PICO*norm2(kfac[k])/(NANO*NANO));
212         fprintf(stdout, "k %6.3f  tau %6.3f  eta %8.5f 10^-3 kg/(m s)\n",
213                 norm(kfac[k]), fitparms[0], eta);
214         fprintf(fp_vk, "%6.3f %g\n", norm(kfac[k]), eta);
215         for (i = 0; i < ncorr; i++)
216         {
217             fprintf(fp, "%g %g\n", i*dt, fit_function(effnVAC, fitparms, i*dt));
218         }
219         fprintf(fp, "&\n");
220     }
221     ffclose(fp);
222     do_view(oenv, fn_tcf, "-nxy");
223
224     if (fn_cub)
225     {
226         fprintf(stdout, "Averaged over k-vectors:\n");
227         fprintf(fp_vk, "&\n");
228         for (k = 0; k < nkc; k++)
229         {
230             tcafc[k][0]  = 1.0;
231             fitparms[0]  = 1;
232             fitparms[1]  = 1;
233             do_lmfit(ncorr, tcafc[k], sig, dt, 0, 0, ncorr*dt,
234                      oenv, bDebugMode(), effnVAC, fitparms, 0);
235             eta = 1000*fitparms[1]*rho/
236                 (4*fitparms[0]*PICO*norm2(kfac[kset_c[k]])/(NANO*NANO));
237             fprintf(stdout,
238                     "k %6.3f  tau %6.3f  Omega %6.3f  eta %8.5f 10^-3 kg/(m s)\n",
239                     norm(kfac[kset_c[k]]), fitparms[0], fitparms[1], eta);
240             fprintf(fp_vk, "%6.3f %g\n", norm(kfac[kset_c[k]]), eta);
241             for (i = 0; i < ncorr; i++)
242             {
243                 fprintf(fp_cub, "%g %g\n", i*dt, fit_function(effnVAC, fitparms, i*dt));
244             }
245             fprintf(fp_cub, "&\n");
246         }
247         fprintf(fp_vk, "&\n");
248         ffclose(fp_cub);
249         do_view(oenv, fn_cub, "-nxy");
250     }
251     ffclose(fp_vk);
252     do_view(oenv, fn_vk, "-nxy");
253 }
254
255
256 int gmx_tcaf(int argc, char *argv[])
257 {
258     const char     *desc[] = {
259         "[THISMODULE] computes tranverse current autocorrelations.",
260         "These are used to estimate the shear viscosity, [GRK]eta[grk].",
261         "For details see: Palmer, Phys. Rev. E 49 (1994) pp 359-366.[PAR]",
262         "Transverse currents are calculated using the",
263         "k-vectors (1,0,0) and (2,0,0) each also in the [IT]y[it]- and [IT]z[it]-direction,",
264         "(1,1,0) and (1,-1,0) each also in the 2 other planes (these vectors",
265         "are not independent) and (1,1,1) and the 3 other box diagonals (also",
266         "not independent). For each k-vector the sine and cosine are used, in",
267         "combination with the velocity in 2 perpendicular directions. This gives",
268         "a total of 16*2*2=64 transverse currents. One autocorrelation is",
269         "calculated fitted for each k-vector, which gives 16 TCAFs. Each of",
270         "these TCAFs is fitted to [MATH]f(t) = [EXP]-v[exp]([COSH]Wv[cosh] + 1/W [SINH]Wv[sinh])[math],",
271         "[MATH]v = -t/(2 [GRK]tau[grk])[math], [MATH]W = [SQRT]1 - 4 [GRK]tau[grk] [GRK]eta[grk]/[GRK]rho[grk] k^2[sqrt][math], which gives 16 values of [GRK]tau[grk]",
272         "and [GRK]eta[grk]. The fit weights decay exponentially with time constant [MATH]w[math] (given with [TT]-wt[tt]) as [MATH][EXP]-t/w[exp][math], and the TCAF and",
273         "fit are calculated up to time [MATH]5*w[math].",
274         "The [GRK]eta[grk] values should be fitted to [MATH]1 - a [GRK]eta[grk](k) k^2[math], from which",
275         "one can estimate the shear viscosity at k=0.[PAR]",
276         "When the box is cubic, one can use the option [TT]-oc[tt], which",
277         "averages the TCAFs over all k-vectors with the same length.",
278         "This results in more accurate TCAFs.",
279         "Both the cubic TCAFs and fits are written to [TT]-oc[tt]",
280         "The cubic [GRK]eta[grk] estimates are also written to [TT]-ov[tt].[PAR]",
281         "With option [TT]-mol[tt], the transverse current is determined of",
282         "molecules instead of atoms. In this case, the index group should",
283         "consist of molecule numbers instead of atom numbers.[PAR]",
284         "The k-dependent viscosities in the [TT]-ov[tt] file should be",
285         "fitted to [MATH][GRK]eta[grk](k) = [GRK]eta[grk][SUB]0[sub] (1 - a k^2)[math] to obtain the viscosity at",
286         "infinite wavelength.[PAR]",
287         "[BB]Note:[bb] make sure you write coordinates and velocities often enough.",
288         "The initial, non-exponential, part of the autocorrelation function",
289         "is very important for obtaining a good fit."
290     };
291
292     static gmx_bool bMol = FALSE, bK34 = FALSE;
293     static real     wt   = 5;
294     t_pargs         pa[] = {
295         { "-mol", FALSE, etBOOL, {&bMol},
296           "Calculate TCAF of molecules" },
297         { "-k34", FALSE, etBOOL, {&bK34},
298           "Also use k=(3,0,0) and k=(4,0,0)" },
299         { "-wt", FALSE, etREAL, {&wt},
300           "Exponential decay time for the TCAF fit weights" }
301     };
302
303     t_topology      top;
304     int             ePBC;
305     t_trxframe      fr;
306     matrix          box;
307     gmx_bool        bTPS, bTop; /* ,bCubic; */
308     int             gnx;
309     atom_id        *index, *atndx = NULL, at;
310     char           *grpname;
311     char            title[256];
312     real            t0, t1, dt, m, mtot, sysmass, rho, sx, cx;
313     t_trxstatus    *status;
314     int             nframes, n_alloc, i, j, k, d;
315     rvec            mv_mol, cm_mol, kfac[NK];
316     int             nkc, nk, ntc;
317     real          **c1, **tc;
318     output_env_t    oenv;
319
320 #define NHISTO 360
321
322     t_filenm  fnm[] = {
323         { efTRN, "-f",    NULL,      ffREAD  },
324         { efTPS, NULL,    NULL,      ffOPTRD },
325         { efNDX, NULL,    NULL,      ffOPTRD },
326         { efXVG, "-ot",  "transcur", ffOPTWR },
327         { efXVG, "-oa",  "tcaf_all", ffWRITE },
328         { efXVG, "-o",   "tcaf",     ffWRITE },
329         { efXVG, "-of",  "tcaf_fit", ffWRITE },
330         { efXVG, "-oc",  "tcaf_cub", ffOPTWR },
331         { efXVG, "-ov",  "visc_k",   ffWRITE }
332     };
333 #define NFILE asize(fnm)
334     int       npargs;
335     t_pargs  *ppa;
336
337     npargs = asize(pa);
338     ppa    = add_acf_pargs(&npargs, pa);
339
340     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
341                            NFILE, fnm, npargs, ppa, asize(desc), desc, 0, NULL, &oenv))
342     {
343         return 0;
344     }
345
346     bTop = read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, NULL, NULL, box,
347                          TRUE);
348     get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
349
350     if (bMol)
351     {
352         if (!bTop)
353         {
354             gmx_fatal(FARGS, "Need a topology to determine the molecules");
355         }
356         atndx = top.mols.index;
357     }
358
359     if (bK34)
360     {
361         nkc = NKC;
362     }
363     else
364     {
365         nkc = NKC0;
366     }
367     nk  = kset_c[nkc];
368     ntc = nk*NPK;
369
370     sprintf(title, "Velocity Autocorrelation Function for %s", grpname);
371
372     sysmass = 0;
373     for (i = 0; i < nk; i++)
374     {
375         if (iprod(v0[i], v1[i]) != 0)
376         {
377             gmx_fatal(FARGS, "DEATH HORROR: vectors not orthogonal");
378         }
379         if (iprod(v0[i], v2[i]) != 0)
380         {
381             gmx_fatal(FARGS, "DEATH HORROR: vectors not orthogonal");
382         }
383         if (iprod(v1[i], v2[i]) != 0)
384         {
385             gmx_fatal(FARGS, "DEATH HORROR: vectors not orthogonal");
386         }
387         unitv(v1[i], v1[i]);
388         unitv(v2[i], v2[i]);
389     }
390     snew(tc, ntc);
391     for (i = 0; i < top.atoms.nr; i++)
392     {
393         sysmass += top.atoms.atom[i].m;
394     }
395
396     read_first_frame(oenv, &status, ftp2fn(efTRN, NFILE, fnm), &fr,
397                      TRX_NEED_X | TRX_NEED_V);
398     t0 = fr.time;
399
400     n_alloc = 0;
401     nframes = 0;
402     rho     = 0;
403     /* bCubic = TRUE; */
404     do
405     {
406         /*
407            bCubic = bCubic && !TRICLINIC(fr.box) &&
408            fabs(fr.box[XX][XX]-fr.box[YY][YY]) < 0.001*fr.box[XX][XX] &&
409            fabs(fr.box[XX][XX]-fr.box[ZZ][ZZ]) < 0.001*fr.box[XX][XX];
410          */
411
412         if (nframes >= n_alloc)
413         {
414             n_alloc += 100;
415             for (i = 0; i < ntc; i++)
416             {
417                 srenew(tc[i], n_alloc);
418             }
419         }
420
421         rho += 1/det(fr.box);
422         for (k = 0; k < nk; k++)
423         {
424             for (d = 0; d < DIM; d++)
425             {
426                 kfac[k][d] = 2*M_PI*v0[k][d]/fr.box[d][d];
427             }
428         }
429         for (i = 0; i < ntc; i++)
430         {
431             tc[i][nframes] = 0;
432         }
433
434         for (i = 0; i < gnx; i++)
435         {
436             if (bMol)
437             {
438                 clear_rvec(mv_mol);
439                 clear_rvec(cm_mol);
440                 mtot = 0;
441                 for (j = 0; j < atndx[index[i]+1] - atndx[index[i]]; j++)
442                 {
443                     at          = atndx[index[i]] + j;
444                     m           = top.atoms.atom[at].m;
445                     mv_mol[XX] += m*fr.v[at][XX];
446                     mv_mol[YY] += m*fr.v[at][YY];
447                     mv_mol[ZZ] += m*fr.v[at][ZZ];
448                     cm_mol[XX] += m*fr.x[at][XX];
449                     cm_mol[YY] += m*fr.x[at][YY];
450                     cm_mol[ZZ] += m*fr.x[at][ZZ];
451                     mtot       += m;
452                 }
453                 svmul(1.0/mtot, cm_mol, cm_mol);
454             }
455             else
456             {
457                 svmul(top.atoms.atom[index[i]].m, fr.v[index[i]], mv_mol);
458             }
459
460             if (!bMol)
461             {
462                 copy_rvec(fr.x[index[i]], cm_mol);
463             }
464             j = 0;
465             for (k = 0; k < nk; k++)
466             {
467                 sx              = sin(iprod(kfac[k], cm_mol));
468                 cx              = cos(iprod(kfac[k], cm_mol));
469                 tc[j][nframes] += sx*iprod(v1[k], mv_mol);
470                 j++;
471                 tc[j][nframes] += cx*iprod(v1[k], mv_mol);
472                 j++;
473                 tc[j][nframes] += sx*iprod(v2[k], mv_mol);
474                 j++;
475                 tc[j][nframes] += cx*iprod(v2[k], mv_mol);
476                 j++;
477             }
478         }
479
480         t1 = fr.time;
481         nframes++;
482     }
483     while (read_next_frame(oenv, status, &fr));
484     close_trj(status);
485
486     dt = (t1-t0)/(nframes-1);
487
488     rho *= sysmass/nframes*AMU/(NANO*NANO*NANO);
489     fprintf(stdout, "Density = %g (kg/m^3)\n", rho);
490     process_tcaf(nframes, dt, nkc, tc, kfac, rho, wt,
491                  opt2fn_null("-ot", NFILE, fnm),
492                  opt2fn("-oa", NFILE, fnm), opt2fn("-o", NFILE, fnm),
493                  opt2fn("-of", NFILE, fnm), opt2fn_null("-oc", NFILE, fnm),
494                  opt2fn("-ov", NFILE, fnm), oenv);
495
496     return 0;
497 }