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.
39 * Implements function to compute many autocorrelation functions
41 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
42 * \ingroup module_correlationfunctions
52 #include "gromacs/correlationfunctions/expfit.h"
53 #include "gromacs/correlationfunctions/integrate.h"
54 #include "gromacs/correlationfunctions/manyautocorrelation.h"
55 #include "gromacs/correlationfunctions/polynomials.h"
56 #include "gromacs/fileio/xvgr.h"
57 #include "gromacs/legacyheaders/macros.h"
58 #include "gromacs/legacyheaders/names.h"
59 #include "gromacs/math/vec.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/futil.h"
62 #include "gromacs/utility/real.h"
63 #include "gromacs/utility/smalloc.h"
65 /*! \brief Shortcut macro to select modes. */
66 #define MODE(x) ((mode & (x)) == (x))
70 int nrestart
, nout
, P
, fitfn
;
71 gmx_bool bFour
, bNormalize
;
72 real tbeginfit
, tendfit
;
75 /*! \brief Global variable set true if initialization routines are called. */
76 static gmx_bool bACFinit
= FALSE
;
78 /*! \brief Data structure for storing command line variables. */
85 /*! \brief Routine to compute ACF using FFT. */
86 static void low_do_four_core(int nfour
, int nframes
, real c1
[], real cfour
[],
97 for (i
= 0; (i
< nframes
); i
++)
104 for (i
= 0; (i
< nframes
); i
++)
106 cfour
[i
] = cos(c1
[i
]);
110 for (i
= 0; (i
< nframes
); i
++)
112 cfour
[i
] = sin(c1
[i
]);
116 gmx_fatal(FARGS
, "nCos = %d, %s %d", nCos
, __FILE__
, __LINE__
);
119 fftcode
= many_auto_correl(1, nframes
, nfour
, &cfour
);
122 /*! \brief Routine to comput ACF without FFT. */
123 static void do_ac_core(int nframes
, int nout
,
124 real corr
[], real c1
[], int nrestart
,
127 int j
, k
, j3
, jk3
, m
, n
;
133 printf("WARNING: setting number of restarts to 1\n");
139 "Starting do_ac_core: nframes=%d, nout=%d, nrestart=%d,mode=%lu\n",
140 nframes
, nout
, nrestart
, mode
);
143 for (j
= 0; (j
< nout
); j
++)
148 /* Loop over starting points. */
149 for (j
= 0; (j
< nframes
); j
+= nrestart
)
153 /* Loop over the correlation length for this starting point */
154 for (k
= 0; (k
< nout
) && (j
+k
< nframes
); k
++)
158 /* Switch over possible ACF types.
159 * It might be more efficient to put the loops inside the switch,
160 * but this is more clear, and save development time!
164 corr
[k
] += c1
[j
]*c1
[j
+k
];
166 else if (MODE(eacCos
))
168 /* Compute the cos (phi(t)-phi(t+dt)) */
169 corr
[k
] += cos(c1
[j
]-c1
[j
+k
]);
171 else if (MODE(eacIden
))
173 /* Check equality (phi(t)==phi(t+dt)) */
174 corr
[k
] += (c1
[j
] == c1
[j
+k
]) ? 1 : 0;
176 else if (MODE(eacP1
) || MODE(eacP2
) || MODE(eacP3
))
180 for (m
= 0; (m
< DIM
); m
++)
185 cth
= cos_angle(xj
, xk
);
187 if (cth
-1.0 > 1.0e-15)
189 printf("j: %d, k: %d, xj:(%g,%g,%g), xk:(%g,%g,%g)\n",
190 j
, k
, xj
[XX
], xj
[YY
], xj
[ZZ
], xk
[XX
], xk
[YY
], xk
[ZZ
]);
197 else if (MODE(eacP3
))
201 corr
[k
] += LegendreP(cth
, mmm
); /* 1.5*cth*cth-0.5; */
203 else if (MODE(eacRcross
))
206 for (m
= 0; (m
< DIM
); m
++)
213 corr
[k
] += iprod(rr
, rr
);
215 else if (MODE(eacVector
))
217 for (m
= 0; (m
< DIM
); m
++)
228 gmx_fatal(FARGS
, "\nInvalid mode (%d) in do_ac_core", mode
);
232 /* Correct for the number of points and copy results to the data array */
233 for (j
= 0; (j
< nout
); j
++)
235 n
= (nframes
-j
+(nrestart
-1))/nrestart
;
240 /*! \brief Routine to normalize ACF, dividing by corr[0]. */
241 void normalize_acf(int nout
, real corr
[])
248 fprintf(debug
, "Before normalization\n");
249 for (j
= 0; (j
< nout
); j
++)
251 fprintf(debug
, "%5d %10f\n", j
, corr
[j
]);
255 /* Normalisation makes that c[0] = 1.0 and that other points are scaled
258 if (fabs(corr
[0]) < 1e-5)
266 for (j
= 0; (j
< nout
); j
++)
273 fprintf(debug
, "After normalization\n");
274 for (j
= 0; (j
< nout
); j
++)
276 fprintf(debug
, "%5d %10f\n", j
, corr
[j
]);
281 /*! \brief Routine that averages ACFs. */
282 void average_acf(gmx_bool bVerbose
, int n
, int nitem
, real
**c1
)
289 printf("Averaging correlation functions\n");
292 for (j
= 0; (j
< n
); j
++)
295 for (i
= 0; (i
< nitem
); i
++)
303 /*! \brief Normalize ACFs. */
304 void norm_and_scale_vectors(int nframes
, real c1
[], real scale
)
309 for (j
= 0; (j
< nframes
); j
++)
313 for (m
= 0; (m
< DIM
); m
++)
320 /*! \brief Debugging */
321 static void dump_tmp(char *s
, int n
, real c
[])
326 fp
= gmx_ffopen(s
, "w");
327 for (i
= 0; (i
< n
); i
++)
329 fprintf(fp
, "%10d %10g\n", i
, c
[i
]);
334 /*! \brief High level ACF routine. */
335 void do_four_core(unsigned long mode
, int nfour
, int nf2
, int nframes
,
336 real c1
[], real csum
[], real ctmp
[])
347 /********************************************
349 ********************************************/
350 low_do_four_core(nfour
, nf2
, c1
, csum
, enNorm
);
352 else if (MODE(eacCos
))
354 /***************************************************
356 ***************************************************/
357 /* Copy the data to temp array. Since we need it twice
358 * we can't overwrite original.
360 for (j
= 0; (j
< nf2
); j
++)
365 /* Cosine term of AC function */
366 low_do_four_core(nfour
, nf2
, ctmp
, cfour
, enCos
);
367 for (j
= 0; (j
< nf2
); j
++)
372 /* Sine term of AC function */
373 low_do_four_core(nfour
, nf2
, ctmp
, cfour
, enSin
);
374 for (j
= 0; (j
< nf2
); j
++)
380 else if (MODE(eacP2
))
382 /***************************************************
383 * Legendre polynomials
384 ***************************************************/
385 /* First normalize the vectors */
386 norm_and_scale_vectors(nframes
, c1
, 1.0);
388 /* For P2 thingies we have to do six FFT based correls
389 * First for XX^2, then for YY^2, then for ZZ^2
390 * Then we have to do XY, YZ and XZ (counting these twice)
391 * After that we sum them and normalise
392 * P2(x) = (3 * cos^2 (x) - 1)/2
393 * for unit vectors u and v we compute the cosine as the inner product
394 * cos(u,v) = uX vX + uY vY + uZ vZ
398 * C(t) = | (3 cos^2(u(t'),u(t'+t)) - 1)/2 dt'
403 * P2(u(0),u(t)) = [3 * (uX(0) uX(t) +
405 * uZ(0) uZ(t))^2 - 1]/2
406 * = [3 * ((uX(0) uX(t))^2 +
409 * 2(uX(0) uY(0) uX(t) uY(t)) +
410 * 2(uX(0) uZ(0) uX(t) uZ(t)) +
411 * 2(uY(0) uZ(0) uY(t) uZ(t))) - 1]/2
413 * = [(3/2) * (<uX^2> + <uY^2> + <uZ^2> +
414 * 2<uXuY> + 2<uXuZ> + 2<uYuZ>) - 0.5]
418 /* Because of normalization the number of -0.5 to subtract
419 * depends on the number of data points!
421 for (j
= 0; (j
< nf2
); j
++)
423 csum
[j
] = -0.5*(nf2
-j
);
426 /***** DIAGONAL ELEMENTS ************/
427 for (m
= 0; (m
< DIM
); m
++)
429 /* Copy the vector data in a linear array */
430 for (j
= 0; (j
< nf2
); j
++)
432 ctmp
[j
] = sqr(c1
[DIM
*j
+m
]);
436 sprintf(buf
, "c1diag%d.xvg", m
);
437 dump_tmp(buf
, nf2
, ctmp
);
440 low_do_four_core(nfour
, nf2
, ctmp
, cfour
, enNorm
);
444 sprintf(buf
, "c1dfout%d.xvg", m
);
445 dump_tmp(buf
, nf2
, cfour
);
448 for (j
= 0; (j
< nf2
); j
++)
450 csum
[j
] += fac
*(cfour
[j
]);
453 /******* OFF-DIAGONAL ELEMENTS **********/
454 for (m
= 0; (m
< DIM
); m
++)
456 /* Copy the vector data in a linear array */
458 for (j
= 0; (j
< nf2
); j
++)
460 ctmp
[j
] = c1
[DIM
*j
+m
]*c1
[DIM
*j
+m1
];
465 sprintf(buf
, "c1off%d.xvg", m
);
466 dump_tmp(buf
, nf2
, ctmp
);
468 low_do_four_core(nfour
, nf2
, ctmp
, cfour
, enNorm
);
471 sprintf(buf
, "c1ofout%d.xvg", m
);
472 dump_tmp(buf
, nf2
, cfour
);
475 for (j
= 0; (j
< nf2
); j
++)
477 csum
[j
] += fac
*cfour
[j
];
481 else if (MODE(eacP1
) || MODE(eacVector
))
483 /***************************************************
485 ***************************************************/
488 /* First normalize the vectors */
489 norm_and_scale_vectors(nframes
, c1
, 1.0);
492 /* For vector thingies we have to do three FFT based correls
493 * First for XX, then for YY, then for ZZ
494 * After that we sum them and normalise
496 for (j
= 0; (j
< nf2
); j
++)
500 for (m
= 0; (m
< DIM
); m
++)
502 /* Copy the vector data in a linear array */
503 for (j
= 0; (j
< nf2
); j
++)
505 ctmp
[j
] = c1
[DIM
*j
+m
];
507 low_do_four_core(nfour
, nf2
, ctmp
, cfour
, enNorm
);
508 for (j
= 0; (j
< nf2
); j
++)
516 gmx_fatal(FARGS
, "\nUnknown mode in do_autocorr (%d)", mode
);
520 for (j
= 0; (j
< nf2
); j
++)
522 c1
[j
] = csum
[j
]/(real
)(nframes
-j
);
526 void low_do_autocorr(const char *fn
, const output_env_t oenv
, const char *title
,
527 int nframes
, int nitem
, int nout
, real
**c1
,
528 real dt
, unsigned long mode
, int nrestart
,
529 gmx_bool bAver
, gmx_bool bNormalize
,
530 gmx_bool bVerbose
, real tbeginfit
, real tendfit
,
533 FILE *fp
, *gp
= NULL
;
537 real c0
, sum
, Ct2av
, Ctav
;
538 gmx_bool bFour
= acf
.bFour
;
540 /* Check flags and parameters */
541 nout
= get_acfnout();
544 nout
= acf
.nout
= (nframes
+1)/2;
546 else if (nout
> nframes
)
551 if (MODE(eacCos
) && MODE(eacVector
))
553 gmx_fatal(FARGS
, "Incompatible options bCos && bVector (%s, %d)",
556 if ((MODE(eacP3
) || MODE(eacRcross
)) && bFour
)
560 fprintf(stderr
, "Can't combine mode %lu with FFT, turning off FFT\n", mode
);
564 if (MODE(eacNormal
) && MODE(eacVector
))
566 gmx_fatal(FARGS
, "Incompatible mode bits: normal and vector (or Legendre)");
569 /* Print flags and parameters */
572 printf("Will calculate %s of %d thingies for %d frames\n",
573 title
? title
: "autocorrelation", nitem
, nframes
);
574 printf("bAver = %s, bFour = %s bNormalize= %s\n",
575 bool_names
[bAver
], bool_names
[bFour
], bool_names
[bNormalize
]);
576 printf("mode = %lu, dt = %g, nrestart = %d\n", mode
, dt
, nrestart
);
580 c0
= log((double)nframes
)/log(2.0);
590 fprintf(debug
, "Using FFT to calculate %s, #points for FFT = %d\n",
594 /* Allocate temp arrays */
600 nfour
= 0; /* To keep the compiler happy */
605 /* Loop over items (e.g. molecules or dihedrals)
606 * In this loop the actual correlation functions are computed, but without
609 k
= max(1, pow(10, (int)(log(nitem
)/log(100))));
610 for (i
= 0; i
< nitem
; i
++)
612 if (bVerbose
&& ((i
%k
== 0 || i
== nitem
-1)))
614 fprintf(stderr
, "\rThingie %d", i
+1);
619 do_four_core(mode
, nfour
, nframes
, nframes
, c1
[i
], csum
, ctmp
);
623 do_ac_core(nframes
, nout
, ctmp
, c1
[i
], nrestart
, mode
);
628 fprintf(stderr
, "\n");
636 fp
= xvgropen(fn
, title
, "Time (ps)", "C(t)", oenv
);
647 average_acf(bVerbose
, nframes
, nitem
, c1
);
652 normalize_acf(nout
, c1
[0]);
655 if (eFitFn
!= effnNONE
)
657 fit_acf(nout
, eFitFn
, oenv
, fn
!= NULL
, tbeginfit
, tendfit
, dt
, c1
[0], fit
);
658 sum
= print_and_integrate(fp
, nout
, dt
, c1
[0], fit
, 1);
662 sum
= print_and_integrate(fp
, nout
, dt
, c1
[0], NULL
, 1);
665 printf("Correlation time (integral over corrfn): %g (ps)\n", sum
);
671 /* Not averaging. Normalize individual ACFs */
675 gp
= xvgropen("ct-distr.xvg", "Correlation times", "item", "time (ps)", oenv
);
677 for (i
= 0; i
< nitem
; i
++)
681 normalize_acf(nout
, c1
[i
]);
683 if (eFitFn
!= effnNONE
)
685 fit_acf(nout
, eFitFn
, oenv
, fn
!= NULL
, tbeginfit
, tendfit
, dt
, c1
[i
], fit
);
686 sum
= print_and_integrate(fp
, nout
, dt
, c1
[i
], fit
, 1);
690 sum
= print_and_integrate(fp
, nout
, dt
, c1
[i
], NULL
, 1);
694 "CORRelation time (integral over corrfn %d): %g (ps)\n",
702 fprintf(gp
, "%5d %.3f\n", i
, sum
);
713 printf("Average correlation time %.3f Std. Dev. %.3f Error %.3f (ps)\n",
714 Ctav
, sqrt((Ct2av
- sqr(Ctav
))),
715 sqrt((Ct2av
- sqr(Ctav
))/(nitem
-1)));
725 /*! \brief Legend for selecting Legendre polynomials. */
726 static const char *Leg
[] = { NULL
, "0", "1", "2", "3", NULL
};
728 t_pargs
*add_acf_pargs(int *npargs
, t_pargs
*pa
)
731 { "-acflen", FALSE
, etINT
, {&acf
.nout
},
732 "Length of the ACF, default is half the number of frames" },
733 { "-normalize", FALSE
, etBOOL
, {&acf
.bNormalize
},
735 { "-fftcorr", FALSE
, etBOOL
, {&acf
.bFour
},
736 "HIDDENUse fast fourier transform for correlation function" },
737 { "-nrestart", FALSE
, etINT
, {&acf
.nrestart
},
738 "HIDDENNumber of frames between time origins for ACF when no FFT is used" },
739 { "-P", FALSE
, etENUM
, {Leg
},
740 "Order of Legendre polynomial for ACF (0 indicates none)" },
741 { "-fitfn", FALSE
, etENUM
, {s_ffn
},
743 { "-beginfit", FALSE
, etREAL
, {&acf
.tbeginfit
},
744 "Time where to begin the exponential fit of the correlation function" },
745 { "-endfit", FALSE
, etREAL
, {&acf
.tendfit
},
746 "Time where to end the exponential fit of the correlation function, -1 is until the end" },
752 snew(ppa
, *npargs
+npa
);
753 for (i
= 0; (i
< *npargs
); i
++)
757 for (i
= 0; (i
< npa
); i
++)
759 ppa
[*npargs
+i
] = acfpa
[i
];
767 acf
.fitfn
= effnEXP1
;
769 acf
.bNormalize
= TRUE
;
778 void do_autocorr(const char *fn
, const output_env_t oenv
, const char *title
,
779 int nframes
, int nitem
, real
**c1
,
780 real dt
, unsigned long mode
, gmx_bool bAver
)
784 printf("ACF data structures have not been initialised. Call add_acf_pargs\n");
787 /* Handle enumerated types */
788 sscanf(Leg
[0], "%d", &acf
.P
);
789 acf
.fitfn
= sffn2effn(s_ffn
);
806 low_do_autocorr(fn
, oenv
, title
, nframes
, nitem
, acf
.nout
, c1
, dt
, mode
,
807 acf
.nrestart
, bAver
, acf
.bNormalize
,
808 bDebugMode(), acf
.tbeginfit
, acf
.tendfit
,
812 int get_acfnout(void)
816 gmx_fatal(FARGS
, "ACF data not initialized yet");
822 int get_acffitfn(void)
826 gmx_fatal(FARGS
, "ACF data not initialized yet");
829 return sffn2effn(s_ffn
);