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