Apply clang-format to source tree
[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,2015,2018,2019, 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 struct gmx_output_env_t;
57
58 /*! \brief Normal correlation f(t)*f(t+dt) */
59 #define eacNormal (1 << 0)
60 /*! \brief Cosine correlation cos(f(t)-f(t+dt)) */
61 #define eacCos (1 << 1)
62 /*! \brief Vector correlation f(t).f(t+dt) */
63 #define eacVector (1 << 2)
64 /*! \brief Norm of cross product |f(t) (x) f(t+dt)| */
65 #define eacRcross (1 << 3 | eacVector)
66 /*! \brief Vector with Legendre polynomial of order 0 (same as vector) */
67 #define eacP0 (1 << 4 | eacVector)
68 /*! \brief Vector with Legendre polynomial of order P_1(f(t).f(t+dt)) */
69 #define eacP1 (1 << 5 | eacVector)
70 /*! \brief Vector with Legendre polynomial of order P_2(f(t).f(t+dt)) */
71 #define eacP2 (1 << 6 | eacVector)
72 /*! \brief Vector with Legendre polynomial of order P_3(f(t).f(t+dt)) */
73 #define eacP3 (1 << 7 | eacVector)
74 /*! \brief Vector with Legendre polynomial of order P_4(f(t).f(t+dt)) */
75 #define eacP4 (1 << 8 | eacVector)
76 /*! \brief Binary identy correlation (f(t) == f(t+dt)) */
77 #define eacIden (1 << 9) // Not supported for multiple cores
78
79 /*! \brief
80  * Add commandline arguments related to autocorrelations to the existing array.
81  * *npargs must be initialised to the number of elements in pa,
82  * it will be incremented appropriately.
83  *
84  * \param npargs The number of arguments before and after (is changed in this function)
85  * \param[in] pa The initial argument list
86  * \return the new array
87  */
88 t_pargs* add_acf_pargs(int* npargs, t_pargs* pa);
89
90 /*! \brief
91  * Returns the number of points to output from a correlation function.
92  * Works only AFTER do_auto_corr has been called!
93  * \return the output length for the correlation function
94  */
95 int get_acfnout();
96
97 /*! \brief
98  * Returns the fitting function selected.
99  * Works only AFTER do_auto_corr has been called!
100  * \return the fit function type.
101  */
102 int get_acffitfn();
103
104 /*! \brief
105  * Calls low_do_autocorr (see below). add_acf_pargs has to be called before this
106  * can be used.
107  * \param[in] fn File name for xvg output (may this be NULL)?
108  * \param[in] oenv The output environment information
109  * \param[in] title is the title in the output file
110  * \param[in] nframes is the number of frames in the time series
111  * \param[in] nitem is the number of items
112  * \param[in,out] c1 is an array of dimension [ 0 .. nitem-1 ] [ 0 .. nframes-1 ]
113  *          on output, this array is filled with the correlation function
114  *          to reduce storage
115  * \param[in] dt is the time between frames
116  * \param[in] mode Different types of ACF can be done, see above
117  * \param[in] bAver    If set, all ndih C(t) functions are averaged into a single
118  *          C(t)
119  */
120 void do_autocorr(const char*             fn,
121                  const gmx_output_env_t* oenv,
122                  const char*             title,
123                  int                     nframes,
124                  int                     nitem,
125                  real**                  c1,
126                  real                    dt,
127                  unsigned long           mode,
128                  gmx_bool                bAver);
129
130 /*! \brief
131  * Low level computation of autocorrelation functions
132  *
133  * do_autocorr calculates autocorrelation functions for many things.
134  * It takes a 2 d array containing nitem arrays of length nframes
135  * for each item the ACF is calculated.
136  *
137  * A number of "modes" exist for computation of the ACF controlled
138  * by variable mode, with the following meaning.
139  *
140  * Mode       | Function
141  * -----------|------------
142  * eacNormal  | C(t) = < X (tau) * X (tau+t) >
143  * eacCos     | C(t) = < cos (X(tau) - X(tau+t)) >
144  * eacIden    | C(t) = < (X(tau) == X(tau+t)) > (not fully supported yet)
145  * eacVector  | C(t) = < X(tau) * X(tau+t)
146  * eacP1      | C(t) = < cos (X(tau) * X(tau+t) >
147  * eacP2      | C(t) = 1/2 * < 3 cos (X(tau) * X(tau+t) - 1 >
148  * eacRcross  | C(t) = < ( X(tau) * X(tau+t) )^2 >
149  *
150  * For modes eacVector, eacP1, eacP2 and eacRcross the input should be
151  * 3 x nframes long, where each triplet is taken as a 3D vector
152  *
153  * For mode eacCos inputdata must be in radians, not degrees!
154  *
155  * Other parameters are:
156  *
157  * \param[in] fn is output filename (.xvg) where the correlation function(s) are printed
158  * \param[in] oenv controls output file properties
159  * \param[in] title is the title in the output file
160  * \param[in] nframes is the number of frames in the time series
161  * \param[in] nitem is the number of items
162  * \param[in] nout
163  * \param[inout] c1 is an array of dimension [ 0 .. nitem-1 ] [ 0 .. nframes-1 ]
164  *          on output, this array is filled with the correlation function
165  *          to reduce storage
166  * \param[in] dt is the time between frames
167  * \param[in] mode Different types of ACF can be done, see above
168  * \param[in] nrestart     is the number of steps between restarts for direct ACFs
169  *              (i.e. without FFT) When set to 1 all points are used as
170  *              time origin for averaging
171  * \param[in] bAver    If set, all ndih C(t) functions are averaged into a single
172  *          C(t)
173  * \param[in] bNormalize   If set, all ACFs will be normalized to start at 0
174  * \param[in] bVerbose If set output to console will be generated
175  * \param[in] tbeginfit Time to start fitting to the ACF
176  * \param[in] tendfit Time to end fitting to the ACF
177  * \param[in] nfitparm Number of fitting parameters in a multi-exponential fit
178  */
179 void low_do_autocorr(const char*             fn,
180                      const gmx_output_env_t* oenv,
181                      const char*             title,
182                      int                     nframes,
183                      int                     nitem,
184                      int                     nout,
185                      real**                  c1,
186                      real                    dt,
187                      unsigned long           mode,
188                      int                     nrestart,
189                      gmx_bool                bAver,
190                      gmx_bool                bNormalize,
191                      gmx_bool                bVerbose,
192                      real                    tbeginfit,
193                      real                    tendfit,
194                      int                     nfitparm);
195
196 #endif