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