Use proper doxygen tags in modular simulator
[gromacs.git] / src / gromacs / correlationfunctions / autocorr.cpp
blob3c6b4b7018dd1c796695ba993535452e19da5115
1 /*
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, 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.
38 /*! \internal \file
39 * \brief
40 * Implements function to compute many autocorrelation functions
42 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
43 * \ingroup module_correlationfunctions
45 #include "gmxpre.h"
47 #include "autocorr.h"
49 #include <cmath>
50 #include <cstdio>
51 #include <cstring>
53 #include <algorithm>
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))
72 typedef struct
74 unsigned long mode;
75 int nrestart, nout, P, fitfn;
76 gmx_bool bFour, bNormalize;
77 real tbeginfit, tendfit;
78 } t_acf;
80 /*! \brief Global variable set true if initialization routines are called. */
81 static gmx_bool bACFinit = FALSE;
83 /*! \brief Data structure for storing command line variables. */
84 static t_acf acf;
86 enum
88 enNorm,
89 enCos,
90 enSin
93 /*! \brief Routine to compute ACF using FFT. */
94 static void low_do_four_core(int nframes, real c1[], real cfour[], int nCos)
96 int i = 0;
97 std::vector<std::vector<real>> data;
98 data.resize(1);
99 data[0].resize(nframes, 0);
100 switch (nCos)
102 case enNorm:
103 for (i = 0; (i < nframes); i++)
105 data[0][i] = c1[i];
107 break;
108 case enCos:
109 for (i = 0; (i < nframes); i++)
111 data[0][i] = cos(c1[i]);
113 break;
114 case enSin:
115 for (i = 0; (i < nframes); i++)
117 data[0][i] = sin(c1[i]);
119 break;
120 default: gmx_fatal(FARGS, "nCos = %d, %s %d", nCos, __FILE__, __LINE__);
123 many_auto_correl(&data);
124 for (i = 0; (i < nframes); i++)
126 cfour[i] = data[0][i];
130 /*! \brief Routine to comput ACF without FFT. */
131 static void do_ac_core(int nframes, int nout, real corr[], real c1[], int nrestart, unsigned long mode)
133 int j, k, j3, jk3, m, n;
134 real ccc, cth;
135 rvec xj, xk;
137 if (nrestart < 1)
139 printf("WARNING: setting number of restarts to 1\n");
140 nrestart = 1;
142 if (debug)
144 fprintf(debug, "Starting do_ac_core: nframes=%d, nout=%d, nrestart=%d,mode=%lu\n", nframes,
145 nout, nrestart, mode);
148 for (j = 0; (j < nout); j++)
150 corr[j] = 0;
153 /* Loop over starting points. */
154 for (j = 0; (j < nframes); j += nrestart)
156 j3 = DIM * j;
158 /* Loop over the correlation length for this starting point */
159 for (k = 0; (k < nout) && (j + k < nframes); k++)
161 jk3 = DIM * (j + 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!
167 if (MODE(eacNormal))
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))
183 unsigned int mmm;
185 for (m = 0; (m < DIM); m++)
187 xj[m] = c1[j3 + m];
188 xk[m] = c1[jk3 + 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", j, k, xj[XX], xj[YY],
195 xj[ZZ], xk[XX], xk[YY], xk[ZZ]);
197 mmm = 1;
198 if (MODE(eacP2))
200 mmm = 2;
202 else if (MODE(eacP3))
204 mmm = 3;
206 corr[k] += LegendreP(cth, mmm); /* 1.5*cth*cth-0.5; */
208 else if (MODE(eacRcross))
210 rvec xj, xk, rr;
211 for (m = 0; (m < DIM); m++)
213 xj[m] = c1[j3 + m];
214 xk[m] = c1[jk3 + m];
216 cprod(xj, xk, rr);
218 corr[k] += iprod(rr, rr);
220 else if (MODE(eacVector))
222 for (m = 0; (m < DIM); m++)
224 xj[m] = c1[j3 + m];
225 xk[m] = c1[jk3 + m];
227 ccc = iprod(xj, xk);
229 corr[k] += ccc;
231 else
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;
241 c1[j] = corr[j] / n;
245 /*! \brief Routine to normalize ACF, dividing by corr[0]. */
246 static void normalize_acf(int nout, real corr[])
248 int j;
249 double c0;
251 if (debug)
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
261 * accordingly.
263 if (fabs(corr[0]) < 1e-5)
265 c0 = 1.0;
267 else
269 c0 = 1.0 / corr[0];
271 for (j = 0; (j < nout); j++)
273 corr[j] *= c0;
276 if (debug)
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)
289 real c0;
290 int i, j;
292 if (bVerbose)
294 printf("Averaging correlation functions\n");
297 for (j = 0; (j < n); j++)
299 c0 = 0;
300 for (i = 0; (i < nitem); i++)
302 c0 += c1[i][j];
304 c1[0][j] = c0 / nitem;
308 /*! \brief Normalize ACFs. */
309 static void norm_and_scale_vectors(int nframes, real c1[], real scale)
311 int j, m;
312 real* rij;
314 for (j = 0; (j < nframes); j++)
316 rij = &(c1[j * DIM]);
317 unitv(rij, rij);
318 for (m = 0; (m < DIM); m++)
320 rij[m] *= scale;
325 /*! \brief Debugging */
326 static void dump_tmp(char* s, int n, real c[])
328 FILE* fp;
329 int i;
331 fp = gmx_ffopen(s, "w");
332 for (i = 0; (i < n); i++)
334 fprintf(fp, "%10d %10g\n", i, c[i]);
336 gmx_ffclose(fp);
339 /*! \brief High level ACF routine. */
340 static void do_four_core(unsigned long mode, int nframes, real c1[], real csum[], real ctmp[])
342 real* cfour;
343 char buf[32];
344 real fac;
345 int j, m, m1;
347 snew(cfour, nframes);
349 if (MODE(eacNormal))
351 /********************************************
352 * N O R M A L
353 ********************************************/
354 low_do_four_core(nframes, c1, csum, enNorm);
356 else if (MODE(eacCos))
358 /***************************************************
359 * C O S I N E
360 ***************************************************/
361 /* Copy the data to temp array. Since we need it twice
362 * we can't overwrite original.
364 for (j = 0; (j < nframes); j++)
366 ctmp[j] = c1[j];
369 /* Cosine term of AC function */
370 low_do_four_core(nframes, ctmp, cfour, enCos);
371 for (j = 0; (j < nframes); j++)
373 c1[j] = cfour[j];
376 /* Sine term of AC function */
377 low_do_four_core(nframes, ctmp, cfour, enSin);
378 for (j = 0; (j < nframes); j++)
380 c1[j] += cfour[j];
381 csum[j] = c1[j];
384 else if (MODE(eacP2))
386 /***************************************************
387 * Legendre polynomials
388 ***************************************************/
389 /* First normalize the vectors */
390 norm_and_scale_vectors(nframes, c1, 1.0);
392 /* For P2 thingies we have to do six FFT based correls
393 * First for XX^2, then for YY^2, then for ZZ^2
394 * Then we have to do XY, YZ and XZ (counting these twice)
395 * After that we sum them and normalise
396 * P2(x) = (3 * cos^2 (x) - 1)/2
397 * for unit vectors u and v we compute the cosine as the inner product
398 * cos(u,v) = uX vX + uY vY + uZ vZ
400 * oo
402 * C(t) = | (3 cos^2(u(t'),u(t'+t)) - 1)/2 dt'
406 * For ACF we need:
407 * P2(u(0),u(t)) = [3 * (uX(0) uX(t) +
408 * uY(0) uY(t) +
409 * uZ(0) uZ(t))^2 - 1]/2
410 * = [3 * ((uX(0) uX(t))^2 +
411 * (uY(0) uY(t))^2 +
412 * (uZ(0) uZ(t))^2 +
413 * 2(uX(0) uY(0) uX(t) uY(t)) +
414 * 2(uX(0) uZ(0) uX(t) uZ(t)) +
415 * 2(uY(0) uZ(0) uY(t) uZ(t))) - 1]/2
417 * = [(3/2) * (<uX^2> + <uY^2> + <uZ^2> +
418 * 2<uXuY> + 2<uXuZ> + 2<uYuZ>) - 0.5]
422 /* Because of normalization the number of -0.5 to subtract
423 * depends on the number of data points!
425 for (j = 0; (j < nframes); j++)
427 csum[j] = -0.5 * (nframes - j);
430 /***** DIAGONAL ELEMENTS ************/
431 for (m = 0; (m < DIM); m++)
433 /* Copy the vector data in a linear array */
434 for (j = 0; (j < nframes); j++)
436 ctmp[j] = gmx::square(c1[DIM * j + m]);
438 if (debug)
440 sprintf(buf, "c1diag%d.xvg", m);
441 dump_tmp(buf, nframes, ctmp);
444 low_do_four_core(nframes, ctmp, cfour, enNorm);
446 if (debug)
448 sprintf(buf, "c1dfout%d.xvg", m);
449 dump_tmp(buf, nframes, cfour);
451 fac = 1.5;
452 for (j = 0; (j < nframes); j++)
454 csum[j] += fac * (cfour[j]);
457 /******* OFF-DIAGONAL ELEMENTS **********/
458 for (m = 0; (m < DIM); m++)
460 /* Copy the vector data in a linear array */
461 m1 = (m + 1) % DIM;
462 for (j = 0; (j < nframes); j++)
464 ctmp[j] = c1[DIM * j + m] * c1[DIM * j + m1];
467 if (debug)
469 sprintf(buf, "c1off%d.xvg", m);
470 dump_tmp(buf, nframes, ctmp);
472 low_do_four_core(nframes, ctmp, cfour, enNorm);
473 if (debug)
475 sprintf(buf, "c1ofout%d.xvg", m);
476 dump_tmp(buf, nframes, cfour);
478 fac = 3.0;
479 for (j = 0; (j < nframes); j++)
481 csum[j] += fac * cfour[j];
485 else if (MODE(eacP1) || MODE(eacVector))
487 /***************************************************
488 * V E C T O R & P1
489 ***************************************************/
490 if (MODE(eacP1))
492 /* First normalize the vectors */
493 norm_and_scale_vectors(nframes, c1, 1.0);
496 /* For vector thingies we have to do three FFT based correls
497 * First for XX, then for YY, then for ZZ
498 * After that we sum them and normalise
500 for (j = 0; (j < nframes); j++)
502 csum[j] = 0.0;
504 for (m = 0; (m < DIM); m++)
506 /* Copy the vector data in a linear array */
507 for (j = 0; (j < nframes); j++)
509 ctmp[j] = c1[DIM * j + m];
511 low_do_four_core(nframes, ctmp, cfour, enNorm);
512 for (j = 0; (j < nframes); j++)
514 csum[j] += cfour[j];
518 else
520 gmx_fatal(FARGS, "\nUnknown mode in do_autocorr (%lu)", mode);
523 sfree(cfour);
524 for (j = 0; (j < nframes); j++)
526 c1[j] = csum[j] / static_cast<real>(nframes - j);
530 void low_do_autocorr(const char* fn,
531 const gmx_output_env_t* oenv,
532 const char* title,
533 int nframes,
534 int nitem,
535 int nout,
536 real** c1,
537 real dt,
538 unsigned long mode,
539 int nrestart,
540 gmx_bool bAver,
541 gmx_bool bNormalize,
542 gmx_bool bVerbose,
543 real tbeginfit,
544 real tendfit,
545 int eFitFn)
547 FILE * fp, *gp = nullptr;
548 int i;
549 real* csum;
550 real * ctmp, *fit;
551 real sum, Ct2av, Ctav;
552 gmx_bool bFour = acf.bFour;
554 /* Check flags and parameters */
555 nout = get_acfnout();
556 if (nout == -1)
558 nout = acf.nout = (nframes + 1) / 2;
560 else if (nout > nframes)
562 nout = nframes;
565 if (MODE(eacCos) && MODE(eacVector))
567 gmx_fatal(FARGS, "Incompatible options bCos && bVector (%s, %d)", __FILE__, __LINE__);
569 if ((MODE(eacP3) || MODE(eacRcross)) && bFour)
571 if (bVerbose)
573 fprintf(stderr, "Can't combine mode %lu with FFT, turning off FFT\n", mode);
575 bFour = FALSE;
577 if (MODE(eacNormal) && MODE(eacVector))
579 gmx_fatal(FARGS, "Incompatible mode bits: normal and vector (or Legendre)");
582 /* Print flags and parameters */
583 if (bVerbose)
585 printf("Will calculate %s of %d thingies for %d frames\n",
586 title ? title : "autocorrelation", nitem, nframes);
587 printf("bAver = %s, bFour = %s bNormalize= %s\n", gmx::boolToString(bAver),
588 gmx::boolToString(bFour), gmx::boolToString(bNormalize));
589 printf("mode = %lu, dt = %g, nrestart = %d\n", mode, dt, nrestart);
591 /* Allocate temp arrays */
592 snew(csum, nframes);
593 snew(ctmp, nframes);
595 /* Loop over items (e.g. molecules or dihedrals)
596 * In this loop the actual correlation functions are computed, but without
597 * normalizing them.
599 for (int i = 0; i < nitem; i++)
601 if (bVerbose && (((i % 100) == 0) || (i == nitem - 1)))
603 fprintf(stderr, "\rThingie %d", i + 1);
604 fflush(stderr);
607 if (bFour)
609 do_four_core(mode, nframes, c1[i], csum, ctmp);
611 else
613 do_ac_core(nframes, nout, ctmp, c1[i], nrestart, mode);
616 if (bVerbose)
618 fprintf(stderr, "\n");
620 sfree(ctmp);
621 sfree(csum);
623 if (fn)
625 snew(fit, nout);
626 fp = xvgropen(fn, title, "Time (ps)", "C(t)", oenv);
628 else
630 fit = nullptr;
631 fp = nullptr;
633 if (bAver)
635 if (nitem > 1)
637 average_acf(bVerbose, nframes, nitem, c1);
640 if (bNormalize)
642 normalize_acf(nout, c1[0]);
645 if (eFitFn != effnNONE)
647 fit_acf(nout, eFitFn, oenv, fn != nullptr, tbeginfit, tendfit, dt, c1[0], fit);
648 sum = print_and_integrate(fp, nout, dt, c1[0], fit, 1);
650 else
652 sum = print_and_integrate(fp, nout, dt, c1[0], nullptr, 1);
654 if (bVerbose)
656 printf("Correlation time (integral over corrfn): %g (ps)\n", sum);
659 else
661 /* Not averaging. Normalize individual ACFs */
662 Ctav = Ct2av = 0;
663 if (debug)
665 gp = xvgropen("ct-distr.xvg", "Correlation times", "item", "time (ps)", oenv);
667 for (i = 0; i < nitem; i++)
669 if (bNormalize)
671 normalize_acf(nout, c1[i]);
673 if (eFitFn != effnNONE)
675 fit_acf(nout, eFitFn, oenv, fn != nullptr, tbeginfit, tendfit, dt, c1[i], fit);
676 sum = print_and_integrate(fp, nout, dt, c1[i], fit, 1);
678 else
680 sum = print_and_integrate(fp, nout, dt, c1[i], nullptr, 1);
681 if (debug)
683 fprintf(debug, "CORRelation time (integral over corrfn %d): %g (ps)\n", i, sum);
686 Ctav += sum;
687 Ct2av += sum * sum;
688 if (debug)
690 fprintf(gp, "%5d %.3f\n", i, sum);
693 if (debug)
695 xvgrclose(gp);
697 if (nitem > 1)
699 Ctav /= nitem;
700 Ct2av /= nitem;
701 printf("Average correlation time %.3f Std. Dev. %.3f Error %.3f (ps)\n", Ctav,
702 std::sqrt((Ct2av - gmx::square(Ctav))),
703 std::sqrt((Ct2av - gmx::square(Ctav)) / (nitem - 1)));
706 if (fp)
708 xvgrclose(fp);
710 sfree(fit);
713 /*! \brief Legend for selecting Legendre polynomials. */
714 static const char* Leg[] = { nullptr, "0", "1", "2", "3", nullptr };
716 t_pargs* add_acf_pargs(int* npargs, t_pargs* pa)
718 t_pargs acfpa[] = {
719 { "-acflen",
720 FALSE,
721 etINT,
722 { &acf.nout },
723 "Length of the ACF, default is half the number of frames" },
724 { "-normalize", FALSE, etBOOL, { &acf.bNormalize }, "Normalize ACF" },
725 { "-fftcorr",
726 FALSE,
727 etBOOL,
728 { &acf.bFour },
729 "HIDDENUse fast fourier transform for correlation function" },
730 { "-nrestart",
731 FALSE,
732 etINT,
733 { &acf.nrestart },
734 "HIDDENNumber of frames between time origins for ACF when no FFT is used" },
735 { "-P", FALSE, etENUM, { Leg }, "Order of Legendre polynomial for ACF (0 indicates none)" },
736 { "-fitfn", FALSE, etENUM, { s_ffn }, "Fit function" },
737 { "-beginfit",
738 FALSE,
739 etREAL,
740 { &acf.tbeginfit },
741 "Time where to begin the exponential fit of the correlation function" },
742 { "-endfit",
743 FALSE,
744 etREAL,
745 { &acf.tendfit },
746 "Time where to end the exponential fit of the correlation function, -1 is until the "
747 "end" },
749 t_pargs* ppa;
750 int i, npa;
752 npa = asize(acfpa);
753 snew(ppa, *npargs + npa);
754 for (i = 0; (i < *npargs); i++)
756 ppa[i] = pa[i];
758 for (i = 0; (i < npa); i++)
760 ppa[*npargs + i] = acfpa[i];
762 (*npargs) += npa;
764 acf.mode = 0;
765 acf.nrestart = 1;
766 acf.nout = -1;
767 acf.P = 0;
768 acf.fitfn = effnEXP1;
769 acf.bFour = TRUE;
770 acf.bNormalize = TRUE;
771 acf.tbeginfit = 0.0;
772 acf.tendfit = -1;
774 bACFinit = TRUE;
776 return ppa;
779 void do_autocorr(const char* fn,
780 const gmx_output_env_t* oenv,
781 const char* title,
782 int nframes,
783 int nitem,
784 real** c1,
785 real dt,
786 unsigned long mode,
787 gmx_bool bAver)
789 if (!bACFinit)
791 printf("ACF data structures have not been initialised. Call add_acf_pargs\n");
794 /* Handle enumerated types */
795 sscanf(Leg[0], "%d", &acf.P);
796 acf.fitfn = sffn2effn(s_ffn);
798 switch (acf.P)
800 case 1: mode = mode | eacP1; break;
801 case 2: mode = mode | eacP2; break;
802 case 3: mode = mode | eacP3; break;
803 default: break;
806 low_do_autocorr(fn, oenv, title, nframes, nitem, acf.nout, c1, dt, mode, acf.nrestart, bAver,
807 acf.bNormalize, bDebugMode(), acf.tbeginfit, acf.tendfit, acf.fitfn);
810 int get_acfnout()
812 if (!bACFinit)
814 gmx_fatal(FARGS, "ACF data not initialized yet");
817 return acf.nout;
820 int get_acffitfn()
822 if (!bACFinit)
824 gmx_fatal(FARGS, "ACF data not initialized yet");
827 return sffn2effn(s_ffn);