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, 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.
43 #include "gromacs/fileio/confio.h"
44 #include "gromacs/utility/fatalerror.h"
45 #include "gromacs/utility/futil.h"
48 #include "gromacs/math/utilities.h"
49 #include "gromacs/math/units.h"
50 #include "gromacs/topology/index.h"
51 #include "gromacs/utility/smalloc.h"
52 #include "gromacs/commandline/pargs.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/fileio/xvgr.h"
58 #include "gromacs/pbcutil/pbc.h"
60 #include "gromacs/fileio/trxio.h"
68 int kset_c
[NKC
+1] = { 0, 3, 9, 13, 16, 19, NK
};
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}};
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
)
80 FILE *fp
, *fp_vk
, *fp_cub
= NULL
;
82 real
**tcaf
, **tcafc
= NULL
, eta
;
85 real fitparms
[3], *sig
;
92 fp
= xvgropen(fn_trans
, "Transverse Current", "Time (ps)", "TC (nm/ps)",
94 for (i
= 0; i
< nframes
; i
++)
96 fprintf(fp
, "%g", i
*dt
);
97 for (j
= 0; j
< ntc
; j
++)
99 fprintf(fp
, " %g", tc
[j
][i
]);
104 do_view(oenv
, fn_trans
, "-nxy");
107 ncorr
= (nframes
+1)/2;
108 if (ncorr
> (int)(5*wt
/dt
+0.5))
110 ncorr
= (int)(5*wt
/dt
+0.5)+1;
113 for (k
= 0; k
< nk
; k
++)
115 snew(tcaf
[k
], ncorr
);
120 for (k
= 0; k
< nkc
; k
++)
122 snew(tcafc
[k
], ncorr
);
126 for (i
= 0; i
< ncorr
; i
++)
128 sig
[i
] = exp(0.5*i
*dt
/wt
);
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);
134 do_view(oenv
, fn_tca
, "-nxy");
136 fp
= xvgropen(fn_tc
, "Transverse Current Autocorrelation Functions",
137 "Time (ps)", "TCAF", oenv
);
138 for (i
= 0; i
< ncorr
; i
++)
141 fprintf(fp
, "%g", i
*dt
);
142 for (k
= 0; k
< nk
; k
++)
144 for (j
= 0; j
< NPK
; j
++)
146 tcaf
[k
][i
] += tc
[NPK
*k
+j
][i
];
150 for (j
= 0; j
< NPK
; j
++)
152 tcafc
[kc
][i
] += tc
[NPK
*k
+j
][i
];
157 fprintf(fp
, " %g", 1.0);
161 tcaf
[k
][i
] /= tcaf
[k
][0];
162 fprintf(fp
, " %g", tcaf
[k
][i
]);
164 if (k
+1 == kset_c
[kc
+1])
172 do_view(oenv
, fn_tc
, "-nxy");
176 fp_cub
= xvgropen(fn_cub
, "TCAFs and fits", "Time (ps)", "TCAF", oenv
);
177 for (kc
= 0; kc
< nkc
; kc
++)
179 fprintf(fp_cub
, "%g %g\n", 0.0, 1.0);
180 for (i
= 1; i
< ncorr
; i
++)
182 tcafc
[kc
][i
] /= tcafc
[kc
][0];
183 fprintf(fp_cub
, "%g %g\n", i
*dt
, tcafc
[kc
][i
]);
185 fprintf(fp_cub
, "%s\n", output_env_get_print_xvgr_codes(oenv
) ? "&" : "");
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 if (output_env_get_print_xvgr_codes(oenv
))
194 fprintf(fp_vk
, "@ s0 symbol 2\n");
195 fprintf(fp_vk
, "@ s0 symbol color 1\n");
196 fprintf(fp_vk
, "@ s0 linestyle 0\n");
199 fprintf(fp_vk
, "@ s1 symbol 3\n");
200 fprintf(fp_vk
, "@ s1 symbol color 2\n");
203 fp
= xvgropen(fn_tcf
, "TCAF Fits", "Time (ps)", "", oenv
);
204 for (k
= 0; k
< nk
; k
++)
209 do_lmfit(ncorr
, tcaf
[k
], sig
, dt
, 0, 0, ncorr
*dt
,
210 oenv
, bDebugMode(), effnVAC
, fitparms
, 0);
211 eta
= 1000*fitparms
[1]*rho
/
212 (4*fitparms
[0]*PICO
*norm2(kfac
[k
])/(NANO
*NANO
));
213 fprintf(stdout
, "k %6.3f tau %6.3f eta %8.5f 10^-3 kg/(m s)\n",
214 norm(kfac
[k
]), fitparms
[0], eta
);
215 fprintf(fp_vk
, "%6.3f %g\n", norm(kfac
[k
]), eta
);
216 for (i
= 0; i
< ncorr
; i
++)
218 fprintf(fp
, "%g %g\n", i
*dt
, fit_function(effnVAC
, fitparms
, i
*dt
));
220 fprintf(fp
, "%s\n", output_env_get_print_xvgr_codes(oenv
) ? "&" : "");
223 do_view(oenv
, fn_tcf
, "-nxy");
227 fprintf(stdout
, "Averaged over k-vectors:\n");
228 fprintf(fp_vk
, "%s\n", output_env_get_print_xvgr_codes(oenv
) ? "&" : "");
229 for (k
= 0; k
< nkc
; k
++)
234 do_lmfit(ncorr
, tcafc
[k
], sig
, dt
, 0, 0, ncorr
*dt
,
235 oenv
, bDebugMode(), effnVAC
, fitparms
, 0);
236 eta
= 1000*fitparms
[1]*rho
/
237 (4*fitparms
[0]*PICO
*norm2(kfac
[kset_c
[k
]])/(NANO
*NANO
));
239 "k %6.3f tau %6.3f Omega %6.3f eta %8.5f 10^-3 kg/(m s)\n",
240 norm(kfac
[kset_c
[k
]]), fitparms
[0], fitparms
[1], eta
);
241 fprintf(fp_vk
, "%6.3f %g\n", norm(kfac
[kset_c
[k
]]), eta
);
242 for (i
= 0; i
< ncorr
; i
++)
244 fprintf(fp_cub
, "%g %g\n", i
*dt
, fit_function(effnVAC
, fitparms
, i
*dt
));
246 fprintf(fp_cub
, "%s\n", output_env_get_print_xvgr_codes(oenv
) ? "&" : "");
248 fprintf(fp_vk
, "%s\n", output_env_get_print_xvgr_codes(oenv
) ? "&" : "");
250 do_view(oenv
, fn_cub
, "-nxy");
253 do_view(oenv
, fn_vk
, "-nxy");
257 int gmx_tcaf(int argc
, char *argv
[])
259 const char *desc
[] = {
260 "[THISMODULE] computes tranverse current autocorrelations.",
261 "These are used to estimate the shear viscosity, [GRK]eta[grk].",
262 "For details see: Palmer, Phys. Rev. E 49 (1994) pp 359-366.[PAR]",
263 "Transverse currents are calculated using the",
264 "k-vectors (1,0,0) and (2,0,0) each also in the [IT]y[it]- and [IT]z[it]-direction,",
265 "(1,1,0) and (1,-1,0) each also in the 2 other planes (these vectors",
266 "are not independent) and (1,1,1) and the 3 other box diagonals (also",
267 "not independent). For each k-vector the sine and cosine are used, in",
268 "combination with the velocity in 2 perpendicular directions. This gives",
269 "a total of 16*2*2=64 transverse currents. One autocorrelation is",
270 "calculated fitted for each k-vector, which gives 16 TCAFs. Each of",
271 "these TCAFs is fitted to [MATH]f(t) = [EXP]-v[exp]([COSH]Wv[cosh] + 1/W [SINH]Wv[sinh])[math],",
272 "[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]",
273 "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",
274 "fit are calculated up to time [MATH]5*w[math].",
275 "The [GRK]eta[grk] values should be fitted to [MATH]1 - a [GRK]eta[grk](k) k^2[math], from which",
276 "one can estimate the shear viscosity at k=0.[PAR]",
277 "When the box is cubic, one can use the option [TT]-oc[tt], which",
278 "averages the TCAFs over all k-vectors with the same length.",
279 "This results in more accurate TCAFs.",
280 "Both the cubic TCAFs and fits are written to [TT]-oc[tt]",
281 "The cubic [GRK]eta[grk] estimates are also written to [TT]-ov[tt].[PAR]",
282 "With option [TT]-mol[tt], the transverse current is determined of",
283 "molecules instead of atoms. In this case, the index group should",
284 "consist of molecule numbers instead of atom numbers.[PAR]",
285 "The k-dependent viscosities in the [TT]-ov[tt] file should be",
286 "fitted to [MATH][GRK]eta[grk](k) = [GRK]eta[grk][SUB]0[sub] (1 - a k^2)[math] to obtain the viscosity at",
287 "infinite wavelength.[PAR]",
288 "[BB]Note:[bb] make sure you write coordinates and velocities often enough.",
289 "The initial, non-exponential, part of the autocorrelation function",
290 "is very important for obtaining a good fit."
293 static gmx_bool bMol
= FALSE
, bK34
= FALSE
;
296 { "-mol", FALSE
, etBOOL
, {&bMol
},
297 "Calculate TCAF of molecules" },
298 { "-k34", FALSE
, etBOOL
, {&bK34
},
299 "Also use k=(3,0,0) and k=(4,0,0)" },
300 { "-wt", FALSE
, etREAL
, {&wt
},
301 "Exponential decay time for the TCAF fit weights" }
308 gmx_bool bTPS
, bTop
; /* ,bCubic; */
310 atom_id
*index
, *atndx
= NULL
, at
;
313 real t0
, t1
, dt
, m
, mtot
, sysmass
, rho
, sx
, cx
;
315 int nframes
, n_alloc
, i
, j
, k
, d
;
316 rvec mv_mol
, cm_mol
, kfac
[NK
];
324 { efTRN
, "-f", NULL
, ffREAD
},
325 { efTPS
, NULL
, NULL
, ffOPTRD
},
326 { efNDX
, NULL
, NULL
, ffOPTRD
},
327 { efXVG
, "-ot", "transcur", ffOPTWR
},
328 { efXVG
, "-oa", "tcaf_all", ffWRITE
},
329 { efXVG
, "-o", "tcaf", ffWRITE
},
330 { efXVG
, "-of", "tcaf_fit", ffWRITE
},
331 { efXVG
, "-oc", "tcaf_cub", ffOPTWR
},
332 { efXVG
, "-ov", "visc_k", ffWRITE
}
334 #define NFILE asize(fnm)
339 ppa
= add_acf_pargs(&npargs
, pa
);
341 if (!parse_common_args(&argc
, argv
, PCA_CAN_VIEW
| PCA_CAN_TIME
| PCA_BE_NICE
,
342 NFILE
, fnm
, npargs
, ppa
, asize(desc
), desc
, 0, NULL
, &oenv
))
347 bTop
= read_tps_conf(ftp2fn(efTPS
, NFILE
, fnm
), title
, &top
, &ePBC
, NULL
, NULL
, box
,
349 get_index(&top
.atoms
, ftp2fn_null(efNDX
, NFILE
, fnm
), 1, &gnx
, &index
, &grpname
);
355 gmx_fatal(FARGS
, "Need a topology to determine the molecules");
357 atndx
= top
.mols
.index
;
371 sprintf(title
, "Velocity Autocorrelation Function for %s", grpname
);
374 for (i
= 0; i
< nk
; i
++)
376 if (iprod(v0
[i
], v1
[i
]) != 0)
378 gmx_fatal(FARGS
, "DEATH HORROR: vectors not orthogonal");
380 if (iprod(v0
[i
], v2
[i
]) != 0)
382 gmx_fatal(FARGS
, "DEATH HORROR: vectors not orthogonal");
384 if (iprod(v1
[i
], v2
[i
]) != 0)
386 gmx_fatal(FARGS
, "DEATH HORROR: vectors not orthogonal");
392 for (i
= 0; i
< top
.atoms
.nr
; i
++)
394 sysmass
+= top
.atoms
.atom
[i
].m
;
397 read_first_frame(oenv
, &status
, ftp2fn(efTRN
, NFILE
, fnm
), &fr
,
398 TRX_NEED_X
| TRX_NEED_V
);
408 bCubic = bCubic && !TRICLINIC(fr.box) &&
409 fabs(fr.box[XX][XX]-fr.box[YY][YY]) < 0.001*fr.box[XX][XX] &&
410 fabs(fr.box[XX][XX]-fr.box[ZZ][ZZ]) < 0.001*fr.box[XX][XX];
413 if (nframes
>= n_alloc
)
416 for (i
= 0; i
< ntc
; i
++)
418 srenew(tc
[i
], n_alloc
);
422 rho
+= 1/det(fr
.box
);
423 for (k
= 0; k
< nk
; k
++)
425 for (d
= 0; d
< DIM
; d
++)
427 kfac
[k
][d
] = 2*M_PI
*v0
[k
][d
]/fr
.box
[d
][d
];
430 for (i
= 0; i
< ntc
; i
++)
435 for (i
= 0; i
< gnx
; i
++)
442 for (j
= 0; j
< atndx
[index
[i
]+1] - atndx
[index
[i
]]; j
++)
444 at
= atndx
[index
[i
]] + j
;
445 m
= top
.atoms
.atom
[at
].m
;
446 mv_mol
[XX
] += m
*fr
.v
[at
][XX
];
447 mv_mol
[YY
] += m
*fr
.v
[at
][YY
];
448 mv_mol
[ZZ
] += m
*fr
.v
[at
][ZZ
];
449 cm_mol
[XX
] += m
*fr
.x
[at
][XX
];
450 cm_mol
[YY
] += m
*fr
.x
[at
][YY
];
451 cm_mol
[ZZ
] += m
*fr
.x
[at
][ZZ
];
454 svmul(1.0/mtot
, cm_mol
, cm_mol
);
458 svmul(top
.atoms
.atom
[index
[i
]].m
, fr
.v
[index
[i
]], mv_mol
);
463 copy_rvec(fr
.x
[index
[i
]], cm_mol
);
466 for (k
= 0; k
< nk
; k
++)
468 sx
= sin(iprod(kfac
[k
], cm_mol
));
469 cx
= cos(iprod(kfac
[k
], cm_mol
));
470 tc
[j
][nframes
] += sx
*iprod(v1
[k
], mv_mol
);
472 tc
[j
][nframes
] += cx
*iprod(v1
[k
], mv_mol
);
474 tc
[j
][nframes
] += sx
*iprod(v2
[k
], mv_mol
);
476 tc
[j
][nframes
] += cx
*iprod(v2
[k
], mv_mol
);
484 while (read_next_frame(oenv
, status
, &fr
));
487 dt
= (t1
-t0
)/(nframes
-1);
489 rho
*= sysmass
/nframes
*AMU
/(NANO
*NANO
*NANO
);
490 fprintf(stdout
, "Density = %g (kg/m^3)\n", rho
);
491 process_tcaf(nframes
, dt
, nkc
, tc
, kfac
, rho
, wt
,
492 opt2fn_null("-ot", NFILE
, fnm
),
493 opt2fn("-oa", NFILE
, fnm
), opt2fn("-o", NFILE
, fnm
),
494 opt2fn("-of", NFILE
, fnm
), opt2fn_null("-oc", NFILE
, fnm
),
495 opt2fn("-ov", NFILE
, fnm
), oenv
);