5f5e719c24c758e176c9b7bf370280a9c00d80cd
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_dielectric.cpp
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,2016,2017,2018,2019, 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 <cmath>
40 #include <cstdio>
41 #include <cstdlib>
42 #include <cstring>
43
44 #include <algorithm>
45
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"
60
61 /* Determines at which point in the array the fit should start */
62 static int calc_nbegin(int nx, real x[], real tbegin)
63 {
64     int nbegin;
65
66     /* Assume input x is sorted */
67     for (nbegin = 0; (nbegin < nx) && (x[nbegin] <= tbegin); nbegin++)
68     {
69         // Deliberately left empty
70     }
71     if ((nbegin == nx) || (nbegin == 0))
72     {
73         gmx_fatal(FARGS, "Begin time %f not in x-domain [%f through %f]\n", tbegin, x[0], x[nx - 1]);
74     }
75
76     /* Take the one closest to tbegin */
77     if (std::abs(x[nbegin] - tbegin) > std::abs(x[nbegin - 1] - tbegin))
78     {
79         nbegin--;
80     }
81
82     printf("nbegin = %d, x[nbegin] = %g, tbegin = %g\n", nbegin, x[nbegin], tbegin);
83
84     return nbegin;
85 }
86
87 static real
88 numerical_deriv(int nx, real x[], const real y[], const real fity[], real combined[], real dy[], real tendInt, int nsmooth)
89 {
90     FILE* tmpfp;
91     int   i, nbegin, i0, i1;
92     real  fac, fx, fy, integralSmth;
93
94     nbegin = calc_nbegin(nx, x, tendInt);
95     if (nsmooth == 0)
96     {
97         for (i = 0; (i < nbegin); i++)
98         {
99             combined[i] = y[i];
100         }
101         fac = y[nbegin] / fity[nbegin];
102         printf("scaling fitted curve by %g\n", fac);
103         for (i = nbegin; (i < nx); i++)
104         {
105             combined[i] = fity[i] * fac;
106         }
107     }
108     else
109     {
110         i0 = std::max(0, nbegin);
111         i1 = std::min(nx - 1, nbegin + nsmooth);
112         printf("Making smooth transition from %d through %d\n", i0, i1);
113         for (i = 0; (i < i0); i++)
114         {
115             combined[i] = y[i];
116         }
117         for (i = i0; (i <= i1); i++)
118         {
119             fx = static_cast<real>(i1 - i) / (i1 - i0);
120             fy = static_cast<real>(i - i0) / (i1 - i0);
121             if (debug)
122             {
123                 fprintf(debug, "x: %g factors for smoothing: %10g %10g\n", x[i], fx, fy);
124             }
125             combined[i] = fx * y[i] + fy * fity[i];
126         }
127         for (i = i1 + 1; (i < nx); i++)
128         {
129             combined[i] = fity[i];
130         }
131     }
132
133     tmpfp        = gmx_ffopen("integral_smth.xvg", "w");
134     integralSmth = print_and_integrate(tmpfp, nx, x[1] - x[0], combined, nullptr, 1);
135     printf("SMOOTH integral = %10.5e\n", integralSmth);
136
137     dy[0] = (combined[1] - combined[0]) / (x[1] - x[0]);
138     for (i = 1; (i < nx - 1); i++)
139     {
140         dy[i] = (combined[i + 1] - combined[i - 1]) / (x[i + 1] - x[i - 1]);
141     }
142     dy[nx - 1] = (combined[nx - 1] - combined[nx - 2]) / (x[nx - 1] - x[nx - 2]);
143
144     for (i = 0; (i < nx); i++)
145     {
146         dy[i] *= -1;
147     }
148
149     return integralSmth;
150 }
151
152 static void do_four(const char*             fn,
153                     const char*             cn,
154                     int                     nx,
155                     const real              x[],
156                     const real              dy[],
157                     real                    eps0,
158                     real                    epsRF,
159                     const gmx_output_env_t* oenv)
160 {
161     FILE *     fp, *cp;
162     t_complex *tmp, gw, hw, kw;
163     int        i, nnx, nxsav;
164     real       fac, nu, dt, maxeps, numax;
165     gmx_fft_t  fft;
166     int        fftcode;
167
168     nxsav = nx;
169     /*while ((dy[nx-1] == 0.0) && (nx > 0))
170        nx--;*/
171     if (nx == 0)
172     {
173         gmx_fatal(FARGS, "Empty dataset in %s, line %d!", __FILE__, __LINE__);
174     }
175
176     nnx = 1;
177     while (nnx < 2 * nx)
178     {
179         nnx *= 2;
180     }
181
182     snew(tmp, 2 * nnx);
183     printf("Doing FFT of %d points\n", nnx);
184     for (i = 0; (i < nx); i++)
185     {
186         tmp[i].re = dy[i];
187     }
188     if ((fftcode = gmx_fft_init_1d_real(&fft, nnx, GMX_FFT_FLAG_NONE)) != 0)
189     {
190         gmx_fatal(FARGS, "gmx_fft_init_1d_real returned %d", fftcode);
191     }
192     if ((fftcode = gmx_fft_1d_real(fft, GMX_FFT_COMPLEX_TO_REAL, tmp, tmp)) != 0)
193     {
194         gmx_fatal(FARGS, "gmx_fft_1d_real returned %d", fftcode);
195     }
196     gmx_fft_destroy(fft);
197     dt = x[1] - x[0];
198     if (epsRF == 0)
199     {
200         fac = (eps0 - 1) / tmp[0].re;
201     }
202     else
203     {
204         fac = ((eps0 - 1) / (2 * epsRF + eps0)) / tmp[0].re;
205     }
206     fp     = xvgropen(fn, "Epsilon(\\8w\\4)", "Freq. (GHz)", "eps", oenv);
207     cp     = xvgropen(cn, "Cole-Cole plot", "Eps'", "Eps''", oenv);
208     maxeps = 0;
209     numax  = 0;
210     for (i = 0; (i < nxsav); i++)
211     {
212         if (epsRF == 0)
213         {
214             kw.re = 1 + fac * tmp[i].re;
215             kw.im = 1 + fac * tmp[i].im;
216         }
217         else
218         {
219             gw = rcmul(fac, tmp[i]);
220             hw = rcmul(2 * epsRF, gw);
221             hw.re += 1.0;
222             gw.re = 1.0 - gw.re;
223             gw.im = -gw.im;
224             kw    = cdiv(hw, gw);
225         }
226         kw.im *= -1;
227
228         nu = (i + 1) * 1000.0 / (nnx * dt);
229         if (kw.im > maxeps)
230         {
231             maxeps = kw.im;
232             numax  = nu;
233         }
234
235         fprintf(fp, "%10.5e  %10.5e  %10.5e\n", nu, kw.re, kw.im);
236         fprintf(cp, "%10.5e  %10.5e\n", kw.re, kw.im);
237     }
238     printf("MAXEPS = %10.5e at frequency %10.5e GHz (tauD = %8.1f ps)\n", maxeps, numax,
239            1000 / (2 * M_PI * numax));
240     xvgrclose(fp);
241     xvgrclose(cp);
242     sfree(tmp);
243 }
244
245 int gmx_dielectric(int argc, char* argv[])
246 {
247     const char* desc[] = {
248         "[THISMODULE] calculates frequency dependent dielectric constants",
249         "from the autocorrelation function of the total dipole moment in",
250         "your simulation. This ACF can be generated by [gmx-dipoles].",
251         "The functional forms of the available functions are:",
252         "",
253         " * One parameter:    y = [EXP]-a[SUB]1[sub] x[exp],",
254         " * Two parameters:   y = a[SUB]2[sub] [EXP]-a[SUB]1[sub] x[exp],",
255         " * Three parameters: y = a[SUB]2[sub] [EXP]-a[SUB]1[sub] x[exp] + (1 - ",
256         "   a[SUB]2[sub]) [EXP]-a[SUB]3[sub] x[exp].",
257         "",
258         "Start values for the fit procedure can be given on the command line.",
259         "It is also possible to fix parameters at their start value, use [TT]-fix[tt]",
260         "with the number of the parameter you want to fix.",
261         "[PAR]",
262         "Three output files are generated, the first contains the ACF,",
263         "an exponential fit to it with 1, 2 or 3 parameters, and the",
264         "numerical derivative of the combination data/fit.",
265         "The second file contains the real and imaginary parts of the",
266         "frequency-dependent dielectric constant, the last gives a plot",
267         "known as the Cole-Cole plot, in which the imaginary",
268         "component is plotted as a function of the real component.",
269         "For a pure exponential relaxation (Debye relaxation) the latter",
270         "plot should be one half of a circle."
271     };
272     t_filenm fnm[] = { { efXVG, "-f", "dipcorr", ffREAD },
273                        { efXVG, "-d", "deriv", ffWRITE },
274                        { efXVG, "-o", "epsw", ffWRITE },
275                        { efXVG, "-c", "cole", ffWRITE } };
276 #define NFILE asize(fnm)
277     gmx_output_env_t* oenv;
278     int               i, j, nx, ny, nxtail, eFitFn, nfitparm;
279     real              dt, integral, fitintegral, fac, rffac;
280     double*           fitparms;
281     double**          yd;
282     real**            y;
283     const char*       legend[] = { "Correlation", "Std. Dev.", "Fit", "Combined", "Derivative" };
284     static int        fix = 0, bX = 1, nsmooth = 3;
285     static real       tendInt = 5.0, tbegin = 5.0, tend = 500.0;
286     static real       A = 0.5, tau1 = 10.0, tau2 = 1.0, eps0 = 80, epsRF = 78.5, tail = 500.0;
287     real              lambda;
288     t_pargs           pa[] = {
289         { "-x1",
290           FALSE,
291           etBOOL,
292           { &bX },
293           "use first column as [IT]x[it]-axis rather than first data set" },
294         { "-eint",
295           FALSE,
296           etREAL,
297           { &tendInt },
298           "Time to end the integration of the data and start to use the fit" },
299         { "-bfit", FALSE, etREAL, { &tbegin }, "Begin time of fit" },
300         { "-efit", FALSE, etREAL, { &tend }, "End time of fit" },
301         { "-tail", FALSE, etREAL, { &tail }, "Length of function including data and tail from fit" },
302         { "-A", FALSE, etREAL, { &A }, "Start value for fit parameter A" },
303         { "-tau1", FALSE, etREAL, { &tau1 }, "Start value for fit parameter [GRK]tau[grk]1" },
304         { "-tau2", FALSE, etREAL, { &tau2 }, "Start value for fit parameter [GRK]tau[grk]2" },
305         { "-eps0", FALSE, etREAL, { &eps0 }, "[GRK]epsilon[grk]0 of your liquid" },
306         { "-epsRF",
307           FALSE,
308           etREAL,
309           { &epsRF },
310           "[GRK]epsilon[grk] of the reaction field used in your simulation. A value of 0 means "
311           "infinity." },
312         { "-fix",
313           FALSE,
314           etINT,
315           { &fix },
316           "Fix parameters at their start values, A (2), tau1 (1), or tau2 (4)" },
317         { "-ffn", FALSE, etENUM, { s_ffn }, "Fit function" },
318         { "-nsmooth", FALSE, etINT, { &nsmooth }, "Number of points for smoothing" }
319     };
320
321     if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW, NFILE, fnm, asize(pa), pa,
322                            asize(desc), desc, 0, nullptr, &oenv))
323     {
324         return 0;
325     }
326     please_cite(stdout, "Spoel98a");
327     printf("WARNING: non-polarizable models can never yield an infinite\n"
328            "dielectric constant that is different from 1. This is incorrect\n"
329            "in the reference given above (Spoel98a).\n\n");
330
331
332     nx     = read_xvg(opt2fn("-f", NFILE, fnm), &yd, &ny);
333     dt     = yd[0][1] - yd[0][0];
334     nxtail = std::min(static_cast<int>(tail / dt), nx);
335
336     printf("Read data set containing %d colums and %d rows\n", ny, nx);
337     printf("Assuming (from data) that timestep is %g, nxtail = %d\n", dt, nxtail);
338     snew(y, 6);
339     for (i = 0; (i < ny); i++)
340     {
341         snew(y[i], std::max(nx, nxtail));
342     }
343     for (i = 0; (i < nx); i++)
344     {
345         y[0][i] = yd[0][i];
346         for (j = 1; (j < ny); j++)
347         {
348             y[j][i] = yd[j][i];
349         }
350     }
351     if (nxtail > nx)
352     {
353         for (i = nx; (i < nxtail); i++)
354         {
355             y[0][i] = dt * i + y[0][0];
356             for (j = 1; (j < ny); j++)
357             {
358                 y[j][i] = 0.0;
359             }
360         }
361         nx = nxtail;
362     }
363
364
365     /* We have read a file WITHOUT standard deviations, so we make our own... */
366     if (ny == 2)
367     {
368         printf("Creating standard deviation numbers ...\n");
369         snew(y[2], nx);
370
371         fac = 1.0 / nx;
372         for (i = 0; (i < nx); i++)
373         {
374             y[2][i] = fac;
375         }
376     }
377
378     eFitFn   = sffn2effn(s_ffn);
379     nfitparm = effnNparams(eFitFn);
380     snew(fitparms, 4);
381     fitparms[0] = tau1;
382     if (nfitparm > 1)
383     {
384         fitparms[1] = A;
385     }
386     if (nfitparm > 2)
387     {
388         fitparms[2] = tau2;
389     }
390
391
392     snew(y[3], nx);
393     snew(y[4], nx);
394     snew(y[5], nx);
395
396     integral = print_and_integrate(nullptr, calc_nbegin(nx, y[0], tbegin), dt, y[1], nullptr, 1);
397     integral += do_lmfit(nx, y[1], y[2], dt, y[0], tbegin, tend, oenv, TRUE, eFitFn, fitparms, fix, nullptr);
398     for (i = 0; i < nx; i++)
399     {
400         y[3][i] = fit_function(eFitFn, fitparms, y[0][i]);
401     }
402
403     if (epsRF == 0)
404     {
405         /* This means infinity! */
406         lambda = 0;
407         rffac  = 1;
408     }
409     else
410     {
411         lambda = (eps0 - 1.0) / (2 * epsRF - 1.0);
412         rffac  = (2 * epsRF + eps0) / (2 * epsRF + 1);
413     }
414     printf("DATA INTEGRAL: %5.1f, tauD(old) = %5.1f ps, "
415            "tau_slope = %5.1f, tau_slope,D = %5.1f ps\n",
416            integral, integral * rffac, fitparms[0], fitparms[0] * rffac);
417
418     printf("tau_D from tau1 = %8.3g , eps(Infty) = %8.3f\n", fitparms[0] * (1 + fitparms[1] * lambda),
419            1 + ((1 - fitparms[1]) * (eps0 - 1)) / (1 + fitparms[1] * lambda));
420
421     fitintegral = numerical_deriv(nx, y[0], y[1], y[3], y[4], y[5], tendInt, nsmooth);
422     printf("FIT INTEGRAL (tau_M): %5.1f, tau_D = %5.1f\n", fitintegral, fitintegral * rffac);
423
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);
426
427     /* Do FFT and analysis */
428     do_four(opt2fn("-o", NFILE, fnm), opt2fn("-c", NFILE, fnm), nx - 1, y[0], y[5], eps0, epsRF, oenv);
429
430     do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy");
431     do_view(oenv, opt2fn("-c", NFILE, fnm), nullptr);
432     do_view(oenv, opt2fn("-d", NFILE, fnm), "-nxy");
433
434     return 0;
435 }