Reorganizing analysis of correlation functions.
[alexxy/gromacs.git] / src / gromacs / correlationfunctions / autocorr.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \libinternal
36  * \defgroup module_correlationfunctions Correlation functions
37  * \ingroup group_analysismodules
38  * \brief
39  * Compute correlation functions and fit analytical functions to the result.
40  */
41 /*! \libinternal
42  * \file
43  * \brief
44  * Declares routine for computing autocorrelation functions
45  *
46  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
47  * \inlibraryapi
48  * \ingroup module_correlationfunctions
49  */
50 #ifndef GMX_AUTOCORR_H
51 #define GMX_AUTOCORR_H
52
53 #include "gromacs/commandline/pargs.h"
54 #include "gromacs/utility/real.h"
55
56 #ifdef __cplusplus
57 extern "C" {
58 #endif
59
60 /*! \brief Normal correlation f(t)*f(t+dt) */
61 #define eacNormal (1<<0)
62 /*! \brief Cosine correlation cos(f(t)-f(t+dt)) */
63 #define eacCos    (1<<1)
64 /*! \brief Vector correlation f(t).f(t+dt) */
65 #define eacVector (1<<2)
66 /*! \brief Norm of cross product |f(t) (x) f(t+dt)| */
67 #define eacRcross (1<<3  | eacVector)
68 /*! \brief Vector with Legendre polynomial of order 0 (same as vector) */
69 #define eacP0     (1<<4  | eacVector)
70 /*! \brief Vector with Legendre polynomial of order P_1(f(t).f(t+dt)) */
71 #define eacP1     (1<<5  | eacVector)
72 /*! \brief Vector with Legendre polynomial of order P_2(f(t).f(t+dt)) */
73 #define eacP2     (1<<6  | eacVector)
74 /*! \brief Vector with Legendre polynomial of order P_3(f(t).f(t+dt)) */
75 #define eacP3     (1<<7  | eacVector)
76 /*! \brief Vector with Legendre polynomial of order P_4(f(t).f(t+dt)) */
77 #define eacP4     (1<<8  | eacVector)
78 /*! \brief Binary identy correlation (f(t) == f(t+dt)) */
79 #define eacIden   (1<<9) //Not supported for multiple cores
80
81 /*! \brief
82  * Add commandline arguments related to autocorrelations to the existing array.
83  * *npargs must be initialised to the number of elements in pa,
84  * it will be incremented appropriately.
85  *
86  * \param npargs The number of arguments before and after (is changed in this function)
87  * \param[in] pa The initial argument list
88  * \return the new array
89  */
90 t_pargs *add_acf_pargs(int *npargs, t_pargs *pa);
91
92 /*! \brief
93  * Returns the number of points to output from a correlation function.
94  * Works only AFTER do_auto_corr has been called!
95  * \return the output length for the correlation function
96  */
97 int get_acfnout(void);
98
99 /*! \brief
100  * Returns the fitting function selected.
101  * Works only AFTER do_auto_corr has been called!
102  * \return the fit function type.
103  */
104 int get_acffitfn(void);
105
106 /*! \brief
107  * Calls low_do_autocorr (see below). add_acf_pargs has to be called before this
108  * can be used.
109  * \param[in] fn File name for xvg output (may this be NULL)?
110  * \param[in] oenv The output environment information
111  * \param[in] title is the title in the output file
112  * \param[in] nframes is the number of frames in the time series
113  * \param[in] nitem is the number of items
114  * \param[inout] c1 is an array of dimension [ 0 .. nitem-1 ] [ 0 .. nframes-1 ]
115  *          on output, this array is filled with the correlation function
116  *          to reduce storage
117  * \param[in] dt is the time between frames
118  * \param[in] mode Different types of ACF can be done, see above
119  * \param[in] bAver    If set, all ndih C(t) functions are averaged into a single
120  *          C(t)
121  */
122 void do_autocorr(const char *fn, const output_env_t oenv,
123                  const char *title,
124                  int nframes, int nitem, real **c1,
125                  real dt, unsigned long mode, gmx_bool bAver);
126
127 /*! \brief
128  * Low level computation of autocorrelation functions
129  *
130  * do_autocorr calculates autocorrelation functions for many things.
131  * It takes a 2 d array containing nitem arrays of length nframes
132  * for each item the ACF is calculated.
133  *
134  * A number of "modes" exist for computation of the ACF controlled
135  * by variable mode, with the following meaning.
136  *
137  * Mode       | Function
138  * -----------|------------
139  * eacNormal  | C(t) = < X (tau) * X (tau+t) >
140  * eacCos     | C(t) = < cos (X(tau) - X(tau+t)) >
141  * eacIden    | C(t) = < (X(tau) == X(tau+t)) > (not fully supported yet)
142  * eacVector  | C(t) = < X(tau) * X(tau+t)
143  * eacP1      | C(t) = < cos (X(tau) * X(tau+t) >
144  * eacP2      | C(t) = 1/2 * < 3 cos (X(tau) * X(tau+t) - 1 >
145  * eacRcross  | C(t) = < ( X(tau) * X(tau+t) )^2 >
146  *
147  * For modes eacVector, eacP1, eacP2 and eacRcross the input should be
148  * 3 x nframes long, where each triplet is taken as a 3D vector
149  *
150  * For mode eacCos inputdata must be in radians, not degrees!
151  *
152  * Other parameters are:
153  *
154  * \param[in] fn is output filename (.xvg) where the correlation function(s) are printed
155  * \param[in] oenv controls output file properties
156  * \param[in] title is the title in the output file
157  * \param[in] nframes is the number of frames in the time series
158  * \param[in] nitem is the number of items
159  * \param[in] nout
160  * \param[inout] c1 is an array of dimension [ 0 .. nitem-1 ] [ 0 .. nframes-1 ]
161  *          on output, this array is filled with the correlation function
162  *          to reduce storage
163  * \param[in] dt is the time between frames
164  * \param[in] mode Different types of ACF can be done, see above
165  * \param[in] nrestart     is the number of steps between restarts for direct ACFs
166  *              (i.e. without FFT) When set to 1 all points are used as
167  *              time origin for averaging
168  * \param[in] bAver    If set, all ndih C(t) functions are averaged into a single
169  *          C(t)
170  * \param[in] bNormalize   If set, all ACFs will be normalized to start at 0
171  * \param[in] bVerbose If set output to console will be generated
172  * \param[in] tbeginfit Time to start fitting to the ACF
173  * \param[in] tendfit Time to end fitting to the ACF
174  * \param[in] nfitparm Number of fitting parameters in a multi-exponential fit
175  */
176 void low_do_autocorr(const char *fn, const output_env_t oenv,
177                      const char *title, int  nframes, int nitem,
178                      int nout, real **c1, real dt, unsigned long mode,
179                      int nrestart, gmx_bool bAver, gmx_bool bNormalize,
180                      gmx_bool bVerbose, real tbeginfit, real tendfit,
181                      int nfitparm);
182
183 #ifdef __cplusplus
184 }
185 #endif
186
187 #endif