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,2015,2016,2017,2018, 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.
39 * Implements function to compute many autocorrelation functions
41 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42 * \ingroup module_correlationfunctions
54 #include "gromacs/correlationfunctions/expfit.h"
55 #include "gromacs/correlationfunctions/integrate.h"
56 #include "gromacs/correlationfunctions/manyautocorrelation.h"
57 #include "gromacs/correlationfunctions/polynomials.h"
58 #include "gromacs/fileio/xvgr.h"
59 #include "gromacs/math/functions.h"
60 #include "gromacs/math/vec.h"
61 #include "gromacs/utility/arraysize.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/futil.h"
64 #include "gromacs/utility/real.h"
65 #include "gromacs/utility/smalloc.h"
66 #include "gromacs/utility/strconvert.h"
68 /*! \brief Shortcut macro to select modes. */
69 #define MODE(x) ((mode & (x)) == (x))
73 int nrestart, nout, P, fitfn;
74 gmx_bool bFour, bNormalize;
75 real tbeginfit, tendfit;
78 /*! \brief Global variable set true if initialization routines are called. */
79 static gmx_bool bACFinit = FALSE;
81 /*! \brief Data structure for storing command line variables. */
88 /*! \brief Routine to compute ACF using FFT. */
89 static void low_do_four_core(int nframes, real c1[], real cfour[],
93 std::vector<std::vector<real> > data;
95 data[0].resize(nframes, 0);
99 for (i = 0; (i < nframes); i++)
105 for (i = 0; (i < nframes); i++)
107 data[0][i] = cos(c1[i]);
111 for (i = 0; (i < nframes); i++)
113 data[0][i] = sin(c1[i]);
117 gmx_fatal(FARGS, "nCos = %d, %s %d", nCos, __FILE__, __LINE__);
120 many_auto_correl(&data);
121 for (i = 0; (i < nframes); i++)
123 cfour[i] = data[0][i];
127 /*! \brief Routine to comput ACF without FFT. */
128 static void do_ac_core(int nframes, int nout,
129 real corr[], real c1[], int nrestart,
132 int j, k, j3, jk3, m, n;
138 printf("WARNING: setting number of restarts to 1\n");
144 "Starting do_ac_core: nframes=%d, nout=%d, nrestart=%d,mode=%lu\n",
145 nframes, nout, nrestart, mode);
148 for (j = 0; (j < nout); j++)
153 /* Loop over starting points. */
154 for (j = 0; (j < nframes); j += nrestart)
158 /* Loop over the correlation length for this starting point */
159 for (k = 0; (k < nout) && (j+k < nframes); k++)
163 /* Switch over possible ACF types.
164 * It might be more efficient to put the loops inside the switch,
165 * but this is more clear, and save development time!
169 corr[k] += c1[j]*c1[j+k];
171 else if (MODE(eacCos))
173 /* Compute the cos (phi(t)-phi(t+dt)) */
174 corr[k] += std::cos(c1[j]-c1[j+k]);
176 else if (MODE(eacIden))
178 /* Check equality (phi(t)==phi(t+dt)) */
179 corr[k] += (c1[j] == c1[j+k]) ? 1 : 0;
181 else if (MODE(eacP1) || MODE(eacP2) || MODE(eacP3))
185 for (m = 0; (m < DIM); m++)
190 cth = cos_angle(xj, xk);
192 if (cth-1.0 > 1.0e-15)
194 printf("j: %d, k: %d, xj:(%g,%g,%g), xk:(%g,%g,%g)\n",
195 j, k, xj[XX], xj[YY], xj[ZZ], xk[XX], xk[YY], xk[ZZ]);
202 else if (MODE(eacP3))
206 corr[k] += LegendreP(cth, mmm); /* 1.5*cth*cth-0.5; */
208 else if (MODE(eacRcross))
211 for (m = 0; (m < DIM); m++)
218 corr[k] += iprod(rr, rr);
220 else if (MODE(eacVector))
222 for (m = 0; (m < DIM); m++)
233 gmx_fatal(FARGS, "\nInvalid mode (%lu) in do_ac_core", mode);
237 /* Correct for the number of points and copy results to the data array */
238 for (j = 0; (j < nout); j++)
240 n = (nframes-j+(nrestart-1))/nrestart;
245 /*! \brief Routine to normalize ACF, dividing by corr[0]. */
246 static void normalize_acf(int nout, real corr[])
253 fprintf(debug, "Before normalization\n");
254 for (j = 0; (j < nout); j++)
256 fprintf(debug, "%5d %10f\n", j, corr[j]);
260 /* Normalisation makes that c[0] = 1.0 and that other points are scaled
263 if (fabs(corr[0]) < 1e-5)
271 for (j = 0; (j < nout); j++)
278 fprintf(debug, "After normalization\n");
279 for (j = 0; (j < nout); j++)
281 fprintf(debug, "%5d %10f\n", j, corr[j]);
286 /*! \brief Routine that averages ACFs. */
287 static void average_acf(gmx_bool bVerbose, int n, int nitem, real **c1)
294 printf("Averaging correlation functions\n");
297 for (j = 0; (j < n); j++)
300 for (i = 0; (i < nitem); i++)
308 /*! \brief Normalize ACFs. */
309 static void norm_and_scale_vectors(int nframes, real c1[], real scale)
314 for (j = 0; (j < nframes); j++)
318 for (m = 0; (m < DIM); m++)
325 /*! \brief Debugging */
326 static void dump_tmp(char *s, int n, real c[])
331 fp = gmx_ffopen(s, "w");
332 for (i = 0; (i < n); i++)
334 fprintf(fp, "%10d %10g\n", i, c[i]);
339 /*! \brief High level ACF routine. */
340 static void do_four_core(unsigned long mode, int nframes,
341 real c1[], real csum[], real ctmp[])
348 snew(cfour, nframes);
352 /********************************************
354 ********************************************/
355 low_do_four_core(nframes, c1, csum, enNorm);
357 else if (MODE(eacCos))
359 /***************************************************
361 ***************************************************/
362 /* Copy the data to temp array. Since we need it twice
363 * we can't overwrite original.
365 for (j = 0; (j < nframes); j++)
370 /* Cosine term of AC function */
371 low_do_four_core(nframes, ctmp, cfour, enCos);
372 for (j = 0; (j < nframes); j++)
377 /* Sine term of AC function */
378 low_do_four_core(nframes, ctmp, cfour, enSin);
379 for (j = 0; (j < nframes); j++)
385 else if (MODE(eacP2))
387 /***************************************************
388 * Legendre polynomials
389 ***************************************************/
390 /* First normalize the vectors */
391 norm_and_scale_vectors(nframes, c1, 1.0);
393 /* For P2 thingies we have to do six FFT based correls
394 * First for XX^2, then for YY^2, then for ZZ^2
395 * Then we have to do XY, YZ and XZ (counting these twice)
396 * After that we sum them and normalise
397 * P2(x) = (3 * cos^2 (x) - 1)/2
398 * for unit vectors u and v we compute the cosine as the inner product
399 * cos(u,v) = uX vX + uY vY + uZ vZ
403 * C(t) = | (3 cos^2(u(t'),u(t'+t)) - 1)/2 dt'
408 * P2(u(0),u(t)) = [3 * (uX(0) uX(t) +
410 * uZ(0) uZ(t))^2 - 1]/2
411 * = [3 * ((uX(0) uX(t))^2 +
414 * 2(uX(0) uY(0) uX(t) uY(t)) +
415 * 2(uX(0) uZ(0) uX(t) uZ(t)) +
416 * 2(uY(0) uZ(0) uY(t) uZ(t))) - 1]/2
418 * = [(3/2) * (<uX^2> + <uY^2> + <uZ^2> +
419 * 2<uXuY> + 2<uXuZ> + 2<uYuZ>) - 0.5]
423 /* Because of normalization the number of -0.5 to subtract
424 * depends on the number of data points!
426 for (j = 0; (j < nframes); j++)
428 csum[j] = -0.5*(nframes-j);
431 /***** DIAGONAL ELEMENTS ************/
432 for (m = 0; (m < DIM); m++)
434 /* Copy the vector data in a linear array */
435 for (j = 0; (j < nframes); j++)
437 ctmp[j] = gmx::square(c1[DIM*j+m]);
441 sprintf(buf, "c1diag%d.xvg", m);
442 dump_tmp(buf, nframes, ctmp);
445 low_do_four_core(nframes, ctmp, cfour, enNorm);
449 sprintf(buf, "c1dfout%d.xvg", m);
450 dump_tmp(buf, nframes, cfour);
453 for (j = 0; (j < nframes); j++)
455 csum[j] += fac*(cfour[j]);
458 /******* OFF-DIAGONAL ELEMENTS **********/
459 for (m = 0; (m < DIM); m++)
461 /* Copy the vector data in a linear array */
463 for (j = 0; (j < nframes); j++)
465 ctmp[j] = c1[DIM*j+m]*c1[DIM*j+m1];
470 sprintf(buf, "c1off%d.xvg", m);
471 dump_tmp(buf, nframes, ctmp);
473 low_do_four_core(nframes, ctmp, cfour, enNorm);
476 sprintf(buf, "c1ofout%d.xvg", m);
477 dump_tmp(buf, nframes, cfour);
480 for (j = 0; (j < nframes); j++)
482 csum[j] += fac*cfour[j];
486 else if (MODE(eacP1) || MODE(eacVector))
488 /***************************************************
490 ***************************************************/
493 /* First normalize the vectors */
494 norm_and_scale_vectors(nframes, c1, 1.0);
497 /* For vector thingies we have to do three FFT based correls
498 * First for XX, then for YY, then for ZZ
499 * After that we sum them and normalise
501 for (j = 0; (j < nframes); j++)
505 for (m = 0; (m < DIM); m++)
507 /* Copy the vector data in a linear array */
508 for (j = 0; (j < nframes); j++)
510 ctmp[j] = c1[DIM*j+m];
512 low_do_four_core(nframes, ctmp, cfour, enNorm);
513 for (j = 0; (j < nframes); j++)
521 gmx_fatal(FARGS, "\nUnknown mode in do_autocorr (%lu)", mode);
525 for (j = 0; (j < nframes); j++)
527 c1[j] = csum[j]/static_cast<real>(nframes-j);
531 void low_do_autocorr(const char *fn, const gmx_output_env_t *oenv, const char *title,
532 int nframes, int nitem, int nout, real **c1,
533 real dt, unsigned long mode, int nrestart,
534 gmx_bool bAver, gmx_bool bNormalize,
535 gmx_bool bVerbose, real tbeginfit, real tendfit,
538 FILE *fp, *gp = nullptr;
542 real sum, Ct2av, Ctav;
543 gmx_bool bFour = acf.bFour;
545 /* Check flags and parameters */
546 nout = get_acfnout();
549 nout = acf.nout = (nframes+1)/2;
551 else if (nout > nframes)
556 if (MODE(eacCos) && MODE(eacVector))
558 gmx_fatal(FARGS, "Incompatible options bCos && bVector (%s, %d)",
561 if ((MODE(eacP3) || MODE(eacRcross)) && bFour)
565 fprintf(stderr, "Can't combine mode %lu with FFT, turning off FFT\n", mode);
569 if (MODE(eacNormal) && MODE(eacVector))
571 gmx_fatal(FARGS, "Incompatible mode bits: normal and vector (or Legendre)");
574 /* Print flags and parameters */
577 printf("Will calculate %s of %d thingies for %d frames\n",
578 title ? title : "autocorrelation", nitem, nframes);
579 printf("bAver = %s, bFour = %s bNormalize= %s\n",
580 gmx::boolToString(bAver), gmx::boolToString(bFour),
581 gmx::boolToString(bNormalize));
582 printf("mode = %lu, dt = %g, nrestart = %d\n", mode, dt, nrestart);
584 /* Allocate temp arrays */
588 /* Loop over items (e.g. molecules or dihedrals)
589 * In this loop the actual correlation functions are computed, but without
592 for (int i = 0; i < nitem; i++)
594 if (bVerbose && (((i % 100) == 0) || (i == nitem-1)))
596 fprintf(stderr, "\rThingie %d", i+1);
602 do_four_core(mode, nframes, c1[i], csum, ctmp);
606 do_ac_core(nframes, nout, ctmp, c1[i], nrestart, mode);
611 fprintf(stderr, "\n");
619 fp = xvgropen(fn, title, "Time (ps)", "C(t)", oenv);
630 average_acf(bVerbose, nframes, nitem, c1);
635 normalize_acf(nout, c1[0]);
638 if (eFitFn != effnNONE)
640 fit_acf(nout, eFitFn, oenv, fn != nullptr, tbeginfit, tendfit, dt, c1[0], fit);
641 sum = print_and_integrate(fp, nout, dt, c1[0], fit, 1);
645 sum = print_and_integrate(fp, nout, dt, c1[0], nullptr, 1);
649 printf("Correlation time (integral over corrfn): %g (ps)\n", sum);
654 /* Not averaging. Normalize individual ACFs */
658 gp = xvgropen("ct-distr.xvg", "Correlation times", "item", "time (ps)", oenv);
660 for (i = 0; i < nitem; i++)
664 normalize_acf(nout, c1[i]);
666 if (eFitFn != effnNONE)
668 fit_acf(nout, eFitFn, oenv, fn != nullptr, tbeginfit, tendfit, dt, c1[i], fit);
669 sum = print_and_integrate(fp, nout, dt, c1[i], fit, 1);
673 sum = print_and_integrate(fp, nout, dt, c1[i], nullptr, 1);
677 "CORRelation time (integral over corrfn %d): %g (ps)\n",
685 fprintf(gp, "%5d %.3f\n", i, sum);
696 printf("Average correlation time %.3f Std. Dev. %.3f Error %.3f (ps)\n",
697 Ctav, std::sqrt((Ct2av - gmx::square(Ctav))),
698 std::sqrt((Ct2av - gmx::square(Ctav))/(nitem-1)));
708 /*! \brief Legend for selecting Legendre polynomials. */
709 static const char *Leg[] = { nullptr, "0", "1", "2", "3", nullptr };
711 t_pargs *add_acf_pargs(int *npargs, t_pargs *pa)
714 { "-acflen", FALSE, etINT, {&acf.nout},
715 "Length of the ACF, default is half the number of frames" },
716 { "-normalize", FALSE, etBOOL, {&acf.bNormalize},
718 { "-fftcorr", FALSE, etBOOL, {&acf.bFour},
719 "HIDDENUse fast fourier transform for correlation function" },
720 { "-nrestart", FALSE, etINT, {&acf.nrestart},
721 "HIDDENNumber of frames between time origins for ACF when no FFT is used" },
722 { "-P", FALSE, etENUM, {Leg},
723 "Order of Legendre polynomial for ACF (0 indicates none)" },
724 { "-fitfn", FALSE, etENUM, {s_ffn},
726 { "-beginfit", FALSE, etREAL, {&acf.tbeginfit},
727 "Time where to begin the exponential fit of the correlation function" },
728 { "-endfit", FALSE, etREAL, {&acf.tendfit},
729 "Time where to end the exponential fit of the correlation function, -1 is until the end" },
735 snew(ppa, *npargs+npa);
736 for (i = 0; (i < *npargs); i++)
740 for (i = 0; (i < npa); i++)
742 ppa[*npargs+i] = acfpa[i];
750 acf.fitfn = effnEXP1;
752 acf.bNormalize = TRUE;
761 void do_autocorr(const char *fn, const gmx_output_env_t *oenv, const char *title,
762 int nframes, int nitem, real **c1,
763 real dt, unsigned long mode, gmx_bool bAver)
767 printf("ACF data structures have not been initialised. Call add_acf_pargs\n");
770 /* Handle enumerated types */
771 sscanf(Leg[0], "%d", &acf.P);
772 acf.fitfn = sffn2effn(s_ffn);
789 low_do_autocorr(fn, oenv, title, nframes, nitem, acf.nout, c1, dt, mode,
790 acf.nrestart, bAver, acf.bNormalize,
791 bDebugMode(), acf.tbeginfit, acf.tendfit,
799 gmx_fatal(FARGS, "ACF data not initialized yet");
809 gmx_fatal(FARGS, "ACF data not initialized yet");
812 return sffn2effn(s_ffn);