int gmx_dielectric(int argc, char *argv[]) { const char *desc[] = { "[THISMODULE] calculates frequency dependent dielectric constants", "from the autocorrelation function of the total dipole moment in", "your simulation. This ACF can be generated by [gmx-dipoles].", "The functional forms of the available functions are:[PAR]", "One parameter: y = [EXP]-a[SUB]1[sub] x[exp],[BR]", "Two parameters: y = a[SUB]2[sub] [EXP]-a[SUB]1[sub] x[exp],[BR]", "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]", "Start values for the fit procedure can be given on the command line.", "It is also possible to fix parameters at their start value, use [TT]-fix[tt]", "with the number of the parameter you want to fix.", "[PAR]", "Three output files are generated, the first contains the ACF,", "an exponential fit to it with 1, 2 or 3 parameters, and the", "numerical derivative of the combination data/fit.", "The second file contains the real and imaginary parts of the", "frequency-dependent dielectric constant, the last gives a plot", "known as the Cole-Cole plot, in which the imaginary", "component is plotted as a function of the real component.", "For a pure exponential relaxation (Debye relaxation) the latter", "plot should be one half of a circle." }; t_filenm fnm[] = { { efXVG, "-f", "dipcorr", ffREAD }, { efXVG, "-d", "deriv", ffWRITE }, { efXVG, "-o", "epsw", ffWRITE }, { efXVG, "-c", "cole", ffWRITE } }; #define NFILE asize(fnm) output_env_t oenv; int i, j, nx, ny, nxtail, eFitFn, nfitparm; real dt, integral, fitintegral, *fitparms, fac, rffac; double **yd; real **y; const char *legend[] = { "Correlation", "Std. Dev.", "Fit", "Combined", "Derivative" }; static int fix = 0, bFour = 0, bX = 1, nsmooth = 3; static real tendInt = 5.0, tbegin = 5.0, tend = 500.0; static real A = 0.5, tau1 = 10.0, tau2 = 1.0, eps0 = 80, epsRF = 78.5, tail = 500.0; real lambda; t_pargs pa[] = { { "-fft", FALSE, etBOOL, {&bFour}, "use fast fourier transform for correlation function" }, { "-x1", FALSE, etBOOL, {&bX}, "use first column as [IT]x[it]-axis rather than first data set" }, { "-eint", FALSE, etREAL, {&tendInt}, "Time to end the integration of the data and start to use the fit"}, { "-bfit", FALSE, etREAL, {&tbegin}, "Begin time of fit" }, { "-efit", FALSE, etREAL, {&tend}, "End time of fit" }, { "-tail", FALSE, etREAL, {&tail}, "Length of function including data and tail from fit" }, { "-A", FALSE, etREAL, {&A}, "Start value for fit parameter A" }, { "-tau1", FALSE, etREAL, {&tau1}, "Start value for fit parameter [GRK]tau[grk]1" }, { "-tau2", FALSE, etREAL, {&tau2}, "Start value for fit parameter [GRK]tau[grk]2" }, { "-eps0", FALSE, etREAL, {&eps0}, "[GRK]epsilon[grk]0 of your liquid" }, { "-epsRF", FALSE, etREAL, {&epsRF}, "[GRK]epsilon[grk] of the reaction field used in your simulation. A value of 0 means infinity." }, { "-fix", FALSE, etINT, {&fix}, "Fix parameters at their start values, A (2), tau1 (1), or tau2 (4)" }, { "-ffn", FALSE, etENUM, {s_ffn}, "Fit function" }, { "-nsmooth", FALSE, etINT, {&nsmooth}, "Number of points for smoothing" } }; if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE, NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv)) { return 0; } please_cite(stdout, "Spoel98a"); printf("WARNING: non-polarizable models can never yield an infinite\n" "dielectric constant that is different from 1. This is incorrect\n" "in the reference given above (Spoel98a).\n\n"); nx = read_xvg(opt2fn("-f", NFILE, fnm), &yd, &ny); dt = yd[0][1] - yd[0][0]; nxtail = min(tail/dt, nx); printf("Read data set containing %d colums and %d rows\n", ny, nx); printf("Assuming (from data) that timestep is %g, nxtail = %d\n", dt, nxtail); snew(y, 6); for (i = 0; (i < ny); i++) { snew(y[i], max(nx, nxtail)); } for (i = 0; (i < nx); i++) { y[0][i] = yd[0][i]; for (j = 1; (j < ny); j++) { y[j][i] = yd[j][i]; } } if (nxtail > nx) { for (i = nx; (i < nxtail); i++) { y[0][i] = dt*i+y[0][0]; for (j = 1; (j < ny); j++) { y[j][i] = 0.0; } } nx = nxtail; } /* We have read a file WITHOUT standard deviations, so we make our own... */ if (ny == 2) { printf("Creating standard deviation numbers ...\n"); srenew(y, 3); snew(y[2], nx); fac = 1.0/((real)nx); for (i = 0; (i < nx); i++) { y[2][i] = fac; } } eFitFn = sffn2effn(s_ffn); nfitparm = nfp_ffn[eFitFn]; snew(fitparms, 4); fitparms[0] = tau1; if (nfitparm > 1) { fitparms[1] = A; } if (nfitparm > 2) { fitparms[2] = tau2; } snew(y[3], nx); snew(y[4], nx); snew(y[5], nx); integral = print_and_integrate(NULL, calc_nbegin(nx, y[0], tbegin), dt, y[1], NULL, 1); integral += do_lmfit(nx, y[1], y[2], dt, y[0], tbegin, tend, oenv, TRUE, eFitFn, fitparms, fix); for (i = 0; i < nx; i++) { y[3][i] = fit_function(eFitFn, fitparms, y[0][i]); } if (epsRF == 0) { /* This means infinity! */ lambda = 0; rffac = 1; } else { lambda = (eps0 - 1.0)/(2*epsRF - 1.0); rffac = (2*epsRF+eps0)/(2*epsRF+1); } printf("DATA INTEGRAL: %5.1f, tauD(old) = %5.1f ps, " "tau_slope = %5.1f, tau_slope,D = %5.1f ps\n", integral, integral*rffac, fitparms[0], fitparms[0]*rffac); printf("tau_D from tau1 = %8.3g , eps(Infty) = %8.3f\n", fitparms[0]*(1 + fitparms[1]*lambda), 1 + ((1 - fitparms[1])*(eps0 - 1))/(1 + fitparms[1]*lambda)); fitintegral = numerical_deriv(nx, y[0], y[1], y[3], y[4], y[5], tendInt, nsmooth); printf("FIT INTEGRAL (tau_M): %5.1f, tau_D = %5.1f\n", fitintegral, fitintegral*rffac); /* Now we have the negative gradient of <Phi(0) Phi(t)> */ write_xvg(opt2fn("-d", NFILE, fnm), "Data", nx-1, 6, y, legend, oenv); /* Do FFT and analysis */ do_four(opt2fn("-o", NFILE, fnm), opt2fn("-c", NFILE, fnm), nx-1, y[0], y[5], eps0, epsRF, oenv); do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy"); do_view(oenv, opt2fn("-c", NFILE, fnm), NULL); do_view(oenv, opt2fn("-d", NFILE, fnm), "-nxy"); return 0; }
static void process_tcaf(int nframes, real dt, int nkc, real **tc, rvec *kfac, real rho, real wt, const char *fn_trans, const char *fn_tca, const char *fn_tc, const char *fn_tcf, const char *fn_cub, const char *fn_vk, const gmx_output_env_t *oenv) { FILE *fp, *fp_vk, *fp_cub = NULL; int nk, ntc; real **tcaf, **tcafc = NULL, eta, *sig; int i, j, k, kc; int ncorr; double fitparms[3]; nk = kset_c[nkc]; ntc = nk*NPK; if (fn_trans) { fp = xvgropen(fn_trans, "Transverse Current", "Time (ps)", "TC (nm/ps)", oenv); for (i = 0; i < nframes; i++) { fprintf(fp, "%g", i*dt); for (j = 0; j < ntc; j++) { fprintf(fp, " %g", tc[j][i]); } fprintf(fp, "\n"); } xvgrclose(fp); do_view(oenv, fn_trans, "-nxy"); } ncorr = (nframes+1)/2; if (ncorr > static_cast<int>(5*wt/dt+0.5)) { ncorr = static_cast<int>(5*wt/dt+0.5)+1; } snew(tcaf, nk); for (k = 0; k < nk; k++) { snew(tcaf[k], ncorr); } if (fn_cub) { snew(tcafc, nkc); for (k = 0; k < nkc; k++) { snew(tcafc[k], ncorr); } } snew(sig, ncorr); for (i = 0; i < ncorr; i++) { sig[i] = std::exp(0.5*i*dt/wt); } low_do_autocorr(fn_tca, oenv, "Transverse Current Autocorrelation Functions", nframes, ntc, ncorr, tc, dt, eacNormal, 1, FALSE, FALSE, FALSE, 0, 0, 0); do_view(oenv, fn_tca, "-nxy"); fp = xvgropen(fn_tc, "Transverse Current Autocorrelation Functions", "Time (ps)", "TCAF", oenv); for (i = 0; i < ncorr; i++) { kc = 0; fprintf(fp, "%g", i*dt); for (k = 0; k < nk; k++) { for (j = 0; j < NPK; j++) { tcaf[k][i] += tc[NPK*k+j][i]; } if (fn_cub) { for (j = 0; j < NPK; j++) { tcafc[kc][i] += tc[NPK*k+j][i]; } } if (i == 0) { fprintf(fp, " %g", 1.0); } else { tcaf[k][i] /= tcaf[k][0]; fprintf(fp, " %g", tcaf[k][i]); } if (k+1 == kset_c[kc+1]) { kc++; } } fprintf(fp, "\n"); } xvgrclose(fp); do_view(oenv, fn_tc, "-nxy"); if (fn_cub) { fp_cub = xvgropen(fn_cub, "TCAFs and fits", "Time (ps)", "TCAF", oenv); for (kc = 0; kc < nkc; kc++) { fprintf(fp_cub, "%g %g\n", 0.0, 1.0); for (i = 1; i < ncorr; i++) { tcafc[kc][i] /= tcafc[kc][0]; fprintf(fp_cub, "%g %g\n", i*dt, tcafc[kc][i]); } fprintf(fp_cub, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : ""); tcafc[kc][0] = 1.0; } } fp_vk = xvgropen(fn_vk, "Fits", "k (nm\\S-1\\N)", "\\8h\\4 (10\\S-3\\N kg m\\S-1\\N s\\S-1\\N)", oenv); if (output_env_get_print_xvgr_codes(oenv)) { fprintf(fp_vk, "@ s0 symbol 2\n"); fprintf(fp_vk, "@ s0 symbol color 1\n"); fprintf(fp_vk, "@ s0 linestyle 0\n"); if (fn_cub) { fprintf(fp_vk, "@ s1 symbol 3\n"); fprintf(fp_vk, "@ s1 symbol color 2\n"); } } fp = xvgropen(fn_tcf, "TCAF Fits", "Time (ps)", "", oenv); for (k = 0; k < nk; k++) { tcaf[k][0] = 1.0; fitparms[0] = 1; fitparms[1] = 1; do_lmfit(ncorr, tcaf[k], sig, dt, 0, 0, ncorr*dt, oenv, bDebugMode(), effnVAC, fitparms, 0, NULL); eta = 1000*fitparms[1]*rho/ (4*fitparms[0]*PICO*norm2(kfac[k])/(NANO*NANO)); fprintf(stdout, "k %6.3f tau %6.3f eta %8.5f 10^-3 kg/(m s)\n", norm(kfac[k]), fitparms[0], eta); fprintf(fp_vk, "%6.3f %g\n", norm(kfac[k]), eta); for (i = 0; i < ncorr; i++) { fprintf(fp, "%g %g\n", i*dt, fit_function(effnVAC, fitparms, i*dt)); } fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : ""); } xvgrclose(fp); do_view(oenv, fn_tcf, "-nxy"); if (fn_cub) { fprintf(stdout, "Averaged over k-vectors:\n"); fprintf(fp_vk, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : ""); for (k = 0; k < nkc; k++) { tcafc[k][0] = 1.0; fitparms[0] = 1; fitparms[1] = 1; do_lmfit(ncorr, tcafc[k], sig, dt, 0, 0, ncorr*dt, oenv, bDebugMode(), effnVAC, fitparms, 0, NULL); eta = 1000*fitparms[1]*rho/ (4*fitparms[0]*PICO*norm2(kfac[kset_c[k]])/(NANO*NANO)); fprintf(stdout, "k %6.3f tau %6.3f Omega %6.3f eta %8.5f 10^-3 kg/(m s)\n", norm(kfac[kset_c[k]]), fitparms[0], fitparms[1], eta); fprintf(fp_vk, "%6.3f %g\n", norm(kfac[kset_c[k]]), eta); for (i = 0; i < ncorr; i++) { fprintf(fp_cub, "%g %g\n", i*dt, fit_function(effnVAC, fitparms, i*dt)); } fprintf(fp_cub, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : ""); } fprintf(fp_vk, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : ""); xvgrclose(fp_cub); do_view(oenv, fn_cub, "-nxy"); } xvgrclose(fp_vk); do_view(oenv, fn_vk, "-nxy"); }
real do_lmfit(int ndata, real c1[], real sig[], real dt, real x0[], real begintimefit, real endtimefit, const output_env_t oenv, gmx_bool bVerbose, int eFitFn, real fitparms[], int fix) { FILE *fp; char buf[32]; int i, j, nparm, nfitpnts; real integral, ttt; real *parm, *dparm; real *x, *y, *dy; real ftol = 1e-4; nparm = nfp_ffn[eFitFn]; if (debug) { fprintf(debug, "There are %d points to fit %d vars!\n", ndata, nparm); fprintf(debug, "Fit to function %d from %g through %g, dt=%g\n", eFitFn, begintimefit, endtimefit, dt); } snew(x, ndata); snew(y, ndata); snew(dy, ndata); j = 0; for (i = 0; i < ndata; i++) { ttt = x0 ? x0[i] : dt*i; if (ttt >= begintimefit && ttt <= endtimefit) { x[j] = ttt; y[j] = c1[i]; /* mrqmin does not like sig to be zero */ if (sig[i] < 1.0e-7) { dy[j] = 1.0e-7; } else { dy[j] = sig[i]; } if (debug) { fprintf(debug, "j= %d, i= %d, x= %g, y= %g, dy= %g\n", j, i, x[j], y[j], dy[j]); } j++; } } nfitpnts = j; integral = 0; if (nfitpnts < nparm) { fprintf(stderr, "Not enough data points for fitting!\n"); } else { snew(parm, nparm); snew(dparm, nparm); if (fitparms) { for (i = 0; (i < nparm); i++) { parm[i] = fitparms[i]; } } if (!lmfit_exp(nfitpnts, x, y, dy, ftol, parm, dparm, bVerbose, eFitFn, fix)) { fprintf(stderr, "Fit failed!\n"); } else if (nparm <= 3) { /* Compute the integral from begintimefit */ if (nparm == 3) { integral = (parm[0]*myexp(begintimefit, parm[1], parm[0]) + parm[2]*myexp(begintimefit, 1-parm[1], parm[2])); } else if (nparm == 2) { integral = parm[0]*myexp(begintimefit, parm[1], parm[0]); } else if (nparm == 1) { integral = parm[0]*myexp(begintimefit, 1, parm[0]); } else { gmx_fatal(FARGS, "nparm = %d in file %s, line %d", nparm, __FILE__, __LINE__); } /* Generate THE output */ if (bVerbose) { fprintf(stderr, "FIT: # points used in fit is: %d\n", nfitpnts); fprintf(stderr, "FIT: %21s%21s%21s\n", "parm0 ", "parm1 (ps) ", "parm2 (ps) "); fprintf(stderr, "FIT: ------------------------------------------------------------\n"); fprintf(stderr, "FIT: %8.3g +/- %8.3g%9.4g +/- %8.3g%8.3g +/- %8.3g\n", parm[0], dparm[0], parm[1], dparm[1], parm[2], dparm[2]); fprintf(stderr, "FIT: Integral (calc with fitted function) from %g ps to inf. is: %g\n", begintimefit, integral); sprintf(buf, "test%d.xvg", nfitpnts); fp = xvgropen(buf, "C(t) + Fit to C(t)", "Time (ps)", "C(t)", oenv); fprintf(fp, "# parm0 = %g, parm1 = %g, parm2 = %g\n", parm[0], parm[1], parm[2]); for (j = 0; (j < nfitpnts); j++) { ttt = x0 ? x0[j] : dt*j; fprintf(fp, "%10.5e %10.5e %10.5e\n", ttt, c1[j], fit_function(eFitFn, parm, ttt)); } xvgrclose(fp); } } if (fitparms) { for (i = 0; (i < nparm); i++) { fitparms[i] = parm[i]; } } sfree(parm); sfree(dparm); } sfree(x); sfree(y); sfree(dy); return integral; }
static void estimate_error(char *eefile,int nb_min,int resol,int n,int nset, double *av,double *sig,real **val,real dt, bool bFitAc,bool bSingleExpFit,bool bAllowNegLTCorr) { FILE *fp; int bs,prev_bs,nbs,nb; real spacing,nbr; int s,i,j; double blav,var; char **leg; real *tbs,*ybs,rtmp,dens,*fitsig,twooe,tau1_est,tau_sig; real fitparm[4]; real ee,a,tau1,tau2; if (n < 4) { fprintf(stdout,"The number of points is smaller than 4, can not make an error estimate\n"); return; } fp = xvgropen(eefile,"Error estimates", "Block size (time)","Error estimate"); if (bPrintXvgrCodes()) { fprintf(fp, "@ subtitle \"using block averaging, total time %g (%d points)\"\n", (n-1)*dt,n); } snew(leg,2*nset); xvgr_legend(fp,2*nset,leg); sfree(leg); spacing = pow(2,1.0/resol); snew(tbs,n); snew(ybs,n); snew(fitsig,n); for(s=0; s<nset; s++) { nbs = 0; prev_bs = 0; nbr = nb_min; while (nbr <= n) { bs = n/(int)nbr; if (bs != prev_bs) { nb = n/bs; var = 0; for(i=0; i<nb; i++) { blav=0; for (j=0; j<bs; j++) { blav += val[s][bs*i+j]; } var += sqr(av[s] - blav/bs); } tbs[nbs] = bs*dt; if (sig[s] == 0) { ybs[nbs] = 0; } else { ybs[nbs] = var/(nb*(nb-1.0))*(n*dt)/(sig[s]*sig[s]); } nbs++; } nbr *= spacing; nb = (int)(nbr+0.5); prev_bs = bs; } if (sig[s] == 0) { ee = 0; a = 1; tau1 = 0; tau2 = 0; } else { for(i=0; i<nbs/2; i++) { rtmp = tbs[i]; tbs[i] = tbs[nbs-1-i]; tbs[nbs-1-i] = rtmp; rtmp = ybs[i]; ybs[i] = ybs[nbs-1-i]; ybs[nbs-1-i] = rtmp; } /* The initial slope of the normalized ybs^2 is 1. * For a single exponential autocorrelation: ybs(tau1) = 2/e tau1 * From this we take our initial guess for tau1. */ twooe = 2/exp(1); i = -1; do { i++; tau1_est = tbs[i]; } while (i < nbs - 1 && (ybs[i] > ybs[i+1] || ybs[i] > twooe*tau1_est)); if (ybs[0] > ybs[1]) { fprintf(stdout,"Data set %d has strange time correlations:\n" "the std. error using single points is larger than that of blocks of 2 points\n" "The error estimate might be inaccurate, check the fit\n", s+1); /* Use the total time as tau for the fitting weights */ tau_sig = (n - 1)*dt; } else { tau_sig = tau1_est; } if (debug) { fprintf(debug,"set %d tau1 estimate %f\n",s+1,tau1_est); } /* Generate more or less appropriate sigma's, * also taking the density of points into account. */ for(i=0; i<nbs; i++) { if (i == 0) { dens = tbs[1]/tbs[0] - 1; } else if (i == nbs-1) { dens = tbs[nbs-1]/tbs[nbs-2] - 1; } else { dens = 0.5*(tbs[i+1]/tbs[i-1] - 1); } fitsig[i] = sqrt((tau_sig + tbs[i])/dens); } if (!bSingleExpFit) { fitparm[0] = tau1_est; fitparm[1] = 0.95; /* We set the initial guess for tau2 * to halfway between tau1_est and the total time (on log scale). */ fitparm[2] = sqrt(tau1_est*(n-1)*dt); do_lmfit(nbs,ybs,fitsig,0,tbs,0,dt*n,bDebugMode(),effnERREST,fitparm,0); fitparm[3] = 1-fitparm[1]; } if (bSingleExpFit || fitparm[0]<0 || fitparm[2]<0 || fitparm[1]<0 || (fitparm[1]>1 && !bAllowNegLTCorr) || fitparm[2]>(n-1)*dt) { if (!bSingleExpFit) { if (fitparm[2] > (n-1)*dt) { fprintf(stdout, "Warning: tau2 is longer than the length of the data (%g)\n" " the statistics might be bad\n", (n-1)*dt); } else { fprintf(stdout,"a fitted parameter is negative\n"); } fprintf(stdout,"invalid fit: e.e. %g a %g tau1 %g tau2 %g\n", sig[s]*anal_ee_inf(fitparm,n*dt), fitparm[1],fitparm[0],fitparm[2]); /* Do a fit with tau2 fixed at the total time. * One could also choose any other large value for tau2. */ fitparm[0] = tau1_est; fitparm[1] = 0.95; fitparm[2] = (n-1)*dt; fprintf(stderr,"Will fix tau2 at the total time: %g\n",fitparm[2]); do_lmfit(nbs,ybs,fitsig,0,tbs,0,dt*n,bDebugMode(),effnERREST,fitparm,4); fitparm[3] = 1-fitparm[1]; } if (bSingleExpFit || fitparm[0]<0 || fitparm[1]<0 || (fitparm[1]>1 && !bAllowNegLTCorr)) { if (!bSingleExpFit) { fprintf(stdout,"a fitted parameter is negative\n"); fprintf(stdout,"invalid fit: e.e. %g a %g tau1 %g tau2 %g\n", sig[s]*anal_ee_inf(fitparm,n*dt), fitparm[1],fitparm[0],fitparm[2]); } /* Do a single exponential fit */ fprintf(stderr,"Will use a single exponential fit for set %d\n",s+1); fitparm[0] = tau1_est; fitparm[1] = 1.0; fitparm[2] = 0.0; do_lmfit(nbs,ybs,fitsig,0,tbs,0,dt*n,bDebugMode(),effnERREST,fitparm,6); fitparm[3] = 1-fitparm[1]; } } ee = sig[s]*anal_ee_inf(fitparm,n*dt); a = fitparm[1]; tau1 = fitparm[0]; tau2 = fitparm[2]; } fprintf(stdout,"Set %3d: err.est. %g a %g tau1 %g tau2 %g\n", s+1,ee,a,tau1,tau2); fprintf(fp,"@ legend string %d \"av %f\"\n",2*s,av[s]); fprintf(fp,"@ legend string %d \"ee %6g\"\n", 2*s+1,sig[s]*anal_ee_inf(fitparm,n*dt)); for(i=0; i<nbs; i++) { fprintf(fp,"%g %g %g\n",tbs[i],sig[s]*sqrt(ybs[i]/(n*dt)), sig[s]*sqrt(fit_function(effnERREST,fitparm,tbs[i])/(n*dt))); } if (bFitAc) { int fitlen; real *ac,acint,ac_fit[4]; snew(ac,n); for(i=0; i<n; i++) { ac[i] = val[s][i] - av[s]; if (i > 0) fitsig[i] = sqrt(i); else fitsig[i] = 1; } low_do_autocorr(NULL,NULL,n,1,-1,&ac, dt,eacNormal,1,FALSE,TRUE, FALSE,0,0,effnNONE,0); fitlen = n/nb_min; /* Integrate ACF only up to fitlen/2 to avoid integrating noise */ acint = 0.5*ac[0]; for(i=1; i<=fitlen/2; i++) { acint += ac[i]; } acint *= dt; /* Generate more or less appropriate sigma's */ for(i=0; i<=fitlen; i++) { fitsig[i] = sqrt(acint + dt*i); } ac_fit[0] = 0.5*acint; ac_fit[1] = 0.95; ac_fit[2] = 10*acint; do_lmfit(n/nb_min,ac,fitsig,dt,0,0,fitlen*dt, bDebugMode(),effnEXP3,ac_fit,0); ac_fit[3] = 1 - ac_fit[1]; fprintf(stdout,"Set %3d: ac erest %g a %g tau1 %g tau2 %g\n", s+1,sig[s]*anal_ee_inf(ac_fit,n*dt), ac_fit[1],ac_fit[0],ac_fit[2]); fprintf(fp,"&\n"); for(i=0; i<nbs; i++) { fprintf(fp,"%g %g\n",tbs[i], sig[s]*sqrt(fit_function(effnERREST,ac_fit,tbs[i]))/(n*dt)); } sfree(ac); } if (s < nset-1) { fprintf(fp,"&\n"); } } sfree(fitsig); sfree(ybs); sfree(tbs); fclose(fp); }
real fit_acf(int ncorr,int fitfn,const output_env_t oenv,gmx_bool bVerbose, real tbeginfit,real tendfit,real dt,real c1[],real *fit) { real fitparm[3]; real tStart,tail_corr,sum,sumtot=0,ct_estimate,*sig; int i,j,jmax,nf_int; gmx_bool bPrint; bPrint = bVerbose || bDebugMode(); if (bPrint) printf("COR:\n"); if (tendfit <= 0) tendfit = ncorr*dt; nf_int = min(ncorr,(int)(tendfit/dt)); sum = print_and_integrate(debug,nf_int,dt,c1,NULL,1); /* Estimate the correlation time for better fitting */ ct_estimate = 0.5*c1[0]; for(i=1; (i<ncorr) && (c1[i]>0); i++) ct_estimate += c1[i]; ct_estimate *= dt/c1[0]; if (bPrint) { printf("COR: Correlation time (plain integral from %6.3f to %6.3f ps) = %8.5f ps\n", 0.0,dt*nf_int,sum); printf("COR: Relaxation times are computed as fit to an exponential:\n"); printf("COR: %s\n",longs_ffn[fitfn]); printf("COR: Fit to correlation function from %6.3f ps to %6.3f ps, results in a\n",tbeginfit,min(ncorr*dt,tendfit)); } tStart = 0; if (bPrint) printf("COR:%11s%11s%11s%11s%11s%11s%11s\n", "Fit from","Integral","Tail Value","Sum (ps)"," a1 (ps)", (nfp_ffn[fitfn]>=2) ? " a2 ()" : "", (nfp_ffn[fitfn]>=3) ? " a3 (ps)" : ""); if (tbeginfit > 0) jmax = 3; else jmax = 1; if (fitfn == effnEXP3) { fitparm[0] = 0.002*ncorr*dt; fitparm[1] = 0.95; fitparm[2] = 0.2*ncorr*dt; } else { /* Good initial guess, this increases the probability of convergence */ fitparm[0] = ct_estimate; fitparm[1] = 1.0; fitparm[2] = 1.0; } /* Generate more or less appropriate sigma's */ snew(sig,ncorr); for(i=0; i<ncorr; i++) sig[i] = sqrt(ct_estimate+dt*i); for(j=0; ((j<jmax) && (tStart < tendfit)); j++) { /* Use the previous fitparm as starting values for the next fit */ nf_int = min(ncorr,(int)((tStart+1e-4)/dt)); sum = print_and_integrate(debug,nf_int,dt,c1,NULL,1); tail_corr = do_lmfit(ncorr,c1,sig,dt,NULL,tStart,tendfit,oenv, bDebugMode(),fitfn,fitparm,0); sumtot = sum+tail_corr; if (fit && ((jmax == 1) || (j == 1))) for(i=0; (i<ncorr); i++) fit[i] = fit_function(fitfn,fitparm,i*dt); if (bPrint) { printf("COR:%11.4e%11.4e%11.4e%11.4e",tStart,sum,tail_corr,sumtot); for(i=0; (i<nfp_ffn[fitfn]); i++) printf(" %11.4e",fitparm[i]); printf("\n"); } tStart += tbeginfit; } sfree(sig); return sumtot; }