Add replacements for pbc enumerations
[gromacs.git] / src / gromacs / gmxana / gmx_tcaf.cpp
blob7e1553e18e45cf397bbced9dc9c4e99b2897d937
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,2018, 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.
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.
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.
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.
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.
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.
37 #include "gmxpre.h"
39 #include <cmath>
40 #include <cstdio>
41 #include <cstring>
43 #include "gromacs/commandline/pargs.h"
44 #include "gromacs/commandline/viewit.h"
45 #include "gromacs/correlationfunctions/autocorr.h"
46 #include "gromacs/correlationfunctions/expfit.h"
47 #include "gromacs/fileio/confio.h"
48 #include "gromacs/fileio/trxio.h"
49 #include "gromacs/fileio/xvgr.h"
50 #include "gromacs/gmxana/gmx_ana.h"
51 #include "gromacs/gmxana/gstat.h"
52 #include "gromacs/math/units.h"
53 #include "gromacs/math/utilities.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/pbcutil/pbc.h"
56 #include "gromacs/topology/index.h"
57 #include "gromacs/topology/topology.h"
58 #include "gromacs/trajectory/trajectoryframe.h"
59 #include "gromacs/utility/arraysize.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/futil.h"
62 #include "gromacs/utility/gmxassert.h"
63 #include "gromacs/utility/smalloc.h"
65 #define NK 24
66 #define NPK 4
68 #define NKC 6
69 #define NKC0 4
70 static const int kset_c[NKC+1] = { 0, 3, 9, 13, 16, 19, NK };
72 static 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 static 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 static 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}};
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 gmx_output_env_t *oenv)
82 FILE *fp, *fp_vk, *fp_cub = nullptr;
83 int nk, ntc;
84 real **tcaf, **tcafc = nullptr, eta, *sig;
85 int i, j, k, kc;
86 int ncorr;
87 double fitparms[3];
89 nk = kset_c[nkc];
90 ntc = nk*NPK;
92 if (fn_trans)
94 fp = xvgropen(fn_trans, "Transverse Current", "Time (ps)", "TC (nm/ps)",
95 oenv);
96 for (i = 0; i < nframes; i++)
98 fprintf(fp, "%g", i*dt);
99 for (j = 0; j < ntc; j++)
101 fprintf(fp, " %g", tc[j][i]);
103 fprintf(fp, "\n");
105 xvgrclose(fp);
106 do_view(oenv, fn_trans, "-nxy");
109 ncorr = (nframes+1)/2;
110 if (ncorr > gmx::roundToInt(5*wt/dt))
112 ncorr = gmx::roundToInt(5*wt/dt)+1;
114 snew(tcaf, nk);
115 for (k = 0; k < nk; k++)
117 snew(tcaf[k], ncorr);
119 if (fn_cub)
121 snew(tcafc, nkc);
122 for (k = 0; k < nkc; k++)
124 snew(tcafc[k], ncorr);
127 snew(sig, ncorr);
128 for (i = 0; i < ncorr; i++)
130 sig[i] = std::exp(0.5*i*dt/wt);
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");
138 fp = xvgropen(fn_tc, "Transverse Current Autocorrelation Functions",
139 "Time (ps)", "TCAF", oenv);
140 for (i = 0; i < ncorr; i++)
142 kc = 0;
143 fprintf(fp, "%g", i*dt);
144 for (k = 0; k < nk; k++)
146 for (j = 0; j < NPK; j++)
148 tcaf[k][i] += tc[NPK*k+j][i];
150 if (fn_cub)
152 for (j = 0; j < NPK; j++)
154 tcafc[kc][i] += tc[NPK*k+j][i];
157 if (i == 0)
159 fprintf(fp, " %g", 1.0);
161 else
163 tcaf[k][i] /= tcaf[k][0];
164 fprintf(fp, " %g", tcaf[k][i]);
166 if (k+1 == kset_c[kc+1])
168 kc++;
171 fprintf(fp, "\n");
173 xvgrclose(fp);
174 do_view(oenv, fn_tc, "-nxy");
176 if (fn_cub)
178 fp_cub = xvgropen(fn_cub, "TCAFs and fits", "Time (ps)", "TCAF", oenv);
179 for (kc = 0; kc < nkc; kc++)
181 fprintf(fp_cub, "%g %g\n", 0.0, 1.0);
182 for (i = 1; i < ncorr; i++)
184 tcafc[kc][i] /= tcafc[kc][0];
185 fprintf(fp_cub, "%g %g\n", i*dt, tcafc[kc][i]);
187 fprintf(fp_cub, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
188 tcafc[kc][0] = 1.0;
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 if (output_env_get_print_xvgr_codes(oenv))
196 fprintf(fp_vk, "@ s0 symbol 2\n");
197 fprintf(fp_vk, "@ s0 symbol color 1\n");
198 fprintf(fp_vk, "@ s0 linestyle 0\n");
199 if (fn_cub)
201 fprintf(fp_vk, "@ s1 symbol 3\n");
202 fprintf(fp_vk, "@ s1 symbol color 2\n");
205 fp = xvgropen(fn_tcf, "TCAF Fits", "Time (ps)", "", oenv);
206 for (k = 0; k < nk; k++)
208 tcaf[k][0] = 1.0;
209 fitparms[0] = 1;
210 fitparms[1] = 1;
211 do_lmfit(ncorr, tcaf[k], sig, dt, nullptr, 0, ncorr*dt,
212 oenv, bDebugMode(), effnVAC, fitparms, 0, nullptr);
213 eta = 1000*fitparms[1]*rho/
214 (4*fitparms[0]*PICO*norm2(kfac[k])/(NANO*NANO));
215 fprintf(stdout, "k %6.3f tau %6.3f eta %8.5f 10^-3 kg/(m s)\n",
216 norm(kfac[k]), fitparms[0], eta);
217 fprintf(fp_vk, "%6.3f %g\n", norm(kfac[k]), eta);
218 for (i = 0; i < ncorr; i++)
220 fprintf(fp, "%g %g\n", i*dt, fit_function(effnVAC, fitparms, i*dt));
222 fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
224 xvgrclose(fp);
225 do_view(oenv, fn_tcf, "-nxy");
227 if (fn_cub)
229 fprintf(stdout, "Averaged over k-vectors:\n");
230 fprintf(fp_vk, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
231 for (k = 0; k < nkc; k++)
233 tcafc[k][0] = 1.0;
234 fitparms[0] = 1;
235 fitparms[1] = 1;
236 do_lmfit(ncorr, tcafc[k], sig, dt, nullptr, 0, ncorr*dt,
237 oenv, bDebugMode(), effnVAC, fitparms, 0, nullptr);
238 eta = 1000*fitparms[1]*rho/
239 (4*fitparms[0]*PICO*norm2(kfac[kset_c[k]])/(NANO*NANO));
240 fprintf(stdout,
241 "k %6.3f tau %6.3f Omega %6.3f eta %8.5f 10^-3 kg/(m s)\n",
242 norm(kfac[kset_c[k]]), fitparms[0], fitparms[1], eta);
243 fprintf(fp_vk, "%6.3f %g\n", norm(kfac[kset_c[k]]), eta);
244 for (i = 0; i < ncorr; i++)
246 fprintf(fp_cub, "%g %g\n", i*dt, fit_function(effnVAC, fitparms, i*dt));
248 fprintf(fp_cub, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
250 fprintf(fp_vk, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
251 xvgrclose(fp_cub);
252 do_view(oenv, fn_cub, "-nxy");
254 xvgrclose(fp_vk);
255 do_view(oenv, fn_vk, "-nxy");
259 int gmx_tcaf(int argc, char *argv[])
261 const char *desc[] = {
262 "[THISMODULE] computes tranverse current autocorrelations.",
263 "These are used to estimate the shear viscosity, [GRK]eta[grk].",
264 "For details see: Palmer, Phys. Rev. E 49 (1994) pp 359-366.[PAR]",
265 "Transverse currents are calculated using the",
266 "k-vectors (1,0,0) and (2,0,0) each also in the [IT]y[it]- and [IT]z[it]-direction,",
267 "(1,1,0) and (1,-1,0) each also in the 2 other planes (these vectors",
268 "are not independent) and (1,1,1) and the 3 other box diagonals (also",
269 "not independent). For each k-vector the sine and cosine are used, in",
270 "combination with the velocity in 2 perpendicular directions. This gives",
271 "a total of 16*2*2=64 transverse currents. One autocorrelation is",
272 "calculated fitted for each k-vector, which gives 16 TCAFs. Each of",
273 "these TCAFs is fitted to [MATH]f(t) = [EXP]-v[exp]([COSH]Wv[cosh] + 1/W [SINH]Wv[sinh])[math],",
274 "[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]",
275 "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",
276 "fit are calculated up to time [MATH]5*w[math].",
277 "The [GRK]eta[grk] values should be fitted to [MATH]1 - a [GRK]eta[grk](k) k^2[math], from which",
278 "one can estimate the shear viscosity at k=0.[PAR]",
279 "When the box is cubic, one can use the option [TT]-oc[tt], which",
280 "averages the TCAFs over all k-vectors with the same length.",
281 "This results in more accurate TCAFs.",
282 "Both the cubic TCAFs and fits are written to [TT]-oc[tt]",
283 "The cubic [GRK]eta[grk] estimates are also written to [TT]-ov[tt].[PAR]",
284 "With option [TT]-mol[tt], the transverse current is determined of",
285 "molecules instead of atoms. In this case, the index group should",
286 "consist of molecule numbers instead of atom numbers.[PAR]",
287 "The k-dependent viscosities in the [TT]-ov[tt] file should be",
288 "fitted to [MATH][GRK]eta[grk](k) = [GRK]eta[grk][SUB]0[sub] (1 - a k^2)[math] to obtain the viscosity at",
289 "infinite wavelength.[PAR]",
290 "[BB]Note:[bb] make sure you write coordinates and velocities often enough.",
291 "The initial, non-exponential, part of the autocorrelation function",
292 "is very important for obtaining a good fit."
295 static gmx_bool bMol = FALSE, bK34 = FALSE;
296 static real wt = 5;
297 t_pargs pa[] = {
298 { "-mol", FALSE, etBOOL, {&bMol},
299 "Calculate TCAF of molecules" },
300 { "-k34", FALSE, etBOOL, {&bK34},
301 "Also use k=(3,0,0) and k=(4,0,0)" },
302 { "-wt", FALSE, etREAL, {&wt},
303 "Exponential decay time for the TCAF fit weights" }
306 t_topology top;
307 int ePBC;
308 t_trxframe fr;
309 matrix box;
310 gmx_bool bTop;
311 int gnx;
312 int *index, *atndx = nullptr, at;
313 char *grpname;
314 char title[256];
315 real t0, t1, dt, m, mtot, sysmass, rho, sx, cx;
316 t_trxstatus *status;
317 int nframes, n_alloc, i, j, k, d;
318 rvec mv_mol, cm_mol, kfac[NK];
319 int nkc, nk, ntc;
320 real **tc;
321 gmx_output_env_t *oenv;
323 #define NHISTO 360
325 t_filenm fnm[] = {
326 { efTRN, "-f", nullptr, ffREAD },
327 { efTPS, nullptr, nullptr, ffOPTRD },
328 { efNDX, nullptr, nullptr, ffOPTRD },
329 { efXVG, "-ot", "transcur", ffOPTWR },
330 { efXVG, "-oa", "tcaf_all", ffWRITE },
331 { efXVG, "-o", "tcaf", ffWRITE },
332 { efXVG, "-of", "tcaf_fit", ffWRITE },
333 { efXVG, "-oc", "tcaf_cub", ffOPTWR },
334 { efXVG, "-ov", "visc_k", ffWRITE }
336 #define NFILE asize(fnm)
337 int npargs;
338 t_pargs *ppa;
340 npargs = asize(pa);
341 ppa = add_acf_pargs(&npargs, pa);
343 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME,
344 NFILE, fnm, npargs, ppa, asize(desc), desc, 0, nullptr, &oenv))
346 sfree(ppa);
347 return 0;
350 bTop = read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &ePBC, nullptr, nullptr, box,
351 TRUE);
352 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
354 if (bMol)
356 if (!bTop)
358 gmx_fatal(FARGS, "Need a topology to determine the molecules");
360 atndx = top.mols.index;
363 if (bK34)
365 nkc = NKC;
367 else
369 nkc = NKC0;
371 nk = kset_c[nkc];
372 GMX_ASSERT(nk >= 16, "Has to be over 16 because nkc is either NKC or NKC0.");
373 ntc = nk*NPK;
375 sprintf(title, "Velocity Autocorrelation Function for %s", grpname);
377 sysmass = 0;
378 for (i = 0; i < nk; i++)
380 if (iprod(v0[i], v1[i]) != 0)
382 gmx_fatal(FARGS, "DEATH HORROR: vectors not orthogonal");
384 if (iprod(v0[i], v2[i]) != 0)
386 gmx_fatal(FARGS, "DEATH HORROR: vectors not orthogonal");
388 if (iprod(v1[i], v2[i]) != 0)
390 gmx_fatal(FARGS, "DEATH HORROR: vectors not orthogonal");
392 unitv(v1[i], v1[i]);
393 unitv(v2[i], v2[i]);
395 snew(tc, ntc);
396 for (i = 0; i < top.atoms.nr; i++)
398 sysmass += top.atoms.atom[i].m;
401 read_first_frame(oenv, &status, ftp2fn(efTRN, NFILE, fnm), &fr,
402 TRX_NEED_X | TRX_NEED_V);
403 t0 = fr.time;
405 n_alloc = 0;
406 nframes = 0;
407 rho = 0;
412 if (nframes >= n_alloc)
414 n_alloc += 100;
415 for (i = 0; i < ntc; i++)
417 srenew(tc[i], n_alloc);
421 rho += 1/det(fr.box);
422 for (k = 0; k < nk; k++)
424 for (d = 0; d < DIM; d++)
426 kfac[k][d] = 2*M_PI*v0[k][d]/fr.box[d][d];
429 for (i = 0; i < ntc; i++)
431 tc[i][nframes] = 0;
434 for (i = 0; i < gnx; i++)
436 if (bMol)
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++)
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;
453 svmul(1.0/mtot, cm_mol, cm_mol);
455 else
457 svmul(top.atoms.atom[index[i]].m, fr.v[index[i]], mv_mol);
460 if (!bMol)
462 copy_rvec(fr.x[index[i]], cm_mol);
464 j = 0;
465 for (k = 0; k < nk; k++)
467 sx = std::sin(iprod(kfac[k], cm_mol));
468 cx = std::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++;
480 t1 = fr.time;
481 nframes++;
483 while (read_next_frame(oenv, status, &fr));
484 close_trx(status);
486 dt = (t1-t0)/(nframes-1);
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);
496 return 0;