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