real fit_acf(int ncorr, int fitfn, const gmx_output_env_t *oenv, gmx_bool bVerbose, real tbeginfit, real tendfit, real dt, real c1[], real *fit) { double fitparm[3]; double tStart, tail_corr, sum, sumtot = 0, c_start, ct_estimate; real *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 = std::min(ncorr, (int)(tendfit/dt)); sum = print_and_integrate(debug, nf_int, dt, c1, NULL, 1); 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", effnDescription(fitfn)); printf("COR: Fit to correlation function from %6.3f ps to %6.3f ps, results in a\n", tbeginfit, std::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)", (effnNparams(fitfn) >= 2) ? " a2 ()" : "", (effnNparams(fitfn) >= 3) ? " a3 (ps)" : ""); } snew(sig, ncorr); if (tbeginfit > 0) { jmax = 3; } else { jmax = 1; } for (j = 0; ((j < jmax) && (tStart < tendfit) && (tStart < ncorr*dt)); j++) { /* Estimate the correlation time for better fitting */ c_start = -1; ct_estimate = 0; for (i = 0; (i < ncorr) && (dt*i < tStart || c1[i] > 0); i++) { if (c_start < 0) { if (dt*i >= tStart) { c_start = c1[i]; ct_estimate = 0.5*c1[i]; } } else { ct_estimate += c1[i]; } } if (c_start > 0) { ct_estimate *= dt/c_start; } else { /* The data is strange, so we need to choose somehting */ ct_estimate = tendfit; } if (debug) { fprintf(debug, "tStart %g ct_estimate: %g\n", tStart, ct_estimate); } if (fitfn == effnEXPEXP) { 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 */ for (i = 0; i < ncorr; i++) { sig[i] = sqrt(ct_estimate+dt*i); } nf_int = std::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, NULL); sumtot = sum+tail_corr; if (fit && ((jmax == 1) || (j == 1))) { double mfp[3]; for (i = 0; (i < 3); i++) { mfp[i] = fitparm[i]; } for (i = 0; (i < ncorr); i++) { fit[i] = lmcurves[fitfn](i*dt, mfp); } } if (bPrint) { printf("COR:%11.4e%11.4e%11.4e%11.4e", tStart, sum, tail_corr, sumtot); for (i = 0; (i < effnNparams(fitfn)); i++) { printf(" %11.4e", fitparm[i]); } printf("\n"); } tStart += tbeginfit; } sfree(sig); return sumtot; }
static void do_fit(FILE *out, int n, gmx_bool bYdy, int ny, real *x0, real **val, int npargs, t_pargs *ppa, const gmx_output_env_t *oenv, const char *fn_fitted) { real *c1 = NULL, *sig = NULL; double *fitparm; real tendfit, tbeginfit; int i, efitfn, nparm; efitfn = get_acffitfn(); nparm = effnNparams(efitfn); fprintf(out, "Will fit to the following function:\n"); fprintf(out, "%s\n", effnDescription(efitfn)); c1 = val[n]; if (bYdy) { c1 = val[n]; sig = val[n+1]; fprintf(out, "Using two columns as y and sigma values\n"); } else { snew(sig, ny); } if (opt2parg_bSet("-beginfit", npargs, ppa)) { tbeginfit = opt2parg_real("-beginfit", npargs, ppa); } else { tbeginfit = x0[0]; } if (opt2parg_bSet("-endfit", npargs, ppa)) { tendfit = opt2parg_real("-endfit", npargs, ppa); } else { tendfit = x0[ny-1]; } snew(fitparm, nparm); switch (efitfn) { case effnEXP1: fitparm[0] = 0.5; break; case effnEXP2: fitparm[0] = 0.5; fitparm[1] = c1[0]; break; case effnEXPEXP: fitparm[0] = 1.0; fitparm[1] = 0.5*c1[0]; fitparm[2] = 10.0; break; case effnEXP5: fitparm[0] = fitparm[2] = 0.5*c1[0]; fitparm[1] = 10; fitparm[3] = 40; fitparm[4] = 0; break; case effnEXP7: fitparm[0] = fitparm[2] = fitparm[4] = 0.33*c1[0]; fitparm[1] = 1; fitparm[3] = 10; fitparm[5] = 100; fitparm[6] = 0; break; case effnEXP9: fitparm[0] = fitparm[2] = fitparm[4] = fitparm[6] = 0.25*c1[0]; fitparm[1] = 0.1; fitparm[3] = 1; fitparm[5] = 10; fitparm[7] = 100; fitparm[8] = 0; break; default: fprintf(out, "Warning: don't know how to initialize the parameters\n"); for (i = 0; (i < nparm); i++) { fitparm[i] = 1; } } fprintf(out, "Starting parameters:\n"); for (i = 0; (i < nparm); i++) { fprintf(out, "a%-2d = %12.5e\n", i+1, fitparm[i]); } if (do_lmfit(ny, c1, sig, 0, x0, tbeginfit, tendfit, oenv, bDebugMode(), efitfn, fitparm, 0, fn_fitted) > 0) { for (i = 0; (i < nparm); i++) { fprintf(out, "a%-2d = %12.5e\n", i+1, fitparm[i]); } } else { fprintf(out, "No solution was found\n"); } }