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