Remove continuation from convert_tpr
[gromacs.git] / src / gromacs / gmxana / gmx_dielectric.cpp
blobc9d4bb53de7afff883d10b10935feced9c8ffaf1
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, 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 <cstdlib>
42 #include <cstring>
44 #include <algorithm>
46 #include "gromacs/commandline/viewit.h"
47 #include "gromacs/correlationfunctions/expfit.h"
48 #include "gromacs/correlationfunctions/integrate.h"
49 #include "gromacs/fft/fft.h"
50 #include "gromacs/fileio/xvgr.h"
51 #include "gromacs/gmxana/gmx_ana.h"
52 #include "gromacs/gmxana/gstat.h"
53 #include "gromacs/math/gmxcomplex.h"
54 #include "gromacs/math/utilities.h"
55 #include "gromacs/utility/arraysize.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/futil.h"
58 #include "gromacs/utility/pleasecite.h"
59 #include "gromacs/utility/smalloc.h"
61 /* Determines at which point in the array the fit should start */
62 int calc_nbegin(int nx, real x[], real tbegin)
64 int nbegin;
66 /* Assume input x is sorted */
67 for (nbegin = 0; (nbegin < nx) && (x[nbegin] <= tbegin); nbegin++)
71 if ((nbegin == nx) || (nbegin == 0))
73 gmx_fatal(FARGS, "Begin time %f not in x-domain [%f through %f]\n",
74 tbegin, x[0], x[nx-1]);
77 /* Take the one closest to tbegin */
78 if (std::abs(x[nbegin]-tbegin) > std::abs(x[nbegin-1]-tbegin))
80 nbegin--;
83 printf("nbegin = %d, x[nbegin] = %g, tbegin = %g\n",
84 nbegin, x[nbegin], tbegin);
86 return nbegin;
89 real numerical_deriv(int nx, real x[], real y[], real fity[], real combined[], real dy[],
90 real tendInt, int nsmooth)
92 FILE *tmpfp;
93 int i, nbegin, i0, i1;
94 real fac, fx, fy, integralSmth;
96 nbegin = calc_nbegin(nx, x, tendInt);
97 if (nsmooth == 0)
99 for (i = 0; (i < nbegin); i++)
101 combined[i] = y[i];
103 fac = y[nbegin]/fity[nbegin];
104 printf("scaling fitted curve by %g\n", fac);
105 for (i = nbegin; (i < nx); i++)
107 combined[i] = fity[i]*fac;
110 else
112 i0 = std::max(0, nbegin);
113 i1 = std::min(nx-1, nbegin+nsmooth);
114 printf("Making smooth transition from %d through %d\n", i0, i1);
115 for (i = 0; (i < i0); i++)
117 combined[i] = y[i];
119 for (i = i0; (i <= i1); i++)
121 fx = static_cast<real>(i1-i)/(i1-i0);
122 fy = static_cast<real>(i-i0)/(i1-i0);
123 if (debug)
125 fprintf(debug, "x: %g factors for smoothing: %10g %10g\n", x[i], fx, fy);
127 combined[i] = fx*y[i] + fy*fity[i];
129 for (i = i1+1; (i < nx); i++)
131 combined[i] = fity[i];
135 tmpfp = gmx_ffopen("integral_smth.xvg", "w");
136 integralSmth = print_and_integrate(tmpfp, nx, x[1]-x[0], combined, NULL, 1);
137 printf("SMOOTH integral = %10.5e\n", integralSmth);
139 dy[0] = (combined[1]-combined[0])/(x[1]-x[0]);
140 for (i = 1; (i < nx-1); i++)
142 dy[i] = (combined[i+1]-combined[i-1])/(x[i+1]-x[i-1]);
144 dy[nx-1] = (combined[nx-1]-combined[nx-2])/(x[nx-1]-x[nx-2]);
146 for (i = 0; (i < nx); i++)
148 dy[i] *= -1;
151 return integralSmth;
154 void do_four(const char *fn, const char *cn, int nx, real x[], real dy[],
155 real eps0, real epsRF, const gmx_output_env_t *oenv)
157 FILE *fp, *cp;
158 t_complex *tmp, gw, hw, kw;
159 int i, nnx, nxsav;
160 real fac, nu, dt, maxeps, numax;
161 gmx_fft_t fft;
162 int fftcode;
164 nxsav = nx;
165 /*while ((dy[nx-1] == 0.0) && (nx > 0))
166 nx--;*/
167 if (nx == 0)
169 gmx_fatal(FARGS, "Empty dataset in %s, line %d!", __FILE__, __LINE__);
172 nnx = 1;
173 while (nnx < 2*nx)
175 nnx *= 2;
178 snew(tmp, 2*nnx);
179 printf("Doing FFT of %d points\n", nnx);
180 for (i = 0; (i < nx); i++)
182 tmp[i].re = dy[i];
184 if ((fftcode = gmx_fft_init_1d_real(&fft, nnx,
185 GMX_FFT_FLAG_NONE)) != 0)
187 gmx_fatal(FARGS, "gmx_fft_init_1d_real returned %d", fftcode);
189 if ((fftcode = gmx_fft_1d_real(fft, GMX_FFT_COMPLEX_TO_REAL,
190 (void *)tmp, (void *)tmp)) != 0)
192 gmx_fatal(FARGS, "gmx_fft_1d_real returned %d", fftcode);
194 gmx_fft_destroy(fft);
195 dt = x[1]-x[0];
196 if (epsRF == 0)
198 fac = (eps0-1)/tmp[0].re;
200 else
202 fac = ((eps0-1)/(2*epsRF+eps0))/tmp[0].re;
204 fp = xvgropen(fn, "Epsilon(\\8w\\4)", "Freq. (GHz)", "eps", oenv);
205 cp = xvgropen(cn, "Cole-Cole plot", "Eps'", "Eps''", oenv);
206 maxeps = 0;
207 numax = 0;
208 for (i = 0; (i < nxsav); i++)
210 if (epsRF == 0)
212 kw.re = 1+fac*tmp[i].re;
213 kw.im = 1+fac*tmp[i].im;
215 else
217 gw = rcmul(fac, tmp[i]);
218 hw = rcmul(2*epsRF, gw);
219 hw.re += 1.0;
220 gw.re = 1.0 - gw.re;
221 gw.im = -gw.im;
222 kw = cdiv(hw, gw);
224 kw.im *= -1;
226 nu = (i+1)*1000.0/(nnx*dt);
227 if (kw.im > maxeps)
229 maxeps = kw.im;
230 numax = nu;
233 fprintf(fp, "%10.5e %10.5e %10.5e\n", nu, kw.re, kw.im);
234 fprintf(cp, "%10.5e %10.5e\n", kw.re, kw.im);
236 printf("MAXEPS = %10.5e at frequency %10.5e GHz (tauD = %8.1f ps)\n",
237 maxeps, numax, 1000/(2*M_PI*numax));
238 xvgrclose(fp);
239 xvgrclose(cp);
240 sfree(tmp);
243 int gmx_dielectric(int argc, char *argv[])
245 const char *desc[] = {
246 "[THISMODULE] calculates frequency dependent dielectric constants",
247 "from the autocorrelation function of the total dipole moment in",
248 "your simulation. This ACF can be generated by [gmx-dipoles].",
249 "The functional forms of the available functions are:",
251 " * One parameter: y = [EXP]-a[SUB]1[sub] x[exp],",
252 " * Two parameters: y = a[SUB]2[sub] [EXP]-a[SUB]1[sub] x[exp],",
253 " * Three parameters: y = a[SUB]2[sub] [EXP]-a[SUB]1[sub] x[exp] + (1 - a[SUB]2[sub]) [EXP]-a[SUB]3[sub] x[exp].",
255 "Start values for the fit procedure can be given on the command line.",
256 "It is also possible to fix parameters at their start value, use [TT]-fix[tt]",
257 "with the number of the parameter you want to fix.",
258 "[PAR]",
259 "Three output files are generated, the first contains the ACF,",
260 "an exponential fit to it with 1, 2 or 3 parameters, and the",
261 "numerical derivative of the combination data/fit.",
262 "The second file contains the real and imaginary parts of the",
263 "frequency-dependent dielectric constant, the last gives a plot",
264 "known as the Cole-Cole plot, in which the imaginary",
265 "component is plotted as a function of the real component.",
266 "For a pure exponential relaxation (Debye relaxation) the latter",
267 "plot should be one half of a circle."
269 t_filenm fnm[] = {
270 { efXVG, "-f", "dipcorr", ffREAD },
271 { efXVG, "-d", "deriv", ffWRITE },
272 { efXVG, "-o", "epsw", ffWRITE },
273 { efXVG, "-c", "cole", ffWRITE }
275 #define NFILE asize(fnm)
276 gmx_output_env_t *oenv;
277 int i, j, nx, ny, nxtail, eFitFn, nfitparm;
278 real dt, integral, fitintegral, fac, rffac;
279 double *fitparms;
280 double **yd;
281 real **y;
282 const char *legend[] = { "Correlation", "Std. Dev.", "Fit", "Combined", "Derivative" };
283 static int fix = 0, bX = 1, nsmooth = 3;
284 static real tendInt = 5.0, tbegin = 5.0, tend = 500.0;
285 static real A = 0.5, tau1 = 10.0, tau2 = 1.0, eps0 = 80, epsRF = 78.5, tail = 500.0;
286 real lambda;
287 t_pargs pa[] = {
288 { "-x1", FALSE, etBOOL, {&bX},
289 "use first column as [IT]x[it]-axis rather than first data set" },
290 { "-eint", FALSE, etREAL, {&tendInt},
291 "Time to end the integration of the data and start to use the fit"},
292 { "-bfit", FALSE, etREAL, {&tbegin},
293 "Begin time of fit" },
294 { "-efit", FALSE, etREAL, {&tend},
295 "End time of fit" },
296 { "-tail", FALSE, etREAL, {&tail},
297 "Length of function including data and tail from fit" },
298 { "-A", FALSE, etREAL, {&A},
299 "Start value for fit parameter A" },
300 { "-tau1", FALSE, etREAL, {&tau1},
301 "Start value for fit parameter [GRK]tau[grk]1" },
302 { "-tau2", FALSE, etREAL, {&tau2},
303 "Start value for fit parameter [GRK]tau[grk]2" },
304 { "-eps0", FALSE, etREAL, {&eps0},
305 "[GRK]epsilon[grk]0 of your liquid" },
306 { "-epsRF", FALSE, etREAL, {&epsRF},
307 "[GRK]epsilon[grk] of the reaction field used in your simulation. A value of 0 means infinity." },
308 { "-fix", FALSE, etINT, {&fix},
309 "Fix parameters at their start values, A (2), tau1 (1), or tau2 (4)" },
310 { "-ffn", FALSE, etENUM, {s_ffn},
311 "Fit function" },
312 { "-nsmooth", FALSE, etINT, {&nsmooth},
313 "Number of points for smoothing" }
316 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
317 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
319 return 0;
321 please_cite(stdout, "Spoel98a");
322 printf("WARNING: non-polarizable models can never yield an infinite\n"
323 "dielectric constant that is different from 1. This is incorrect\n"
324 "in the reference given above (Spoel98a).\n\n");
327 nx = read_xvg(opt2fn("-f", NFILE, fnm), &yd, &ny);
328 dt = yd[0][1] - yd[0][0];
329 nxtail = std::min(static_cast<int>(tail/dt), nx);
331 printf("Read data set containing %d colums and %d rows\n", ny, nx);
332 printf("Assuming (from data) that timestep is %g, nxtail = %d\n",
333 dt, nxtail);
334 snew(y, 6);
335 for (i = 0; (i < ny); i++)
337 snew(y[i], std::max(nx, nxtail));
339 for (i = 0; (i < nx); i++)
341 y[0][i] = yd[0][i];
342 for (j = 1; (j < ny); j++)
344 y[j][i] = yd[j][i];
347 if (nxtail > nx)
349 for (i = nx; (i < nxtail); i++)
351 y[0][i] = dt*i+y[0][0];
352 for (j = 1; (j < ny); j++)
354 y[j][i] = 0.0;
357 nx = nxtail;
361 /* We have read a file WITHOUT standard deviations, so we make our own... */
362 if (ny == 2)
364 printf("Creating standard deviation numbers ...\n");
365 snew(y[2], nx);
367 fac = 1.0/nx;
368 for (i = 0; (i < nx); i++)
370 y[2][i] = fac;
374 eFitFn = sffn2effn(s_ffn);
375 nfitparm = effnNparams(eFitFn);
376 snew(fitparms, 4);
377 fitparms[0] = tau1;
378 if (nfitparm > 1)
380 fitparms[1] = A;
382 if (nfitparm > 2)
384 fitparms[2] = tau2;
388 snew(y[3], nx);
389 snew(y[4], nx);
390 snew(y[5], nx);
392 integral = print_and_integrate(NULL, calc_nbegin(nx, y[0], tbegin),
393 dt, y[1], NULL, 1);
394 integral += do_lmfit(nx, y[1], y[2], dt, y[0], tbegin, tend,
395 oenv, TRUE, eFitFn, fitparms, fix, NULL);
396 for (i = 0; i < nx; i++)
398 y[3][i] = fit_function(eFitFn, fitparms, y[0][i]);
401 if (epsRF == 0)
403 /* This means infinity! */
404 lambda = 0;
405 rffac = 1;
407 else
409 lambda = (eps0 - 1.0)/(2*epsRF - 1.0);
410 rffac = (2*epsRF+eps0)/(2*epsRF+1);
412 printf("DATA INTEGRAL: %5.1f, tauD(old) = %5.1f ps, "
413 "tau_slope = %5.1f, tau_slope,D = %5.1f ps\n",
414 integral, integral*rffac, fitparms[0], fitparms[0]*rffac);
416 printf("tau_D from tau1 = %8.3g , eps(Infty) = %8.3f\n",
417 fitparms[0]*(1 + fitparms[1]*lambda),
418 1 + ((1 - fitparms[1])*(eps0 - 1))/(1 + fitparms[1]*lambda));
420 fitintegral = numerical_deriv(nx, y[0], y[1], y[3], y[4], y[5], tendInt, nsmooth);
421 printf("FIT INTEGRAL (tau_M): %5.1f, tau_D = %5.1f\n",
422 fitintegral, fitintegral*rffac);
424 /* Now we have the negative gradient of <Phi(0) Phi(t)> */
425 write_xvg(opt2fn("-d", NFILE, fnm), "Data", nx-1, 6, y, legend, oenv);
427 /* Do FFT and analysis */
428 do_four(opt2fn("-o", NFILE, fnm), opt2fn("-c", NFILE, fnm),
429 nx-1, y[0], y[5], eps0, epsRF, oenv);
431 do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy");
432 do_view(oenv, opt2fn("-c", NFILE, fnm), NULL);
433 do_view(oenv, opt2fn("-d", NFILE, fnm), "-nxy");
435 return 0;