Update instructions in containers.rst
[gromacs.git] / src / gromacs / gmxana / gmx_tcaf.cpp
blob03e9f2c6ad4d7f0022bcbbf6bf941b21563e0ad9
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 #include "gmxpre.h"
40 #include <cmath>
41 #include <cstdio>
42 #include <cstring>
44 #include "gromacs/commandline/pargs.h"
45 #include "gromacs/commandline/viewit.h"
46 #include "gromacs/correlationfunctions/autocorr.h"
47 #include "gromacs/correlationfunctions/expfit.h"
48 #include "gromacs/fileio/confio.h"
49 #include "gromacs/fileio/trxio.h"
50 #include "gromacs/fileio/xvgr.h"
51 #include "gromacs/gmxana/gmx_ana.h"
52 #include "gromacs/gmxana/gstat.h"
53 #include "gromacs/math/units.h"
54 #include "gromacs/math/utilities.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/pbcutil/pbc.h"
57 #include "gromacs/topology/index.h"
58 #include "gromacs/topology/topology.h"
59 #include "gromacs/trajectory/trajectoryframe.h"
60 #include "gromacs/utility/arraysize.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/futil.h"
63 #include "gromacs/utility/gmxassert.h"
64 #include "gromacs/utility/smalloc.h"
66 #define NK 24
67 #define NPK 4
69 #define NKC 6
70 #define NKC0 4
71 static const int kset_c[NKC + 1] = { 0, 3, 9, 13, 16, 19, NK };
73 static rvec v0[NK] = { { 1, 0, 0 }, { 0, 1, 0 }, { 0, 0, 1 }, { 1, 1, 0 }, { 1, -1, 0 },
74 { 1, 0, 1 }, { 1, 0, -1 }, { 0, 1, 1 }, { 0, 1, -1 }, { 1, 1, 1 },
75 { 1, 1, -1 }, { 1, -1, 1 }, { -1, 1, 1 }, { 2, 0, 0 }, { 0, 2, 0 },
76 { 0, 0, 2 }, { 3, 0, 0 }, { 0, 3, 0 }, { 0, 0, 3 }, { 4, 0, 0 },
77 { 0, 4, 0 }, { 0, 0, 4 } };
78 static rvec v1[NK] = { { 0, 1, 0 }, { 0, 0, 1 }, { 1, 0, 0 }, { 0, 0, 1 }, { 0, 0, 1 },
79 { 0, 1, 0 }, { 0, 1, 0 }, { 1, 0, 0 }, { 1, 0, 0 }, { 1, -1, 0 },
80 { 1, -1, 0 }, { 1, 0, -1 }, { 0, 1, -1 }, { 0, 1, 0 }, { 0, 0, 1 },
81 { 1, 0, 0 }, { 0, 1, 0 }, { 0, 0, 1 }, { 1, 0, 0 }, { 0, 1, 0 },
82 { 0, 0, 1 }, { 1, 0, 0 } };
83 static rvec v2[NK] = { { 0, 0, 1 }, { 1, 0, 0 }, { 0, 1, 0 }, { 1, -1, 0 }, { 1, 1, 0 },
84 { 1, 0, -1 }, { 1, 0, 1 }, { 0, 1, -1 }, { 0, 1, 1 }, { 1, 1, -2 },
85 { 1, 1, 2 }, { 1, 2, 1 }, { 2, 1, 1 }, { 0, 0, 1 }, { 1, 0, 0 },
86 { 0, 1, 0 }, { 0, 0, 1 }, { 1, 0, 0 }, { 0, 1, 0 }, { 0, 0, 1 },
87 { 1, 0, 0 }, { 0, 1, 0 } };
89 static void process_tcaf(int nframes,
90 real dt,
91 int nkc,
92 real** tc,
93 rvec* kfac,
94 real rho,
95 real wt,
96 const char* fn_trans,
97 const char* fn_tca,
98 const char* fn_tc,
99 const char* fn_tcf,
100 const char* fn_cub,
101 const char* fn_vk,
102 const gmx_output_env_t* oenv)
104 FILE * fp, *fp_vk, *fp_cub = nullptr;
105 int nk, ntc;
106 real **tcaf, **tcafc = nullptr, eta, *sig;
107 int i, j, k, kc;
108 int ncorr;
109 double fitparms[3];
111 nk = kset_c[nkc];
112 ntc = nk * NPK;
114 if (fn_trans)
116 fp = xvgropen(fn_trans, "Transverse Current", "Time (ps)", "TC (nm/ps)", oenv);
117 for (i = 0; i < nframes; i++)
119 fprintf(fp, "%g", i * dt);
120 for (j = 0; j < ntc; j++)
122 fprintf(fp, " %g", tc[j][i]);
124 fprintf(fp, "\n");
126 xvgrclose(fp);
127 do_view(oenv, fn_trans, "-nxy");
130 ncorr = (nframes + 1) / 2;
131 if (ncorr > gmx::roundToInt(5 * wt / dt))
133 ncorr = gmx::roundToInt(5 * wt / dt) + 1;
135 snew(tcaf, nk);
136 for (k = 0; k < nk; k++)
138 snew(tcaf[k], ncorr);
140 if (fn_cub)
142 snew(tcafc, nkc);
143 for (k = 0; k < nkc; k++)
145 snew(tcafc[k], ncorr);
148 snew(sig, ncorr);
149 for (i = 0; i < ncorr; i++)
151 sig[i] = std::exp(0.5 * i * dt / wt);
154 low_do_autocorr(fn_tca, oenv, "Transverse Current Autocorrelation Functions", nframes, ntc,
155 ncorr, tc, dt, eacNormal, 1, FALSE, FALSE, FALSE, 0, 0, 0);
156 do_view(oenv, fn_tca, "-nxy");
158 fp = xvgropen(fn_tc, "Transverse Current Autocorrelation Functions", "Time (ps)", "TCAF", oenv);
159 for (i = 0; i < ncorr; i++)
161 kc = 0;
162 fprintf(fp, "%g", i * dt);
163 for (k = 0; k < nk; k++)
165 for (j = 0; j < NPK; j++)
167 tcaf[k][i] += tc[NPK * k + j][i];
169 if (fn_cub)
171 for (j = 0; j < NPK; j++)
173 tcafc[kc][i] += tc[NPK * k + j][i];
176 if (i == 0)
178 fprintf(fp, " %g", 1.0);
180 else
182 tcaf[k][i] /= tcaf[k][0];
183 fprintf(fp, " %g", tcaf[k][i]);
185 if (k + 1 == kset_c[kc + 1])
187 kc++;
190 fprintf(fp, "\n");
192 xvgrclose(fp);
193 do_view(oenv, fn_tc, "-nxy");
195 if (fn_cub)
197 fp_cub = xvgropen(fn_cub, "TCAFs and fits", "Time (ps)", "TCAF", oenv);
198 for (kc = 0; kc < nkc; kc++)
200 fprintf(fp_cub, "%g %g\n", 0.0, 1.0);
201 for (i = 1; i < ncorr; i++)
203 tcafc[kc][i] /= tcafc[kc][0];
204 fprintf(fp_cub, "%g %g\n", i * dt, tcafc[kc][i]);
206 fprintf(fp_cub, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
207 tcafc[kc][0] = 1.0;
211 fp_vk = xvgropen(fn_vk, "Fits", "k (nm\\S-1\\N)", "\\8h\\4 (10\\S-3\\N kg m\\S-1\\N s\\S-1\\N)", oenv);
212 if (output_env_get_print_xvgr_codes(oenv))
214 fprintf(fp_vk, "@ s0 symbol 2\n");
215 fprintf(fp_vk, "@ s0 symbol color 1\n");
216 fprintf(fp_vk, "@ s0 linestyle 0\n");
217 if (fn_cub)
219 fprintf(fp_vk, "@ s1 symbol 3\n");
220 fprintf(fp_vk, "@ s1 symbol color 2\n");
223 fp = xvgropen(fn_tcf, "TCAF Fits", "Time (ps)", "", oenv);
224 for (k = 0; k < nk; k++)
226 tcaf[k][0] = 1.0;
227 fitparms[0] = 1;
228 fitparms[1] = 1;
229 do_lmfit(ncorr, tcaf[k], sig, dt, nullptr, 0, ncorr * dt, oenv, bDebugMode(), effnVAC,
230 fitparms, 0, nullptr);
231 eta = 1000 * fitparms[1] * rho / (4 * fitparms[0] * PICO * norm2(kfac[k]) / (NANO * NANO));
232 fprintf(stdout, "k %6.3f tau %6.3f eta %8.5f 10^-3 kg/(m s)\n", norm(kfac[k]), fitparms[0], eta);
233 fprintf(fp_vk, "%6.3f %g\n", norm(kfac[k]), eta);
234 for (i = 0; i < ncorr; i++)
236 fprintf(fp, "%g %g\n", i * dt, fit_function(effnVAC, fitparms, i * dt));
238 fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
240 xvgrclose(fp);
241 do_view(oenv, fn_tcf, "-nxy");
243 if (fn_cub)
245 fprintf(stdout, "Averaged over k-vectors:\n");
246 fprintf(fp_vk, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
247 for (k = 0; k < nkc; k++)
249 tcafc[k][0] = 1.0;
250 fitparms[0] = 1;
251 fitparms[1] = 1;
252 do_lmfit(ncorr, tcafc[k], sig, dt, nullptr, 0, ncorr * dt, oenv, bDebugMode(), effnVAC,
253 fitparms, 0, nullptr);
254 eta = 1000 * fitparms[1] * rho
255 / (4 * fitparms[0] * PICO * norm2(kfac[kset_c[k]]) / (NANO * NANO));
256 fprintf(stdout, "k %6.3f tau %6.3f Omega %6.3f eta %8.5f 10^-3 kg/(m s)\n",
257 norm(kfac[kset_c[k]]), fitparms[0], fitparms[1], eta);
258 fprintf(fp_vk, "%6.3f %g\n", norm(kfac[kset_c[k]]), eta);
259 for (i = 0; i < ncorr; i++)
261 fprintf(fp_cub, "%g %g\n", i * dt, fit_function(effnVAC, fitparms, i * dt));
263 fprintf(fp_cub, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
265 fprintf(fp_vk, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
266 xvgrclose(fp_cub);
267 do_view(oenv, fn_cub, "-nxy");
269 xvgrclose(fp_vk);
270 do_view(oenv, fn_vk, "-nxy");
274 int gmx_tcaf(int argc, char* argv[])
276 const char* desc[] = {
277 "[THISMODULE] computes tranverse current autocorrelations.",
278 "These are used to estimate the shear viscosity, [GRK]eta[grk].",
279 "For details see: Palmer, Phys. Rev. E 49 (1994) pp 359-366.[PAR]",
280 "Transverse currents are calculated using the",
281 "k-vectors (1,0,0) and (2,0,0) each also in the [IT]y[it]- and [IT]z[it]-direction,",
282 "(1,1,0) and (1,-1,0) each also in the 2 other planes (these vectors",
283 "are not independent) and (1,1,1) and the 3 other box diagonals (also",
284 "not independent). For each k-vector the sine and cosine are used, in",
285 "combination with the velocity in 2 perpendicular directions. This gives",
286 "a total of 16*2*2=64 transverse currents. One autocorrelation is",
287 "calculated fitted for each k-vector, which gives 16 TCAFs. Each of",
288 "these TCAFs is fitted to [MATH]f(t) = [EXP]-v[exp]([COSH]Wv[cosh] + 1/W ",
289 "[SINH]Wv[sinh])[math],",
290 "[MATH]v = -t/(2 [GRK]tau[grk])[math], [MATH]W = [SQRT]1 - 4 [GRK]tau[grk] ",
291 "[GRK]eta[grk]/[GRK]rho[grk] k^2[sqrt][math], which gives 16 values of [GRK]tau[grk]",
292 "and [GRK]eta[grk]. The fit weights decay exponentially with time constant [MATH]w[math] ",
293 "(given with [TT]-wt[tt]) as [MATH][EXP]-t/w[exp][math], and the TCAF and",
294 "fit are calculated up to time [MATH]5*w[math].",
295 "The [GRK]eta[grk] values should be fitted to [MATH]1 - a [GRK]eta[grk](k) k^2[math], ",
296 "from which one can estimate the shear viscosity at k=0.[PAR]",
297 "When the box is cubic, one can use the option [TT]-oc[tt], which",
298 "averages the TCAFs over all k-vectors with the same length.",
299 "This results in more accurate TCAFs.",
300 "Both the cubic TCAFs and fits are written to [TT]-oc[tt]",
301 "The cubic [GRK]eta[grk] estimates are also written to [TT]-ov[tt].[PAR]",
302 "With option [TT]-mol[tt], the transverse current is determined of",
303 "molecules instead of atoms. In this case, the index group should",
304 "consist of molecule numbers instead of atom numbers.[PAR]",
305 "The k-dependent viscosities in the [TT]-ov[tt] file should be",
306 "fitted to [MATH][GRK]eta[grk](k) = [GRK]eta[grk][SUB]0[sub] (1 - a k^2)[math] to obtain ",
307 "the viscosity at",
308 "infinite wavelength.[PAR]",
309 "[BB]Note:[bb] make sure you write coordinates and velocities often enough.",
310 "The initial, non-exponential, part of the autocorrelation function",
311 "is very important for obtaining a good fit."
314 static gmx_bool bMol = FALSE, bK34 = FALSE;
315 static real wt = 5;
316 t_pargs pa[] = {
317 { "-mol", FALSE, etBOOL, { &bMol }, "Calculate TCAF of molecules" },
318 { "-k34", FALSE, etBOOL, { &bK34 }, "Also use k=(3,0,0) and k=(4,0,0)" },
319 { "-wt", FALSE, etREAL, { &wt }, "Exponential decay time for the TCAF fit weights" }
322 t_topology top;
323 PbcType pbcType;
324 t_trxframe fr;
325 matrix box;
326 gmx_bool bTop;
327 int gnx;
328 int * index, *atndx = nullptr, at;
329 char* grpname;
330 char title[256];
331 real t0, t1, dt, m, mtot, sysmass, rho, sx, cx;
332 t_trxstatus* status;
333 int nframes, n_alloc, i, j, k, d;
334 rvec mv_mol, cm_mol, kfac[NK];
335 int nkc, nk, ntc;
336 real** tc;
337 gmx_output_env_t* oenv;
339 #define NHISTO 360
341 t_filenm fnm[] = { { efTRN, "-f", nullptr, ffREAD }, { efTPS, nullptr, nullptr, ffOPTRD },
342 { efNDX, nullptr, nullptr, ffOPTRD }, { efXVG, "-ot", "transcur", ffOPTWR },
343 { efXVG, "-oa", "tcaf_all", ffWRITE }, { efXVG, "-o", "tcaf", ffWRITE },
344 { efXVG, "-of", "tcaf_fit", ffWRITE }, { efXVG, "-oc", "tcaf_cub", ffOPTWR },
345 { efXVG, "-ov", "visc_k", ffWRITE } };
346 #define NFILE asize(fnm)
347 int npargs;
348 t_pargs* ppa;
350 npargs = asize(pa);
351 ppa = add_acf_pargs(&npargs, pa);
353 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME, NFILE, fnm, npargs, ppa,
354 asize(desc), desc, 0, nullptr, &oenv))
356 sfree(ppa);
357 return 0;
360 bTop = read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &pbcType, nullptr, nullptr, box, TRUE);
361 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
363 if (bMol)
365 if (!bTop)
367 gmx_fatal(FARGS, "Need a topology to determine the molecules");
369 atndx = top.mols.index;
372 if (bK34)
374 nkc = NKC;
376 else
378 nkc = NKC0;
380 nk = kset_c[nkc];
381 GMX_ASSERT(nk >= 16, "Has to be over 16 because nkc is either NKC or NKC0.");
382 ntc = nk * NPK;
384 sprintf(title, "Velocity Autocorrelation Function for %s", grpname);
386 sysmass = 0;
387 for (i = 0; i < nk; i++)
389 if (iprod(v0[i], v1[i]) != 0)
391 gmx_fatal(FARGS, "DEATH HORROR: vectors not orthogonal");
393 if (iprod(v0[i], v2[i]) != 0)
395 gmx_fatal(FARGS, "DEATH HORROR: vectors not orthogonal");
397 if (iprod(v1[i], v2[i]) != 0)
399 gmx_fatal(FARGS, "DEATH HORROR: vectors not orthogonal");
401 unitv(v1[i], v1[i]);
402 unitv(v2[i], v2[i]);
404 snew(tc, ntc);
405 for (i = 0; i < top.atoms.nr; i++)
407 sysmass += top.atoms.atom[i].m;
410 read_first_frame(oenv, &status, ftp2fn(efTRN, NFILE, fnm), &fr, TRX_NEED_X | TRX_NEED_V);
411 t0 = fr.time;
413 n_alloc = 0;
414 nframes = 0;
415 rho = 0;
420 if (nframes >= n_alloc)
422 n_alloc += 100;
423 for (i = 0; i < ntc; i++)
425 srenew(tc[i], n_alloc);
429 rho += 1 / det(fr.box);
430 for (k = 0; k < nk; k++)
432 for (d = 0; d < DIM; d++)
434 kfac[k][d] = 2 * M_PI * v0[k][d] / fr.box[d][d];
437 for (i = 0; i < ntc; i++)
439 tc[i][nframes] = 0;
442 for (i = 0; i < gnx; i++)
444 if (bMol)
446 clear_rvec(mv_mol);
447 clear_rvec(cm_mol);
448 mtot = 0;
449 for (j = 0; j < atndx[index[i] + 1] - atndx[index[i]]; j++)
451 at = atndx[index[i]] + j;
452 m = top.atoms.atom[at].m;
453 mv_mol[XX] += m * fr.v[at][XX];
454 mv_mol[YY] += m * fr.v[at][YY];
455 mv_mol[ZZ] += m * fr.v[at][ZZ];
456 cm_mol[XX] += m * fr.x[at][XX];
457 cm_mol[YY] += m * fr.x[at][YY];
458 cm_mol[ZZ] += m * fr.x[at][ZZ];
459 mtot += m;
461 svmul(1.0 / mtot, cm_mol, cm_mol);
463 else
465 svmul(top.atoms.atom[index[i]].m, fr.v[index[i]], mv_mol);
468 if (!bMol)
470 copy_rvec(fr.x[index[i]], cm_mol);
472 j = 0;
473 for (k = 0; k < nk; k++)
475 sx = std::sin(iprod(kfac[k], cm_mol));
476 cx = std::cos(iprod(kfac[k], cm_mol));
477 tc[j][nframes] += sx * iprod(v1[k], mv_mol);
478 j++;
479 tc[j][nframes] += cx * iprod(v1[k], mv_mol);
480 j++;
481 tc[j][nframes] += sx * iprod(v2[k], mv_mol);
482 j++;
483 tc[j][nframes] += cx * iprod(v2[k], mv_mol);
484 j++;
488 t1 = fr.time;
489 nframes++;
490 } while (read_next_frame(oenv, status, &fr));
491 close_trx(status);
493 dt = (t1 - t0) / (nframes - 1);
495 rho *= sysmass / nframes * AMU / (NANO * NANO * NANO);
496 fprintf(stdout, "Density = %g (kg/m^3)\n", rho);
497 process_tcaf(nframes, dt, nkc, tc, kfac, rho, wt, opt2fn_null("-ot", NFILE, fnm),
498 opt2fn("-oa", NFILE, fnm), opt2fn("-o", NFILE, fnm), opt2fn("-of", NFILE, fnm),
499 opt2fn_null("-oc", NFILE, fnm), opt2fn("-ov", NFILE, fnm), oenv);
501 return 0;