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