2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
48 #include "gromacs/utility/smalloc.h"
50 #include "gromacs/fileio/futil.h"
53 #include "gmx_fatal.h"
57 #define MODE(x) ((mode & (x)) == (x))
61 int nrestart, nout, P, fitfn;
62 gmx_bool bFour, bNormalize;
63 real tbeginfit, tendfit;
66 static gmx_bool bACFinit = FALSE;
73 int sffn2effn(const char **sffn)
78 for (i = 0; i < effnNR; i++)
80 if (sffn[i+1] && strcmp(sffn[0], sffn[i+1]) == 0)
89 static void low_do_four_core(int nfour, int nframes, real c1[], real cfour[],
90 int nCos, gmx_bool bPadding)
99 for (i = 0; (i < nframes); i++)
106 for (i = 0; (i < nframes); i++)
108 cfour[i] = cos(c1[i]);
112 for (i = 0; (i < nframes); i++)
114 cfour[i] = sin(c1[i]);
118 gmx_fatal(FARGS, "nCos = %d, %s %d", nCos, __FILE__, __LINE__);
120 for (; (i < nfour); i++)
128 /* printf("aver = %g\n",aver); */
129 for (i = 0; (i < nframes); i++)
136 correl(cfour-1, cfour-1, nfour, ans-1);
140 for (i = 0; (i < nfour); i++)
142 cfour[i] = ans[i]+sqr(aver);
147 for (i = 0; (i < nfour); i++)
156 static void do_ac_core(int nframes, int nout,
157 real corr[], real c1[], int nrestart,
160 int j, k, j3, jk3, m, n;
166 printf("WARNING: setting number of restarts to 1\n");
172 "Starting do_ac_core: nframes=%d, nout=%d, nrestart=%d,mode=%lu\n",
173 nframes, nout, nrestart, mode);
176 for (j = 0; (j < nout); j++)
181 /* Loop over starting points. */
182 for (j = 0; (j < nframes); j += nrestart)
186 /* Loop over the correlation length for this starting point */
187 for (k = 0; (k < nout) && (j+k < nframes); k++)
191 /* Switch over possible ACF types.
192 * It might be more efficient to put the loops inside the switch,
193 * but this is more clear, and save development time!
197 corr[k] += c1[j]*c1[j+k];
199 else if (MODE(eacCos))
201 /* Compute the cos (phi(t)-phi(t+dt)) */
202 corr[k] += cos(c1[j]-c1[j+k]);
204 else if (MODE(eacIden))
206 /* Check equality (phi(t)==phi(t+dt)) */
207 corr[k] += (c1[j] == c1[j+k]) ? 1 : 0;
209 else if (MODE(eacP1) || MODE(eacP2) || MODE(eacP3))
211 for (m = 0; (m < DIM); m++)
216 cth = cos_angle(xj, xk);
218 if (cth-1.0 > 1.0e-15)
220 printf("j: %d, k: %d, xj:(%g,%g,%g), xk:(%g,%g,%g)\n",
221 j, k, xj[XX], xj[YY], xj[ZZ], xk[XX], xk[YY], xk[ZZ]);
224 corr[k] += LegendreP(cth, mode); /* 1.5*cth*cth-0.5; */
226 else if (MODE(eacRcross))
228 for (m = 0; (m < DIM); m++)
235 corr[k] += iprod(rr, rr);
237 else if (MODE(eacVector))
239 for (m = 0; (m < DIM); m++)
250 gmx_fatal(FARGS, "\nInvalid mode (%d) in do_ac_core", mode);
254 /* Correct for the number of points and copy results to the data array */
255 for (j = 0; (j < nout); j++)
257 n = (nframes-j+(nrestart-1))/nrestart;
262 void normalize_acf(int nout, real corr[])
269 fprintf(debug, "Before normalization\n");
270 for (j = 0; (j < nout); j++)
272 fprintf(debug, "%5d %10f\n", j, corr[j]);
276 /* Normalisation makes that c[0] = 1.0 and that other points are scaled
287 for (j = 0; (j < nout); j++)
294 fprintf(debug, "After normalization\n");
295 for (j = 0; (j < nout); j++)
297 fprintf(debug, "%5d %10f\n", j, corr[j]);
302 void average_acf(gmx_bool bVerbose, int n, int nitem, real **c1)
309 printf("Averaging correlation functions\n");
312 for (j = 0; (j < n); j++)
315 for (i = 0; (i < nitem); i++)
323 void norm_and_scale_vectors(int nframes, real c1[], real scale)
328 for (j = 0; (j < nframes); j++)
332 for (m = 0; (m < DIM); m++)
339 void dump_tmp(char *s, int n, real c[])
344 fp = gmx_ffopen(s, "w");
345 for (i = 0; (i < n); i++)
347 fprintf(fp, "%10d %10g\n", i, c[i]);
352 real print_and_integrate(FILE *fp, int n, real dt, real c[], real *fit, int nskip)
357 /* Use trapezoidal rule for calculating integral */
359 for (j = 0; (j < n); j++)
362 if (fp && (nskip == 0 || j % nskip == 0))
364 fprintf(fp, "%10.3f %10.5f\n", j*dt, c0);
368 sum += dt*(c0+c[j-1]);
376 for (j = 0; (j < n); j++)
378 if (nskip == 0 || j % nskip == 0)
380 fprintf(fp, "%10.3f %10.5f\n", j*dt, fit[j]);
389 real evaluate_integral(int n, real x[], real y[], real dy[], real aver_start,
392 double sum, sum_var, w;
393 double sum_tail = 0, sum2_tail = 0;
394 int j, nsum_tail = 0;
396 /* Use trapezoidal rule for calculating integral */
399 gmx_fatal(FARGS, "Evaluating integral: n = %d (file %s, line %d)",
400 n, __FILE__, __LINE__);
405 for (j = 0; (j < n); j++)
410 w += 0.5*(x[j] - x[j-1]);
414 w += 0.5*(x[j+1] - x[j]);
419 /* Assume all errors are uncorrelated */
420 sum_var += sqr(w*dy[j]);
423 if ((aver_start > 0) && (x[j] >= aver_start))
426 sum2_tail += sqrt(sum_var);
433 sum = sum_tail/nsum_tail;
434 /* This is a worst case estimate, assuming all stddev's are correlated. */
435 *stddev = sum2_tail/nsum_tail;
439 *stddev = sqrt(sum_var);
445 void do_four_core(unsigned long mode, int nfour, int nf2, int nframes,
446 real c1[], real csum[], real ctmp[])
457 /********************************************
459 ********************************************/
460 low_do_four_core(nfour, nf2, c1, csum, enNorm, FALSE);
462 else if (MODE(eacCos))
464 /***************************************************
466 ***************************************************/
467 /* Copy the data to temp array. Since we need it twice
468 * we can't overwrite original.
470 for (j = 0; (j < nf2); j++)
475 /* Cosine term of AC function */
476 low_do_four_core(nfour, nf2, ctmp, cfour, enCos, FALSE);
477 for (j = 0; (j < nf2); j++)
482 /* Sine term of AC function */
483 low_do_four_core(nfour, nf2, ctmp, cfour, enSin, FALSE);
484 for (j = 0; (j < nf2); j++)
490 else if (MODE(eacP2))
492 /***************************************************
493 * Legendre polynomials
494 ***************************************************/
495 /* First normalize the vectors */
496 norm_and_scale_vectors(nframes, c1, 1.0);
498 /* For P2 thingies we have to do six FFT based correls
499 * First for XX^2, then for YY^2, then for ZZ^2
500 * Then we have to do XY, YZ and XZ (counting these twice)
501 * After that we sum them and normalise
502 * P2(x) = (3 * cos^2 (x) - 1)/2
503 * for unit vectors u and v we compute the cosine as the inner product
504 * cos(u,v) = uX vX + uY vY + uZ vZ
508 * C(t) = | (3 cos^2(u(t'),u(t'+t)) - 1)/2 dt'
513 * P2(u(0),u(t)) = [3 * (uX(0) uX(t) +
515 * uZ(0) uZ(t))^2 - 1]/2
516 * = [3 * ((uX(0) uX(t))^2 +
519 * 2(uX(0) uY(0) uX(t) uY(t)) +
520 * 2(uX(0) uZ(0) uX(t) uZ(t)) +
521 * 2(uY(0) uZ(0) uY(t) uZ(t))) - 1]/2
523 * = [(3/2) * (<uX^2> + <uY^2> + <uZ^2> +
524 * 2<uXuY> + 2<uXuZ> + 2<uYuZ>) - 0.5]
528 /* Because of normalization the number of -0.5 to subtract
529 * depends on the number of data points!
531 for (j = 0; (j < nf2); j++)
533 csum[j] = -0.5*(nf2-j);
536 /***** DIAGONAL ELEMENTS ************/
537 for (m = 0; (m < DIM); m++)
539 /* Copy the vector data in a linear array */
540 for (j = 0; (j < nf2); j++)
542 ctmp[j] = sqr(c1[DIM*j+m]);
546 sprintf(buf, "c1diag%d.xvg", m);
547 dump_tmp(buf, nf2, ctmp);
550 low_do_four_core(nfour, nf2, ctmp, cfour, enNorm, FALSE);
554 sprintf(buf, "c1dfout%d.xvg", m);
555 dump_tmp(buf, nf2, cfour);
558 for (j = 0; (j < nf2); j++)
560 csum[j] += fac*(cfour[j]);
563 /******* OFF-DIAGONAL ELEMENTS **********/
564 for (m = 0; (m < DIM); m++)
566 /* Copy the vector data in a linear array */
568 for (j = 0; (j < nf2); j++)
570 ctmp[j] = c1[DIM*j+m]*c1[DIM*j+m1];
575 sprintf(buf, "c1off%d.xvg", m);
576 dump_tmp(buf, nf2, ctmp);
578 low_do_four_core(nfour, nf2, ctmp, cfour, enNorm, FALSE);
581 sprintf(buf, "c1ofout%d.xvg", m);
582 dump_tmp(buf, nf2, cfour);
585 for (j = 0; (j < nf2); j++)
587 csum[j] += fac*cfour[j];
591 else if (MODE(eacP1) || MODE(eacVector))
593 /***************************************************
595 ***************************************************/
598 /* First normalize the vectors */
599 norm_and_scale_vectors(nframes, c1, 1.0);
602 /* For vector thingies we have to do three FFT based correls
603 * First for XX, then for YY, then for ZZ
604 * After that we sum them and normalise
606 for (j = 0; (j < nf2); j++)
610 for (m = 0; (m < DIM); m++)
612 /* Copy the vector data in a linear array */
613 for (j = 0; (j < nf2); j++)
615 ctmp[j] = c1[DIM*j+m];
617 low_do_four_core(nfour, nf2, ctmp, cfour, enNorm, FALSE);
618 for (j = 0; (j < nf2); j++)
626 gmx_fatal(FARGS, "\nUnknown mode in do_autocorr (%d)", mode);
630 for (j = 0; (j < nf2); j++)
632 c1[j] = csum[j]/(real)(nframes-j);
636 real fit_acf(int ncorr, int fitfn, const output_env_t oenv, gmx_bool bVerbose,
637 real tbeginfit, real tendfit, real dt, real c1[], real *fit)
640 real tStart, tail_corr, sum, sumtot = 0, c_start, ct_estimate, *sig;
641 int i, j, jmax, nf_int;
644 bPrint = bVerbose || bDebugMode();
655 nf_int = min(ncorr, (int)(tendfit/dt));
656 sum = print_and_integrate(debug, nf_int, dt, c1, NULL, 1);
658 /* Estimate the correlation time for better fitting */
659 ct_estimate = 0.5*c1[0];
660 for (i = 1; (i < ncorr) && (c1[i] > 0); i++)
662 ct_estimate += c1[i];
664 ct_estimate *= dt/c1[0];
668 printf("COR: Correlation time (plain integral from %6.3f to %6.3f ps) = %8.5f ps\n",
669 0.0, dt*nf_int, sum);
670 printf("COR: Relaxation times are computed as fit to an exponential:\n");
671 printf("COR: %s\n", longs_ffn[fitfn]);
672 printf("COR: Fit to correlation function from %6.3f ps to %6.3f ps, results in a\n", tbeginfit, min(ncorr*dt, tendfit));
678 printf("COR:%11s%11s%11s%11s%11s%11s%11s\n",
679 "Fit from", "Integral", "Tail Value", "Sum (ps)", " a1 (ps)",
680 (nfp_ffn[fitfn] >= 2) ? " a2 ()" : "",
681 (nfp_ffn[fitfn] >= 3) ? " a3 (ps)" : "");
694 for (j = 0; ((j < jmax) && (tStart < tendfit) && (tStart < ncorr*dt)); j++)
696 /* Estimate the correlation time for better fitting */
699 for (i = 0; (i < ncorr) && (dt*i < tStart || c1[i] > 0); i++)
706 ct_estimate = 0.5*c1[i];
711 ct_estimate += c1[i];
716 ct_estimate *= dt/c_start;
720 /* The data is strange, so we need to choose somehting */
721 ct_estimate = tendfit;
725 fprintf(debug, "tStart %g ct_estimate: %g\n", tStart, ct_estimate);
728 if (fitfn == effnEXP3)
730 fitparm[0] = 0.002*ncorr*dt;
732 fitparm[2] = 0.2*ncorr*dt;
736 /* Good initial guess, this increases the probability of convergence */
737 fitparm[0] = ct_estimate;
742 /* Generate more or less appropriate sigma's */
743 for (i = 0; i < ncorr; i++)
745 sig[i] = sqrt(ct_estimate+dt*i);
748 nf_int = min(ncorr, (int)((tStart+1e-4)/dt));
749 sum = print_and_integrate(debug, nf_int, dt, c1, NULL, 1);
750 tail_corr = do_lmfit(ncorr, c1, sig, dt, NULL, tStart, tendfit, oenv,
751 bDebugMode(), fitfn, fitparm, 0);
752 sumtot = sum+tail_corr;
753 if (fit && ((jmax == 1) || (j == 1)))
755 for (i = 0; (i < ncorr); i++)
757 fit[i] = fit_function(fitfn, fitparm, i*dt);
762 printf("COR:%11.4e%11.4e%11.4e%11.4e", tStart, sum, tail_corr, sumtot);
763 for (i = 0; (i < nfp_ffn[fitfn]); i++)
765 printf(" %11.4e", fitparm[i]);
776 void low_do_autocorr(const char *fn, const output_env_t oenv, const char *title,
777 int nframes, int nitem, int nout, real **c1,
778 real dt, unsigned long mode, int nrestart,
779 gmx_bool bAver, gmx_bool bNormalize,
780 gmx_bool bVerbose, real tbeginfit, real tendfit,
783 FILE *fp, *gp = NULL;
787 real c0, sum, Ct2av, Ctav;
788 gmx_bool bFour = acf.bFour;
790 /* Check flags and parameters */
791 nout = get_acfnout();
794 nout = acf.nout = (nframes+1)/2;
796 else if (nout > nframes)
801 if (MODE(eacCos) && MODE(eacVector))
803 gmx_fatal(FARGS, "Incompatible options bCos && bVector (%s, %d)",
806 if ((MODE(eacP3) || MODE(eacRcross)) && bFour)
808 fprintf(stderr, "Can't combine mode %lu with FFT, turning off FFT\n", mode);
811 if (MODE(eacNormal) && MODE(eacVector))
813 gmx_fatal(FARGS, "Incompatible mode bits: normal and vector (or Legendre)");
816 /* Print flags and parameters */
819 printf("Will calculate %s of %d thingies for %d frames\n",
820 title ? title : "autocorrelation", nitem, nframes);
821 printf("bAver = %s, bFour = %s bNormalize= %s\n",
822 bool_names[bAver], bool_names[bFour], bool_names[bNormalize]);
823 printf("mode = %lu, dt = %g, nrestart = %d\n", mode, dt, nrestart);
827 c0 = log((double)nframes)/log(2.0);
837 fprintf(debug, "Using FFT to calculate %s, #points for FFT = %d\n",
841 /* Allocate temp arrays */
847 nfour = 0; /* To keep the compiler happy */
852 /* Loop over items (e.g. molecules or dihedrals)
853 * In this loop the actual correlation functions are computed, but without
856 k = max(1, pow(10, (int)(log(nitem)/log(100))));
857 for (i = 0; i < nitem; i++)
859 if (bVerbose && ((i%k == 0 || i == nitem-1)))
861 fprintf(stderr, "\rThingie %d", i+1);
866 do_four_core(mode, nfour, nframes, nframes, c1[i], csum, ctmp);
870 do_ac_core(nframes, nout, ctmp, c1[i], nrestart, mode);
875 fprintf(stderr, "\n");
883 fp = xvgropen(fn, title, "Time (ps)", "C(t)", oenv);
894 average_acf(bVerbose, nframes, nitem, c1);
899 normalize_acf(nout, c1[0]);
902 if (eFitFn != effnNONE)
904 fit_acf(nout, eFitFn, oenv, fn != NULL, tbeginfit, tendfit, dt, c1[0], fit);
905 sum = print_and_integrate(fp, nout, dt, c1[0], fit, 1);
909 sum = print_and_integrate(fp, nout, dt, c1[0], NULL, 1);
912 printf("Correlation time (integral over corrfn): %g (ps)\n", sum);
918 /* Not averaging. Normalize individual ACFs */
922 gp = xvgropen("ct-distr.xvg", "Correlation times", "item", "time (ps)", oenv);
924 for (i = 0; i < nitem; i++)
928 normalize_acf(nout, c1[i]);
930 if (eFitFn != effnNONE)
932 fit_acf(nout, eFitFn, oenv, fn != NULL, tbeginfit, tendfit, dt, c1[i], fit);
933 sum = print_and_integrate(fp, nout, dt, c1[i], fit, 1);
937 sum = print_and_integrate(fp, nout, dt, c1[i], NULL, 1);
941 "CORRelation time (integral over corrfn %d): %g (ps)\n",
949 fprintf(gp, "%5d %.3f\n", i, sum);
960 printf("Average correlation time %.3f Std. Dev. %.3f Error %.3f (ps)\n",
961 Ctav, sqrt((Ct2av - sqr(Ctav))),
962 sqrt((Ct2av - sqr(Ctav))/(nitem-1)));
972 static const char *Leg[] = { NULL, "0", "1", "2", "3", NULL };
974 t_pargs *add_acf_pargs(int *npargs, t_pargs *pa)
977 { "-acflen", FALSE, etINT, {&acf.nout},
978 "Length of the ACF, default is half the number of frames" },
979 { "-normalize", FALSE, etBOOL, {&acf.bNormalize},
981 { "-fftcorr", FALSE, etBOOL, {&acf.bFour},
982 "HIDDENUse fast fourier transform for correlation function" },
983 { "-nrestart", FALSE, etINT, {&acf.nrestart},
984 "HIDDENNumber of frames between time origins for ACF when no FFT is used" },
985 { "-P", FALSE, etENUM, {Leg},
986 "Order of Legendre polynomial for ACF (0 indicates none)" },
987 { "-fitfn", FALSE, etENUM, {s_ffn},
989 { "-beginfit", FALSE, etREAL, {&acf.tbeginfit},
990 "Time where to begin the exponential fit of the correlation function" },
991 { "-endfit", FALSE, etREAL, {&acf.tendfit},
992 "Time where to end the exponential fit of the correlation function, -1 is until the end" },
994 #define NPA asize(acfpa)
998 snew(ppa, *npargs+NPA);
999 for (i = 0; (i < *npargs); i++)
1003 for (i = 0; (i < NPA); i++)
1005 ppa[*npargs+i] = acfpa[i];
1013 acf.fitfn = effnEXP1;
1015 acf.bNormalize = TRUE;
1016 acf.tbeginfit = 0.0;
1024 void do_autocorr(const char *fn, const output_env_t oenv, const char *title,
1025 int nframes, int nitem, real **c1,
1026 real dt, unsigned long mode, gmx_bool bAver)
1030 printf("ACF data structures have not been initialised. Call add_acf_pargs\n");
1033 /* Handle enumerated types */
1034 sscanf(Leg[0], "%d", &acf.P);
1035 acf.fitfn = sffn2effn(s_ffn);
1040 mode = mode | eacP1;
1043 mode = mode | eacP2;
1046 mode = mode | eacP3;
1052 low_do_autocorr(fn, oenv, title, nframes, nitem, acf.nout, c1, dt, mode,
1053 acf.nrestart, bAver, acf.bNormalize,
1054 bDebugMode(), acf.tbeginfit, acf.tendfit,
1058 int get_acfnout(void)
1062 gmx_fatal(FARGS, "ACF data not initialized yet");
1068 int get_acffitfn(void)
1072 gmx_fatal(FARGS, "ACF data not initialized yet");
1075 return sffn2effn(s_ffn);
1078 void cross_corr(int n, real f[], real g[], real corr[])
1084 correl(f-1, g-1, n, corr-1);
1088 for (i = 0; (i < n); i++)
1090 for (j = i; (j < n); j++)
1092 corr[j-i] += f[i]*g[j];
1094 corr[i] /= (real)(n-i);