3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
52 #include "gmxcomplex.h"
55 #include "gmx_fatal.h"
58 /* Determines at which point in the array the fit should start */
59 int calc_nbegin(int nx,real x[],real tbegin)
64 /* Assume input x is sorted */
65 for(nbegin=0; (nbegin < nx) && (x[nbegin] <= tbegin); nbegin++)
67 if ((nbegin == nx) || (nbegin == 0))
68 gmx_fatal(FARGS,"Begin time %f not in x-domain [%f through %f]\n",
71 /* Take the one closest to tbegin */
72 if (fabs(x[nbegin]-tbegin) > fabs(x[nbegin-1]-tbegin))
75 printf("nbegin = %d, x[nbegin] = %g, tbegin = %g\n",
76 nbegin,x[nbegin],tbegin);
81 real numerical_deriv(int nx,real x[],real y[],real fity[],real combined[],real dy[],
82 real tendInt,int nsmooth)
86 real fac,fx,fy,integralSmth;
88 nbegin = calc_nbegin(nx,x,tendInt);
90 for(i=0; (i<nbegin); i++)
92 fac = y[nbegin]/fity[nbegin];
93 printf("scaling fitted curve by %g\n",fac);
94 for(i=nbegin; (i<nx); i++)
95 combined[i]=fity[i]*fac;
99 i1 = min(nx-1,nbegin+nsmooth);
100 printf("Making smooth transition from %d through %d\n",i0,i1);
101 for(i=0; (i<i0); i++)
103 for(i=i0; (i<=i1); i++) {
104 fx = (i1-i)/(real)(i1-i0);
105 fy = (i-i0)/(real)(i1-i0);
107 fprintf(debug,"x: %g factors for smoothing: %10g %10g\n",x[i],fx,fy);
108 combined[i] = fx*y[i] + fy*fity[i];
110 for(i=i1+1; (i<nx); i++)
114 tmpfp = ffopen("integral_smth.xvg","w");
115 integralSmth=print_and_integrate(tmpfp,nx,x[1]-x[0],combined,NULL,1);
116 printf("SMOOTH integral = %10.5e\n",integralSmth);
118 dy[0] = (combined[1]-combined[0])/(x[1]-x[0]);
119 for(i=1; (i<nx-1); i++) {
120 dy[i] = (combined[i+1]-combined[i-1])/(x[i+1]-x[i-1]);
122 dy[nx-1] = (combined[nx-1]-combined[nx-2])/(x[nx-1]-x[nx-2]);
124 for(i=0; (i<nx); i++)
130 void do_four(const char *fn,const char *cn,int nx,real x[],real dy[],
131 real eps0,real epsRF, const output_env_t oenv)
134 t_complex *tmp,gw,hw,kw;
136 real fac,nu,dt,*ptr,maxeps,numax;
139 /*while ((dy[nx-1] == 0.0) && (nx > 0))
142 gmx_fatal(FARGS,"Empty dataset in %s, line %d!",__FILE__,__LINE__);
149 printf("Doing FFT of %d points\n",nnx);
150 for(i=0; (i<nx); i++)
157 fac = (eps0-1)/tmp[0].re;
159 fac=((eps0-1)/(2*epsRF+eps0))/tmp[0].re;
160 fp=xvgropen(fn,"Epsilon(\\8w\\4)","Freq. (GHz)","eps",oenv);
161 cp=xvgropen(cn,"Cole-Cole plot","Eps'","Eps''",oenv);
164 for(i=0; (i<nxsav); i++) {
166 kw.re = 1+fac*tmp[i].re;
167 kw.im = 1+fac*tmp[i].im;
170 gw = rcmul(fac,tmp[i]);
171 hw = rcmul(2*epsRF,gw);
179 nu = (i+1)*1000.0/(nnx*dt);
180 if (kw.im > maxeps) {
185 fprintf(fp,"%10.5e %10.5e %10.5e\n",nu,kw.re,kw.im);
186 fprintf(cp,"%10.5e %10.5e\n",kw.re,kw.im);
188 printf("MAXEPS = %10.5e at frequency %10.5e GHz (tauD = %8.1f ps)\n",
189 maxeps,numax,1000/(2*M_PI*numax));
195 int gmx_dielectric(int argc,char *argv[])
197 const char *desc[] = {
198 "[TT]g_dielectric[tt] calculates frequency dependent dielectric constants",
199 "from the autocorrelation function of the total dipole moment in",
200 "your simulation. This ACF can be generated by [TT]g_dipoles[tt].",
201 "The functional forms of the available functions are:[PAR]",
202 "One parameter: y = [EXP]-a[SUB]1[sub] x[exp],[BR]",
203 "Two parameters: y = a[SUB]2[sub] [EXP]-a[SUB]1[sub] x[exp],[BR]",
204 "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]",
205 "Start values for the fit procedure can be given on the command line.",
206 "It is also possible to fix parameters at their start value, use [TT]-fix[tt]",
207 "with the number of the parameter you want to fix.",
209 "Three output files are generated, the first contains the ACF,",
210 "an exponential fit to it with 1, 2 or 3 parameters, and the",
211 "numerical derivative of the combination data/fit.",
212 "The second file contains the real and imaginary parts of the",
213 "frequency-dependent dielectric constant, the last gives a plot",
214 "known as the Cole-Cole plot, in which the imaginary",
215 "component is plotted as a function of the real component.",
216 "For a pure exponential relaxation (Debye relaxation) the latter",
217 "plot should be one half of a circle."
220 { efXVG, "-f", "dipcorr",ffREAD },
221 { efXVG, "-d", "deriv", ffWRITE },
222 { efXVG, "-o", "epsw", ffWRITE },
223 { efXVG, "-c", "cole", ffWRITE }
225 #define NFILE asize(fnm)
227 int i,j,nx,ny,nxtail,eFitFn,nfitparm;
228 real dt,integral,fitintegral,*fitparms,fac,rffac;
231 const char *legend[] = { "Correlation", "Std. Dev.", "Fit", "Combined", "Derivative" };
232 static int fix=0,bFour = 0,bX = 1,nsmooth=3;
233 static real tendInt=5.0,tbegin=5.0,tend=500.0;
234 static real A=0.5,tau1=10.0,tau2=1.0,eps0=80,epsRF=78.5,tail=500.0;
237 { "-fft", FALSE, etBOOL, {&bFour},
238 "use fast fourier transform for correlation function" },
239 { "-x1", FALSE, etBOOL, {&bX},
240 "use first column as [IT]x[it]-axis rather than first data set" },
241 { "-eint", FALSE, etREAL, {&tendInt},
242 "Time to end the integration of the data and start to use the fit"},
243 { "-bfit", FALSE, etREAL, {&tbegin},
244 "Begin time of fit" },
245 { "-efit", FALSE, etREAL, {&tend},
247 { "-tail", FALSE, etREAL, {&tail},
248 "Length of function including data and tail from fit" },
249 { "-A", FALSE, etREAL, {&A},
250 "Start value for fit parameter A" },
251 { "-tau1", FALSE, etREAL, {&tau1},
252 "Start value for fit parameter [GRK]tau[grk]1" },
253 { "-tau2", FALSE, etREAL, {&tau2},
254 "Start value for fit parameter [GRK]tau[grk]2" },
255 { "-eps0", FALSE, etREAL, {&eps0},
256 "[GRK]epsilon[grk]0 of your liquid" },
257 { "-epsRF", FALSE, etREAL, {&epsRF},
258 "[GRK]epsilon[grk] of the reaction field used in your simulation. A value of 0 means infinity." },
259 { "-fix", FALSE, etINT, {&fix},
260 "Fix parameters at their start values, A (2), tau1 (1), or tau2 (4)" },
261 { "-ffn", FALSE, etENUM, {s_ffn},
263 { "-nsmooth", FALSE, etINT, {&nsmooth},
264 "Number of points for smoothing" }
267 CopyRight(stderr,argv[0]);
268 parse_common_args(&argc,argv,PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE,
269 NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL,&oenv);
270 please_cite(stdout,"Spoel98a");
271 printf("WARNING: non-polarizable models can never yield an infinite\n"
272 "dielectric constant that is different from 1. This is incorrect\n"
273 "in the reference given above (Spoel98a).\n\n");
276 nx = read_xvg(opt2fn("-f",NFILE,fnm),&yd,&ny);
277 dt = yd[0][1] - yd[0][0];
278 nxtail = min(tail/dt,nx);
280 printf("Read data set containing %d colums and %d rows\n",ny,nx);
281 printf("Assuming (from data) that timestep is %g, nxtail = %d\n",
284 for(i=0; (i<ny); i++)
285 snew(y[i],max(nx,nxtail));
286 for(i=0; (i<nx); i++) {
288 for(j=1; (j<ny); j++)
292 for(i=nx; (i<nxtail); i++) {
293 y[0][i] = dt*i+y[0][0];
294 for(j=1; (j<ny); j++)
301 /* We have read a file WITHOUT standard deviations, so we make our own... */
303 printf("Creating standard deviation numbers ...\n");
308 for(i=0; (i<nx); i++)
312 eFitFn = sffn2effn(s_ffn);
313 nfitparm = nfp_ffn[eFitFn];
326 integral = print_and_integrate(NULL,calc_nbegin(nx,y[0],tbegin),
328 integral += do_lmfit(nx,y[1],y[2],dt,y[0],tbegin,tend,
329 oenv,TRUE,eFitFn,fitparms,fix);
331 y[3][i] = fit_function(eFitFn,fitparms,y[0][i]);
334 /* This means infinity! */
339 lambda = (eps0 - 1.0)/(2*epsRF - 1.0);
340 rffac = (2*epsRF+eps0)/(2*epsRF+1);
342 printf("DATA INTEGRAL: %5.1f, tauD(old) = %5.1f ps, "
343 "tau_slope = %5.1f, tau_slope,D = %5.1f ps\n",
344 integral,integral*rffac,fitparms[0],fitparms[0]*rffac);
346 printf("tau_D from tau1 = %8.3g , eps(Infty) = %8.3f\n",
347 fitparms[0]*(1 + fitparms[1]*lambda),
348 1 + ((1 - fitparms[1])*(eps0 - 1))/(1 + fitparms[1]*lambda));
350 fitintegral=numerical_deriv(nx,y[0],y[1],y[3],y[4],y[5],tendInt,nsmooth);
351 printf("FIT INTEGRAL (tau_M): %5.1f, tau_D = %5.1f\n",
352 fitintegral,fitintegral*rffac);
354 /* Now we have the negative gradient of <Phi(0) Phi(t)> */
355 write_xvg(opt2fn("-d",NFILE,fnm),"Data",nx-1,6,y,legend,oenv);
357 /* Do FFT and analysis */
358 do_four(opt2fn("-o",NFILE,fnm),opt2fn("-c",NFILE,fnm),
359 nx-1,y[0],y[5],eps0,epsRF,oenv);
361 do_view(oenv,opt2fn("-o",NFILE,fnm),"-nxy");
362 do_view(oenv,opt2fn("-c",NFILE,fnm),NULL);
363 do_view(oenv,opt2fn("-d",NFILE,fnm),"-nxy");