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