Implements function to compute many autocorrelation functions.
- Author
- David van der Spoel david.nosp@m..van.nosp@m.dersp.nosp@m.oel@.nosp@m.icm.u.nosp@m.u.se
|
static void | low_do_four_core (int nframes, real c1[], real cfour[], int nCos) |
| Routine to compute ACF using FFT.
|
|
static void | do_ac_core (int nframes, int nout, real corr[], real c1[], int nrestart, unsigned long mode) |
| Routine to comput ACF without FFT.
|
|
static void | normalize_acf (int nout, real corr[]) |
| Routine to normalize ACF, dividing by corr[0].
|
|
static void | average_acf (gmx_bool bVerbose, int n, int nitem, real **c1) |
| Routine that averages ACFs.
|
|
static void | norm_and_scale_vectors (int nframes, real c1[], real scale) |
| Normalize ACFs.
|
|
static void | dump_tmp (char *s, int n, real c[]) |
| Debugging.
|
|
static void | do_four_core (unsigned long mode, int nframes, real c1[], real csum[], real ctmp[]) |
| High level ACF routine.
|
|
void | low_do_autocorr (const char *fn, const gmx_output_env_t *oenv, const char *title, int nframes, int nitem, int nout, real **c1, real dt, unsigned long mode, int nrestart, gmx_bool bAver, gmx_bool bNormalize, gmx_bool bVerbose, real tbeginfit, real tendfit, int eFitFn) |
| Low level computation of autocorrelation functions. More...
|
|
t_pargs * | add_acf_pargs (int *npargs, t_pargs *pa) |
| Add commandline arguments related to autocorrelations to the existing array. *npargs must be initialised to the number of elements in pa, it will be incremented appropriately. More...
|
|
void | do_autocorr (const char *fn, const gmx_output_env_t *oenv, const char *title, int nframes, int nitem, real **c1, real dt, unsigned long mode, gmx_bool bAver) |
| Calls low_do_autocorr (see below). add_acf_pargs has to be called before this can be used. More...
|
|
int | get_acfnout () |
| Returns the number of points to output from a correlation function. Works only AFTER do_auto_corr has been called! More...
|
|
int | get_acffitfn () |
| Returns the fitting function selected. Works only AFTER do_auto_corr has been called! More...
|
|
void low_do_autocorr |
( |
const char * |
fn, |
|
|
const gmx_output_env_t * |
oenv, |
|
|
const char * |
title, |
|
|
int |
nframes, |
|
|
int |
nitem, |
|
|
int |
nout, |
|
|
real ** |
c1, |
|
|
real |
dt, |
|
|
unsigned long |
mode, |
|
|
int |
nrestart, |
|
|
gmx_bool |
bAver, |
|
|
gmx_bool |
bNormalize, |
|
|
gmx_bool |
bVerbose, |
|
|
real |
tbeginfit, |
|
|
real |
tendfit, |
|
|
int |
nfitparm |
|
) |
| |
Low level computation of autocorrelation functions.
do_autocorr calculates autocorrelation functions for many things. It takes a 2 d array containing nitem arrays of length nframes for each item the ACF is calculated.
A number of "modes" exist for computation of the ACF controlled by variable mode, with the following meaning.
Mode | Function |
eacNormal | C(t) = < X (tau) * X (tau+t) > |
eacCos | C(t) = < cos (X(tau) - X(tau+t)) > |
eacIden | C(t) = < (X(tau) == X(tau+t)) > (not fully supported yet) |
eacVector | C(t) = < X(tau) * X(tau+t) |
eacP1 | C(t) = < cos (X(tau) * X(tau+t) > |
eacP2 | C(t) = 1/2 * < 3 cos (X(tau) * X(tau+t) - 1 > |
eacRcross | C(t) = < ( X(tau) * X(tau+t) )^2 > |
For modes eacVector, eacP1, eacP2 and eacRcross the input should be 3 x nframes long, where each triplet is taken as a 3D vector
For mode eacCos inputdata must be in radians, not degrees!
Other parameters are:
- Parameters
-
[in] | fn | is output filename (.xvg) where the correlation function(s) are printed |
[in] | oenv | controls output file properties |
[in] | title | is the title in the output file |
[in] | nframes | is the number of frames in the time series |
[in] | nitem | is the number of items |
[in] | nout | |
[in,out] | c1 | is an array of dimension [ 0 .. nitem-1 ] [ 0 .. nframes-1 ] on output, this array is filled with the correlation function to reduce storage |
[in] | dt | is the time between frames |
[in] | mode | Different types of ACF can be done, see above |
[in] | nrestart | is the number of steps between restarts for direct ACFs (i.e. without FFT) When set to 1 all points are used as time origin for averaging |
[in] | bAver | If set, all ndih C(t) functions are averaged into a single C(t) |
[in] | bNormalize | If set, all ACFs will be normalized to start at 0 |
[in] | bVerbose | If set output to console will be generated |
[in] | tbeginfit | Time to start fitting to the ACF |
[in] | tendfit | Time to end fitting to the ACF |
[in] | nfitparm | Number of fitting parameters in a multi-exponential fit |