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