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 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
40 * Implements function to compute many autocorrelation functions
42 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
43 * \ingroup module_correlationfunctions
55 #include "gromacs/correlationfunctions/expfit.h"
56 #include "gromacs/correlationfunctions/integrate.h"
57 #include "gromacs/correlationfunctions/manyautocorrelation.h"
58 #include "gromacs/correlationfunctions/polynomials.h"
59 #include "gromacs/fileio/xvgr.h"
60 #include "gromacs/math/functions.h"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/utility/arraysize.h"
63 #include "gromacs/utility/fatalerror.h"
64 #include "gromacs/utility/futil.h"
65 #include "gromacs/utility/real.h"
66 #include "gromacs/utility/smalloc.h"
67 #include "gromacs/utility/strconvert.h"
69 /*! \brief Shortcut macro to select modes. */
70 #define MODE(x) ((mode & (x)) == (x))
75 int nrestart, nout, P, fitfn;
76 gmx_bool bFour, bNormalize;
77 real tbeginfit, tendfit;
80 /*! \brief Global variable set true if initialization routines are called. */
81 // NOLINTNEXTLINE(cppcoreguidelines-avoid-non-const-global-variables)
82 static gmx_bool bACFinit = FALSE;
84 /*! \brief Data structure for storing command line variables. */
85 // NOLINTNEXTLINE(cppcoreguidelines-avoid-non-const-global-variables)
95 /*! \brief Routine to compute ACF using FFT. */
96 static void low_do_four_core(int nframes, real c1[], real cfour[], int nCos)
99 std::vector<std::vector<real>> data;
101 data[0].resize(nframes, 0);
105 for (i = 0; (i < nframes); i++)
111 for (i = 0; (i < nframes); i++)
113 data[0][i] = cos(c1[i]);
117 for (i = 0; (i < nframes); i++)
119 data[0][i] = sin(c1[i]);
122 default: gmx_fatal(FARGS, "nCos = %d, %s %d", nCos, __FILE__, __LINE__);
125 many_auto_correl(&data);
126 for (i = 0; (i < nframes); i++)
128 cfour[i] = data[0][i];
132 /*! \brief Routine to comput ACF without FFT. */
133 static void do_ac_core(int nframes, int nout, real corr[], real c1[], int nrestart, unsigned long mode)
135 int j, k, j3, jk3, m, n;
141 printf("WARNING: setting number of restarts to 1\n");
146 fprintf(debug, "Starting do_ac_core: nframes=%d, nout=%d, nrestart=%d,mode=%lu\n", nframes, nout, nrestart, mode);
149 for (j = 0; (j < nout); j++)
154 /* Loop over starting points. */
155 for (j = 0; (j < nframes); j += nrestart)
159 /* Loop over the correlation length for this starting point */
160 for (k = 0; (k < nout) && (j + k < nframes); k++)
164 /* Switch over possible ACF types.
165 * It might be more efficient to put the loops inside the switch,
166 * but this is more clear, and save development time!
170 corr[k] += c1[j] * c1[j + k];
172 else if (MODE(eacCos))
174 /* Compute the cos (phi(t)-phi(t+dt)) */
175 corr[k] += std::cos(c1[j] - c1[j + k]);
177 else if (MODE(eacIden))
179 /* Check equality (phi(t)==phi(t+dt)) */
180 corr[k] += (c1[j] == c1[j + k]) ? 1 : 0;
182 else if (MODE(eacP1) || MODE(eacP2) || MODE(eacP3))
186 for (m = 0; (m < DIM); m++)
191 cth = cos_angle(xj, xk);
193 if (cth - 1.0 > 1.0e-15)
195 printf("j: %d, k: %d, xj:(%g,%g,%g), xk:(%g,%g,%g)\n",
210 else if (MODE(eacP3))
214 corr[k] += LegendreP(cth, mmm); /* 1.5*cth*cth-0.5; */
216 else if (MODE(eacRcross))
219 for (m = 0; (m < DIM); m++)
226 corr[k] += iprod(rr, rr);
228 else if (MODE(eacVector))
230 for (m = 0; (m < DIM); m++)
241 gmx_fatal(FARGS, "\nInvalid mode (%lu) in do_ac_core", mode);
245 /* Correct for the number of points and copy results to the data array */
246 for (j = 0; (j < nout); j++)
248 n = (nframes - j + (nrestart - 1)) / nrestart;
253 /*! \brief Routine to normalize ACF, dividing by corr[0]. */
254 static void normalize_acf(int nout, real corr[])
261 fprintf(debug, "Before normalization\n");
262 for (j = 0; (j < nout); j++)
264 fprintf(debug, "%5d %10f\n", j, corr[j]);
268 /* Normalisation makes that c[0] = 1.0 and that other points are scaled
271 if (fabs(corr[0]) < 1e-5)
279 for (j = 0; (j < nout); j++)
286 fprintf(debug, "After normalization\n");
287 for (j = 0; (j < nout); j++)
289 fprintf(debug, "%5d %10f\n", j, corr[j]);
294 /*! \brief Routine that averages ACFs. */
295 static void average_acf(gmx_bool bVerbose, int n, int nitem, real** c1)
302 printf("Averaging correlation functions\n");
305 for (j = 0; (j < n); j++)
308 for (i = 0; (i < nitem); i++)
312 c1[0][j] = c0 / nitem;
316 /*! \brief Normalize ACFs. */
317 static void norm_and_scale_vectors(int nframes, real c1[], real scale)
322 for (j = 0; (j < nframes); j++)
324 rij = &(c1[j * DIM]);
326 for (m = 0; (m < DIM); m++)
333 /*! \brief Debugging */
334 static void dump_tmp(char* s, int n, real c[])
339 fp = gmx_ffopen(s, "w");
340 for (i = 0; (i < n); i++)
342 fprintf(fp, "%10d %10g\n", i, c[i]);
347 /*! \brief High level ACF routine. */
348 static void do_four_core(unsigned long mode, int nframes, real c1[], real csum[], real ctmp[])
355 snew(cfour, nframes);
359 /********************************************
361 ********************************************/
362 low_do_four_core(nframes, c1, csum, enNorm);
364 else if (MODE(eacCos))
366 /***************************************************
368 ***************************************************/
369 /* Copy the data to temp array. Since we need it twice
370 * we can't overwrite original.
372 for (j = 0; (j < nframes); j++)
377 /* Cosine term of AC function */
378 low_do_four_core(nframes, ctmp, cfour, enCos);
379 for (j = 0; (j < nframes); j++)
384 /* Sine term of AC function */
385 low_do_four_core(nframes, ctmp, cfour, enSin);
386 for (j = 0; (j < nframes); j++)
392 else if (MODE(eacP2))
394 /***************************************************
395 * Legendre polynomials
396 ***************************************************/
397 /* First normalize the vectors */
398 norm_and_scale_vectors(nframes, c1, 1.0);
400 /* For P2 thingies we have to do six FFT based correls
401 * First for XX^2, then for YY^2, then for ZZ^2
402 * Then we have to do XY, YZ and XZ (counting these twice)
403 * After that we sum them and normalise
404 * P2(x) = (3 * cos^2 (x) - 1)/2
405 * for unit vectors u and v we compute the cosine as the inner product
406 * cos(u,v) = uX vX + uY vY + uZ vZ
410 * C(t) = | (3 cos^2(u(t'),u(t'+t)) - 1)/2 dt'
415 * P2(u(0),u(t)) = [3 * (uX(0) uX(t) +
417 * uZ(0) uZ(t))^2 - 1]/2
418 * = [3 * ((uX(0) uX(t))^2 +
421 * 2(uX(0) uY(0) uX(t) uY(t)) +
422 * 2(uX(0) uZ(0) uX(t) uZ(t)) +
423 * 2(uY(0) uZ(0) uY(t) uZ(t))) - 1]/2
425 * = [(3/2) * (<uX^2> + <uY^2> + <uZ^2> +
426 * 2<uXuY> + 2<uXuZ> + 2<uYuZ>) - 0.5]
430 /* Because of normalization the number of -0.5 to subtract
431 * depends on the number of data points!
433 for (j = 0; (j < nframes); j++)
435 csum[j] = -0.5 * (nframes - j);
438 /***** DIAGONAL ELEMENTS ************/
439 for (m = 0; (m < DIM); m++)
441 /* Copy the vector data in a linear array */
442 for (j = 0; (j < nframes); j++)
444 ctmp[j] = gmx::square(c1[DIM * j + m]);
448 sprintf(buf, "c1diag%d.xvg", m);
449 dump_tmp(buf, nframes, ctmp);
452 low_do_four_core(nframes, ctmp, cfour, enNorm);
456 sprintf(buf, "c1dfout%d.xvg", m);
457 dump_tmp(buf, nframes, cfour);
460 for (j = 0; (j < nframes); j++)
462 csum[j] += fac * (cfour[j]);
465 /******* OFF-DIAGONAL ELEMENTS **********/
466 for (m = 0; (m < DIM); m++)
468 /* Copy the vector data in a linear array */
470 for (j = 0; (j < nframes); j++)
472 ctmp[j] = c1[DIM * j + m] * c1[DIM * j + m1];
477 sprintf(buf, "c1off%d.xvg", m);
478 dump_tmp(buf, nframes, ctmp);
480 low_do_four_core(nframes, ctmp, cfour, enNorm);
483 sprintf(buf, "c1ofout%d.xvg", m);
484 dump_tmp(buf, nframes, cfour);
487 for (j = 0; (j < nframes); j++)
489 csum[j] += fac * cfour[j];
493 else if (MODE(eacP1) || MODE(eacVector))
495 /***************************************************
497 ***************************************************/
500 /* First normalize the vectors */
501 norm_and_scale_vectors(nframes, c1, 1.0);
504 /* For vector thingies we have to do three FFT based correls
505 * First for XX, then for YY, then for ZZ
506 * After that we sum them and normalise
508 for (j = 0; (j < nframes); j++)
512 for (m = 0; (m < DIM); m++)
514 /* Copy the vector data in a linear array */
515 for (j = 0; (j < nframes); j++)
517 ctmp[j] = c1[DIM * j + m];
519 low_do_four_core(nframes, ctmp, cfour, enNorm);
520 for (j = 0; (j < nframes); j++)
528 gmx_fatal(FARGS, "\nUnknown mode in do_autocorr (%lu)", mode);
532 for (j = 0; (j < nframes); j++)
534 c1[j] = csum[j] / static_cast<real>(nframes - j);
538 void low_do_autocorr(const char* fn,
539 const gmx_output_env_t* oenv,
555 FILE * fp, *gp = nullptr;
559 real sum, Ct2av, Ctav;
560 gmx_bool bFour = acf.bFour;
562 /* Check flags and parameters */
563 nout = get_acfnout();
566 nout = acf.nout = (nframes + 1) / 2;
568 else if (nout > nframes)
573 if (MODE(eacCos) && MODE(eacVector))
575 gmx_fatal(FARGS, "Incompatible options bCos && bVector (%s, %d)", __FILE__, __LINE__);
577 if ((MODE(eacP3) || MODE(eacRcross)) && bFour)
581 fprintf(stderr, "Can't combine mode %lu with FFT, turning off FFT\n", mode);
585 if (MODE(eacNormal) && MODE(eacVector))
587 gmx_fatal(FARGS, "Incompatible mode bits: normal and vector (or Legendre)");
590 /* Print flags and parameters */
593 printf("Will calculate %s of %d thingies for %d frames\n", title ? title : "autocorrelation", nitem, nframes);
594 printf("bAver = %s, bFour = %s bNormalize= %s\n",
595 gmx::boolToString(bAver),
596 gmx::boolToString(bFour),
597 gmx::boolToString(bNormalize));
598 printf("mode = %lu, dt = %g, nrestart = %d\n", mode, dt, nrestart);
600 /* Allocate temp arrays */
604 /* Loop over items (e.g. molecules or dihedrals)
605 * In this loop the actual correlation functions are computed, but without
608 for (int i = 0; i < nitem; i++)
610 if (bVerbose && (((i % 100) == 0) || (i == nitem - 1)))
612 fprintf(stderr, "\rThingie %d", i + 1);
618 do_four_core(mode, nframes, c1[i], csum, ctmp);
622 do_ac_core(nframes, nout, ctmp, c1[i], nrestart, mode);
627 fprintf(stderr, "\n");
635 fp = xvgropen(fn, title, "Time (ps)", "C(t)", oenv);
646 average_acf(bVerbose, nframes, nitem, c1);
651 normalize_acf(nout, c1[0]);
654 if (eFitFn != effnNONE)
656 fit_acf(nout, eFitFn, oenv, fn != nullptr, tbeginfit, tendfit, dt, c1[0], fit);
657 sum = print_and_integrate(fp, nout, dt, c1[0], fit, 1);
661 sum = print_and_integrate(fp, nout, dt, c1[0], nullptr, 1);
665 printf("Correlation time (integral over corrfn): %g (ps)\n", sum);
670 /* Not averaging. Normalize individual ACFs */
674 gp = xvgropen("ct-distr.xvg", "Correlation times", "item", "time (ps)", oenv);
676 for (i = 0; i < nitem; i++)
680 normalize_acf(nout, c1[i]);
682 if (eFitFn != effnNONE)
684 fit_acf(nout, eFitFn, oenv, fn != nullptr, tbeginfit, tendfit, dt, c1[i], fit);
685 sum = print_and_integrate(fp, nout, dt, c1[i], fit, 1);
689 sum = print_and_integrate(fp, nout, dt, c1[i], nullptr, 1);
692 fprintf(debug, "CORRelation time (integral over corrfn %d): %g (ps)\n", i, sum);
699 fprintf(gp, "%5d %.3f\n", i, sum);
710 printf("Average correlation time %.3f Std. Dev. %.3f Error %.3f (ps)\n",
712 std::sqrt((Ct2av - gmx::square(Ctav))),
713 std::sqrt((Ct2av - gmx::square(Ctav)) / (nitem - 1)));
723 /*! \brief Legend for selecting Legendre polynomials. */
724 // NOLINTNEXTLINE(cppcoreguidelines-avoid-non-const-global-variables)
725 static const char* Leg[] = { nullptr, "0", "1", "2", "3", nullptr };
727 t_pargs* add_acf_pargs(int* npargs, t_pargs* pa)
734 "Length of the ACF, default is half the number of frames" },
735 { "-normalize", FALSE, etBOOL, { &acf.bNormalize }, "Normalize ACF" },
740 "HIDDENUse fast fourier transform for correlation function" },
745 "HIDDENNumber of frames between time origins for ACF when no FFT is used" },
746 { "-P", FALSE, etENUM, { Leg }, "Order of Legendre polynomial for ACF (0 indicates none)" },
747 { "-fitfn", FALSE, etENUM, { s_ffn }, "Fit function" },
752 "Time where to begin the exponential fit of the correlation function" },
757 "Time where to end the exponential fit of the correlation function, -1 is until the "
764 snew(ppa, *npargs + npa);
765 for (i = 0; (i < *npargs); i++)
769 for (i = 0; (i < npa); i++)
771 ppa[*npargs + i] = acfpa[i];
779 acf.fitfn = effnEXP1;
781 acf.bNormalize = TRUE;
790 void do_autocorr(const char* fn,
791 const gmx_output_env_t* oenv,
802 printf("ACF data structures have not been initialised. Call add_acf_pargs\n");
805 /* Handle enumerated types */
806 sscanf(Leg[0], "%d", &acf.P);
807 acf.fitfn = sffn2effn(s_ffn);
811 case 1: mode = mode | eacP1; break;
812 case 2: mode = mode | eacP2; break;
813 case 3: mode = mode | eacP3; break;
839 gmx_fatal(FARGS, "ACF data not initialized yet");
849 gmx_fatal(FARGS, "ACF data not initialized yet");
852 return sffn2effn(s_ffn);