Gromacs
2018.8
|
Declares routine for computing autocorrelation functions.
Macros | |
#define | eacNormal (1<<0) |
Normal correlation f(t)*f(t+dt) | |
#define | eacCos (1<<1) |
Cosine correlation cos(f(t)-f(t+dt)) | |
#define | eacVector (1<<2) |
Vector correlation f(t).f(t+dt) | |
#define | eacRcross (1<<3 | eacVector) |
Norm of cross product |f(t) (x) f(t+dt)|. | |
#define | eacP0 (1<<4 | eacVector) |
Vector with Legendre polynomial of order 0 (same as vector) | |
#define | eacP1 (1<<5 | eacVector) |
Vector with Legendre polynomial of order P_1(f(t).f(t+dt)) | |
#define | eacP2 (1<<6 | eacVector) |
Vector with Legendre polynomial of order P_2(f(t).f(t+dt)) | |
#define | eacP3 (1<<7 | eacVector) |
Vector with Legendre polynomial of order P_3(f(t).f(t+dt)) | |
#define | eacP4 (1<<8 | eacVector) |
Vector with Legendre polynomial of order P_4(f(t).f(t+dt)) | |
#define | eacIden (1<<9) |
Binary identy correlation (f(t) == f(t+dt)) | |
Functions | |
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... | |
int | get_acfnout (void) |
Returns the number of points to output from a correlation function. Works only AFTER do_auto_corr has been called! More... | |
int | get_acffitfn (void) |
Returns the fitting function selected. Works only AFTER do_auto_corr has been called! 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... | |
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. More... | |
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.
npargs | The number of arguments before and after (is changed in this function) | |
[in] | pa | The initial argument list |
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.
[in] | fn | File name for xvg output (may this be NULL)? |
[in] | oenv | The output environment information |
[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,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] | bAver | If set, all ndih C(t) functions are averaged into a single C(t) |
int get_acffitfn | ( | void | ) |
Returns the fitting function selected. Works only AFTER do_auto_corr has been called!
int get_acfnout | ( | void | ) |
Returns the number of points to output from a correlation function. Works only AFTER do_auto_corr has been called!
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:
[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 |