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.
45 #include "gromacs/legacyheaders/macros.h"
46 #include "gromacs/legacyheaders/typedefs.h"
47 #include "gromacs/utility/smalloc.h"
48 #include "gromacs/fileio/xvgr.h"
49 #include "gromacs/utility/futil.h"
51 #include "gromacs/legacyheaders/names.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/math/vec.h"
56 #define MODE(x) ((mode & (x)) == (x))
60 int nrestart, nout, P, fitfn;
61 gmx_bool bFour, bNormalize;
62 real tbeginfit, tendfit;
65 static gmx_bool bACFinit = FALSE;
72 int sffn2effn(const char **sffn)
77 for (i = 0; i < effnNR; i++)
79 if (sffn[i+1] && strcmp(sffn[0], sffn[i+1]) == 0)
88 static void low_do_four_core(int nfour, int nframes, real c1[], real cfour[],
89 int nCos, gmx_bool bPadding)
98 for (i = 0; (i < nframes); i++)
105 for (i = 0; (i < nframes); i++)
107 cfour[i] = cos(c1[i]);
111 for (i = 0; (i < nframes); i++)
113 cfour[i] = sin(c1[i]);
117 gmx_fatal(FARGS, "nCos = %d, %s %d", nCos, __FILE__, __LINE__);
119 for (; (i < nfour); i++)
127 /* printf("aver = %g\n",aver); */
128 for (i = 0; (i < nframes); i++)
135 correl(cfour-1, cfour-1, nfour, ans-1);
139 for (i = 0; (i < nfour); i++)
141 cfour[i] = ans[i]+sqr(aver);
146 for (i = 0; (i < nfour); i++)
155 static void do_ac_core(int nframes, int nout,
156 real corr[], real c1[], int nrestart,
159 int j, k, j3, jk3, m, n;
165 printf("WARNING: setting number of restarts to 1\n");
171 "Starting do_ac_core: nframes=%d, nout=%d, nrestart=%d,mode=%lu\n",
172 nframes, nout, nrestart, mode);
175 for (j = 0; (j < nout); j++)
180 /* Loop over starting points. */
181 for (j = 0; (j < nframes); j += nrestart)
185 /* Loop over the correlation length for this starting point */
186 for (k = 0; (k < nout) && (j+k < nframes); k++)
190 /* Switch over possible ACF types.
191 * It might be more efficient to put the loops inside the switch,
192 * but this is more clear, and save development time!
196 corr[k] += c1[j]*c1[j+k];
198 else if (MODE(eacCos))
200 /* Compute the cos (phi(t)-phi(t+dt)) */
201 corr[k] += cos(c1[j]-c1[j+k]);
203 else if (MODE(eacIden))
205 /* Check equality (phi(t)==phi(t+dt)) */
206 corr[k] += (c1[j] == c1[j+k]) ? 1 : 0;
208 else if (MODE(eacP1) || MODE(eacP2) || MODE(eacP3))
210 for (m = 0; (m < DIM); m++)
215 cth = cos_angle(xj, xk);
217 if (cth-1.0 > 1.0e-15)
219 printf("j: %d, k: %d, xj:(%g,%g,%g), xk:(%g,%g,%g)\n",
220 j, k, xj[XX], xj[YY], xj[ZZ], xk[XX], xk[YY], xk[ZZ]);
223 corr[k] += LegendreP(cth, mode); /* 1.5*cth*cth-0.5; */
225 else if (MODE(eacRcross))
227 for (m = 0; (m < DIM); m++)
234 corr[k] += iprod(rr, rr);
236 else if (MODE(eacVector))
238 for (m = 0; (m < DIM); m++)
249 gmx_fatal(FARGS, "\nInvalid mode (%d) in do_ac_core", mode);
253 /* Correct for the number of points and copy results to the data array */
254 for (j = 0; (j < nout); j++)
256 n = (nframes-j+(nrestart-1))/nrestart;
261 void normalize_acf(int nout, real corr[])
268 fprintf(debug, "Before normalization\n");
269 for (j = 0; (j < nout); j++)
271 fprintf(debug, "%5d %10f\n", j, corr[j]);
275 /* Normalisation makes that c[0] = 1.0 and that other points are scaled
286 for (j = 0; (j < nout); j++)
293 fprintf(debug, "After normalization\n");
294 for (j = 0; (j < nout); j++)
296 fprintf(debug, "%5d %10f\n", j, corr[j]);
301 void average_acf(gmx_bool bVerbose, int n, int nitem, real **c1)
308 printf("Averaging correlation functions\n");
311 for (j = 0; (j < n); j++)
314 for (i = 0; (i < nitem); i++)
322 void norm_and_scale_vectors(int nframes, real c1[], real scale)
327 for (j = 0; (j < nframes); j++)
331 for (m = 0; (m < DIM); m++)
338 void dump_tmp(char *s, int n, real c[])
343 fp = gmx_ffopen(s, "w");
344 for (i = 0; (i < n); i++)
346 fprintf(fp, "%10d %10g\n", i, c[i]);
351 real print_and_integrate(FILE *fp, int n, real dt, real c[], real *fit, int nskip)
356 /* Use trapezoidal rule for calculating integral */
358 for (j = 0; (j < n); j++)
361 if (fp && (nskip == 0 || j % nskip == 0))
363 fprintf(fp, "%10.3f %10.5f\n", j*dt, c0);
367 sum += dt*(c0+c[j-1]);
375 for (j = 0; (j < n); j++)
377 if (nskip == 0 || j % nskip == 0)
379 fprintf(fp, "%10.3f %10.5f\n", j*dt, fit[j]);
388 real evaluate_integral(int n, real x[], real y[], real dy[], real aver_start,
391 double sum, sum_var, w;
392 double sum_tail = 0, sum2_tail = 0;
393 int j, nsum_tail = 0;
395 /* Use trapezoidal rule for calculating integral */
398 gmx_fatal(FARGS, "Evaluating integral: n = %d (file %s, line %d)",
399 n, __FILE__, __LINE__);
404 for (j = 0; (j < n); j++)
409 w += 0.5*(x[j] - x[j-1]);
413 w += 0.5*(x[j+1] - x[j]);
418 /* Assume all errors are uncorrelated */
419 sum_var += sqr(w*dy[j]);
422 if ((aver_start > 0) && (x[j] >= aver_start))
425 sum2_tail += sqrt(sum_var);
432 sum = sum_tail/nsum_tail;
433 /* This is a worst case estimate, assuming all stddev's are correlated. */
434 *stddev = sum2_tail/nsum_tail;
438 *stddev = sqrt(sum_var);
444 void do_four_core(unsigned long mode, int nfour, int nf2, int nframes,
445 real c1[], real csum[], real ctmp[])
456 /********************************************
458 ********************************************/
459 low_do_four_core(nfour, nf2, c1, csum, enNorm, FALSE);
461 else if (MODE(eacCos))
463 /***************************************************
465 ***************************************************/
466 /* Copy the data to temp array. Since we need it twice
467 * we can't overwrite original.
469 for (j = 0; (j < nf2); j++)
474 /* Cosine term of AC function */
475 low_do_four_core(nfour, nf2, ctmp, cfour, enCos, FALSE);
476 for (j = 0; (j < nf2); j++)
481 /* Sine term of AC function */
482 low_do_four_core(nfour, nf2, ctmp, cfour, enSin, FALSE);
483 for (j = 0; (j < nf2); j++)
489 else if (MODE(eacP2))
491 /***************************************************
492 * Legendre polynomials
493 ***************************************************/
494 /* First normalize the vectors */
495 norm_and_scale_vectors(nframes, c1, 1.0);
497 /* For P2 thingies we have to do six FFT based correls
498 * First for XX^2, then for YY^2, then for ZZ^2
499 * Then we have to do XY, YZ and XZ (counting these twice)
500 * After that we sum them and normalise
501 * P2(x) = (3 * cos^2 (x) - 1)/2
502 * for unit vectors u and v we compute the cosine as the inner product
503 * cos(u,v) = uX vX + uY vY + uZ vZ
507 * C(t) = | (3 cos^2(u(t'),u(t'+t)) - 1)/2 dt'
512 * P2(u(0),u(t)) = [3 * (uX(0) uX(t) +
514 * uZ(0) uZ(t))^2 - 1]/2
515 * = [3 * ((uX(0) uX(t))^2 +
518 * 2(uX(0) uY(0) uX(t) uY(t)) +
519 * 2(uX(0) uZ(0) uX(t) uZ(t)) +
520 * 2(uY(0) uZ(0) uY(t) uZ(t))) - 1]/2
522 * = [(3/2) * (<uX^2> + <uY^2> + <uZ^2> +
523 * 2<uXuY> + 2<uXuZ> + 2<uYuZ>) - 0.5]
527 /* Because of normalization the number of -0.5 to subtract
528 * depends on the number of data points!
530 for (j = 0; (j < nf2); j++)
532 csum[j] = -0.5*(nf2-j);
535 /***** DIAGONAL ELEMENTS ************/
536 for (m = 0; (m < DIM); m++)
538 /* Copy the vector data in a linear array */
539 for (j = 0; (j < nf2); j++)
541 ctmp[j] = sqr(c1[DIM*j+m]);
545 sprintf(buf, "c1diag%d.xvg", m);
546 dump_tmp(buf, nf2, ctmp);
549 low_do_four_core(nfour, nf2, ctmp, cfour, enNorm, FALSE);
553 sprintf(buf, "c1dfout%d.xvg", m);
554 dump_tmp(buf, nf2, cfour);
557 for (j = 0; (j < nf2); j++)
559 csum[j] += fac*(cfour[j]);
562 /******* OFF-DIAGONAL ELEMENTS **********/
563 for (m = 0; (m < DIM); m++)
565 /* Copy the vector data in a linear array */
567 for (j = 0; (j < nf2); j++)
569 ctmp[j] = c1[DIM*j+m]*c1[DIM*j+m1];
574 sprintf(buf, "c1off%d.xvg", m);
575 dump_tmp(buf, nf2, ctmp);
577 low_do_four_core(nfour, nf2, ctmp, cfour, enNorm, FALSE);
580 sprintf(buf, "c1ofout%d.xvg", m);
581 dump_tmp(buf, nf2, cfour);
584 for (j = 0; (j < nf2); j++)
586 csum[j] += fac*cfour[j];
590 else if (MODE(eacP1) || MODE(eacVector))
592 /***************************************************
594 ***************************************************/
597 /* First normalize the vectors */
598 norm_and_scale_vectors(nframes, c1, 1.0);
601 /* For vector thingies we have to do three FFT based correls
602 * First for XX, then for YY, then for ZZ
603 * After that we sum them and normalise
605 for (j = 0; (j < nf2); j++)
609 for (m = 0; (m < DIM); m++)
611 /* Copy the vector data in a linear array */
612 for (j = 0; (j < nf2); j++)
614 ctmp[j] = c1[DIM*j+m];
616 low_do_four_core(nfour, nf2, ctmp, cfour, enNorm, FALSE);
617 for (j = 0; (j < nf2); j++)
625 gmx_fatal(FARGS, "\nUnknown mode in do_autocorr (%d)", mode);
629 for (j = 0; (j < nf2); j++)
631 c1[j] = csum[j]/(real)(nframes-j);
635 real fit_acf(int ncorr, int fitfn, const output_env_t oenv, gmx_bool bVerbose,
636 real tbeginfit, real tendfit, real dt, real c1[], real *fit)
639 real tStart, tail_corr, sum, sumtot = 0, c_start, ct_estimate, *sig;
640 int i, j, jmax, nf_int;
643 bPrint = bVerbose || bDebugMode();
654 nf_int = min(ncorr, (int)(tendfit/dt));
655 sum = print_and_integrate(debug, nf_int, dt, c1, NULL, 1);
657 /* Estimate the correlation time for better fitting */
658 ct_estimate = 0.5*c1[0];
659 for (i = 1; (i < ncorr) && (c1[i] > 0); i++)
661 ct_estimate += c1[i];
663 ct_estimate *= dt/c1[0];
667 printf("COR: Correlation time (plain integral from %6.3f to %6.3f ps) = %8.5f ps\n",
668 0.0, dt*nf_int, sum);
669 printf("COR: Relaxation times are computed as fit to an exponential:\n");
670 printf("COR: %s\n", longs_ffn[fitfn]);
671 printf("COR: Fit to correlation function from %6.3f ps to %6.3f ps, results in a\n", tbeginfit, min(ncorr*dt, tendfit));
677 printf("COR:%11s%11s%11s%11s%11s%11s%11s\n",
678 "Fit from", "Integral", "Tail Value", "Sum (ps)", " a1 (ps)",
679 (nfp_ffn[fitfn] >= 2) ? " a2 ()" : "",
680 (nfp_ffn[fitfn] >= 3) ? " a3 (ps)" : "");
693 for (j = 0; ((j < jmax) && (tStart < tendfit) && (tStart < ncorr*dt)); j++)
695 /* Estimate the correlation time for better fitting */
698 for (i = 0; (i < ncorr) && (dt*i < tStart || c1[i] > 0); i++)
705 ct_estimate = 0.5*c1[i];
710 ct_estimate += c1[i];
715 ct_estimate *= dt/c_start;
719 /* The data is strange, so we need to choose somehting */
720 ct_estimate = tendfit;
724 fprintf(debug, "tStart %g ct_estimate: %g\n", tStart, ct_estimate);
727 if (fitfn == effnEXP3)
729 fitparm[0] = 0.002*ncorr*dt;
731 fitparm[2] = 0.2*ncorr*dt;
735 /* Good initial guess, this increases the probability of convergence */
736 fitparm[0] = ct_estimate;
741 /* Generate more or less appropriate sigma's */
742 for (i = 0; i < ncorr; i++)
744 sig[i] = sqrt(ct_estimate+dt*i);
747 nf_int = min(ncorr, (int)((tStart+1e-4)/dt));
748 sum = print_and_integrate(debug, nf_int, dt, c1, NULL, 1);
749 tail_corr = do_lmfit(ncorr, c1, sig, dt, NULL, tStart, tendfit, oenv,
750 bDebugMode(), fitfn, fitparm, 0);
751 sumtot = sum+tail_corr;
752 if (fit && ((jmax == 1) || (j == 1)))
754 for (i = 0; (i < ncorr); i++)
756 fit[i] = fit_function(fitfn, fitparm, i*dt);
761 printf("COR:%11.4e%11.4e%11.4e%11.4e", tStart, sum, tail_corr, sumtot);
762 for (i = 0; (i < nfp_ffn[fitfn]); i++)
764 printf(" %11.4e", fitparm[i]);
775 void low_do_autocorr(const char *fn, const output_env_t oenv, const char *title,
776 int nframes, int nitem, int nout, real **c1,
777 real dt, unsigned long mode, int nrestart,
778 gmx_bool bAver, gmx_bool bNormalize,
779 gmx_bool bVerbose, real tbeginfit, real tendfit,
782 FILE *fp, *gp = NULL;
786 real c0, sum, Ct2av, Ctav;
787 gmx_bool bFour = acf.bFour;
789 /* Check flags and parameters */
790 nout = get_acfnout();
793 nout = acf.nout = (nframes+1)/2;
795 else if (nout > nframes)
800 if (MODE(eacCos) && MODE(eacVector))
802 gmx_fatal(FARGS, "Incompatible options bCos && bVector (%s, %d)",
805 if ((MODE(eacP3) || MODE(eacRcross)) && bFour)
807 fprintf(stderr, "Can't combine mode %lu with FFT, turning off FFT\n", mode);
810 if (MODE(eacNormal) && MODE(eacVector))
812 gmx_fatal(FARGS, "Incompatible mode bits: normal and vector (or Legendre)");
815 /* Print flags and parameters */
818 printf("Will calculate %s of %d thingies for %d frames\n",
819 title ? title : "autocorrelation", nitem, nframes);
820 printf("bAver = %s, bFour = %s bNormalize= %s\n",
821 bool_names[bAver], bool_names[bFour], bool_names[bNormalize]);
822 printf("mode = %lu, dt = %g, nrestart = %d\n", mode, dt, nrestart);
826 c0 = log((double)nframes)/log(2.0);
836 fprintf(debug, "Using FFT to calculate %s, #points for FFT = %d\n",
840 /* Allocate temp arrays */
846 nfour = 0; /* To keep the compiler happy */
851 /* Loop over items (e.g. molecules or dihedrals)
852 * In this loop the actual correlation functions are computed, but without
855 k = max(1, pow(10, (int)(log(nitem)/log(100))));
856 for (i = 0; i < nitem; i++)
858 if (bVerbose && ((i%k == 0 || i == nitem-1)))
860 fprintf(stderr, "\rThingie %d", i+1);
865 do_four_core(mode, nfour, nframes, nframes, c1[i], csum, ctmp);
869 do_ac_core(nframes, nout, ctmp, c1[i], nrestart, mode);
874 fprintf(stderr, "\n");
882 fp = xvgropen(fn, title, "Time (ps)", "C(t)", oenv);
893 average_acf(bVerbose, nframes, nitem, c1);
898 normalize_acf(nout, c1[0]);
901 if (eFitFn != effnNONE)
903 fit_acf(nout, eFitFn, oenv, fn != NULL, tbeginfit, tendfit, dt, c1[0], fit);
904 sum = print_and_integrate(fp, nout, dt, c1[0], fit, 1);
908 sum = print_and_integrate(fp, nout, dt, c1[0], NULL, 1);
911 printf("Correlation time (integral over corrfn): %g (ps)\n", sum);
917 /* Not averaging. Normalize individual ACFs */
921 gp = xvgropen("ct-distr.xvg", "Correlation times", "item", "time (ps)", oenv);
923 for (i = 0; i < nitem; i++)
927 normalize_acf(nout, c1[i]);
929 if (eFitFn != effnNONE)
931 fit_acf(nout, eFitFn, oenv, fn != NULL, tbeginfit, tendfit, dt, c1[i], fit);
932 sum = print_and_integrate(fp, nout, dt, c1[i], fit, 1);
936 sum = print_and_integrate(fp, nout, dt, c1[i], NULL, 1);
940 "CORRelation time (integral over corrfn %d): %g (ps)\n",
948 fprintf(gp, "%5d %.3f\n", i, sum);
959 printf("Average correlation time %.3f Std. Dev. %.3f Error %.3f (ps)\n",
960 Ctav, sqrt((Ct2av - sqr(Ctav))),
961 sqrt((Ct2av - sqr(Ctav))/(nitem-1)));
971 static const char *Leg[] = { NULL, "0", "1", "2", "3", NULL };
973 t_pargs *add_acf_pargs(int *npargs, t_pargs *pa)
976 { "-acflen", FALSE, etINT, {&acf.nout},
977 "Length of the ACF, default is half the number of frames" },
978 { "-normalize", FALSE, etBOOL, {&acf.bNormalize},
980 { "-fftcorr", FALSE, etBOOL, {&acf.bFour},
981 "HIDDENUse fast fourier transform for correlation function" },
982 { "-nrestart", FALSE, etINT, {&acf.nrestart},
983 "HIDDENNumber of frames between time origins for ACF when no FFT is used" },
984 { "-P", FALSE, etENUM, {Leg},
985 "Order of Legendre polynomial for ACF (0 indicates none)" },
986 { "-fitfn", FALSE, etENUM, {s_ffn},
988 { "-beginfit", FALSE, etREAL, {&acf.tbeginfit},
989 "Time where to begin the exponential fit of the correlation function" },
990 { "-endfit", FALSE, etREAL, {&acf.tendfit},
991 "Time where to end the exponential fit of the correlation function, -1 is until the end" },
993 #define NPA asize(acfpa)
997 snew(ppa, *npargs+NPA);
998 for (i = 0; (i < *npargs); i++)
1002 for (i = 0; (i < NPA); i++)
1004 ppa[*npargs+i] = acfpa[i];
1012 acf.fitfn = effnEXP1;
1014 acf.bNormalize = TRUE;
1015 acf.tbeginfit = 0.0;
1023 void do_autocorr(const char *fn, const output_env_t oenv, const char *title,
1024 int nframes, int nitem, real **c1,
1025 real dt, unsigned long mode, gmx_bool bAver)
1029 printf("ACF data structures have not been initialised. Call add_acf_pargs\n");
1032 /* Handle enumerated types */
1033 sscanf(Leg[0], "%d", &acf.P);
1034 acf.fitfn = sffn2effn(s_ffn);
1039 mode = mode | eacP1;
1042 mode = mode | eacP2;
1045 mode = mode | eacP3;
1051 low_do_autocorr(fn, oenv, title, nframes, nitem, acf.nout, c1, dt, mode,
1052 acf.nrestart, bAver, acf.bNormalize,
1053 bDebugMode(), acf.tbeginfit, acf.tendfit,
1057 int get_acfnout(void)
1061 gmx_fatal(FARGS, "ACF data not initialized yet");
1067 int get_acffitfn(void)
1071 gmx_fatal(FARGS, "ACF data not initialized yet");
1074 return sffn2effn(s_ffn);
1077 void cross_corr(int n, real f[], real g[], real corr[])
1083 correl(f-1, g-1, n, corr-1);
1087 for (i = 0; (i < n); i++)
1089 for (j = i; (j < n); j++)
1091 corr[j-i] += f[i]*g[j];
1093 corr[i] /= (real)(n-i);