Split off the NMR related analyses from gmx energy.
[gromacs.git] / src / gromacs / correlationfunctions / autocorr.h
blob8e5c67c1d9d5452c7206e9ba0b4a808d4d81f793
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015, 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.
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.
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.
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.
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.
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.
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.
41 /*! \libinternal
42 * \file
43 * \brief
44 * Declares routine for computing autocorrelation functions
46 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
47 * \inlibraryapi
48 * \ingroup module_correlationfunctions
50 #ifndef GMX_AUTOCORR_H
51 #define GMX_AUTOCORR_H
53 #include "gromacs/commandline/pargs.h"
54 #include "gromacs/utility/real.h"
56 struct gmx_output_env_t;
58 #ifdef __cplusplus
59 extern "C" {
60 #endif
62 /*! \brief Normal correlation f(t)*f(t+dt) */
63 #define eacNormal (1<<0)
64 /*! \brief Cosine correlation cos(f(t)-f(t+dt)) */
65 #define eacCos (1<<1)
66 /*! \brief Vector correlation f(t).f(t+dt) */
67 #define eacVector (1<<2)
68 /*! \brief Norm of cross product |f(t) (x) f(t+dt)| */
69 #define eacRcross (1<<3 | eacVector)
70 /*! \brief Vector with Legendre polynomial of order 0 (same as vector) */
71 #define eacP0 (1<<4 | eacVector)
72 /*! \brief Vector with Legendre polynomial of order P_1(f(t).f(t+dt)) */
73 #define eacP1 (1<<5 | eacVector)
74 /*! \brief Vector with Legendre polynomial of order P_2(f(t).f(t+dt)) */
75 #define eacP2 (1<<6 | eacVector)
76 /*! \brief Vector with Legendre polynomial of order P_3(f(t).f(t+dt)) */
77 #define eacP3 (1<<7 | eacVector)
78 /*! \brief Vector with Legendre polynomial of order P_4(f(t).f(t+dt)) */
79 #define eacP4 (1<<8 | eacVector)
80 /*! \brief Binary identy correlation (f(t) == f(t+dt)) */
81 #define eacIden (1<<9) //Not supported for multiple cores
83 /*! \brief
84 * Add commandline arguments related to autocorrelations to the existing array.
85 * *npargs must be initialised to the number of elements in pa,
86 * it will be incremented appropriately.
88 * \param npargs The number of arguments before and after (is changed in this function)
89 * \param[in] pa The initial argument list
90 * \return the new array
92 t_pargs *add_acf_pargs(int *npargs, t_pargs *pa);
94 /*! \brief
95 * Returns the number of points to output from a correlation function.
96 * Works only AFTER do_auto_corr has been called!
97 * \return the output length for the correlation function
99 int get_acfnout(void);
101 /*! \brief
102 * Returns the fitting function selected.
103 * Works only AFTER do_auto_corr has been called!
104 * \return the fit function type.
106 int get_acffitfn(void);
108 /*! \brief
109 * Calls low_do_autocorr (see below). add_acf_pargs has to be called before this
110 * can be used.
111 * \param[in] fn File name for xvg output (may this be NULL)?
112 * \param[in] oenv The output environment information
113 * \param[in] title is the title in the output file
114 * \param[in] nframes is the number of frames in the time series
115 * \param[in] nitem is the number of items
116 * \param[inout] c1 is an array of dimension [ 0 .. nitem-1 ] [ 0 .. nframes-1 ]
117 * on output, this array is filled with the correlation function
118 * to reduce storage
119 * \param[in] dt is the time between frames
120 * \param[in] mode Different types of ACF can be done, see above
121 * \param[in] bAver If set, all ndih C(t) functions are averaged into a single
122 * C(t)
124 void do_autocorr(const char *fn, const gmx_output_env_t *oenv,
125 const char *title,
126 int nframes, int nitem, real **c1,
127 real dt, unsigned long mode, gmx_bool bAver);
129 /*! \brief
130 * Low level computation of autocorrelation functions
132 * do_autocorr calculates autocorrelation functions for many things.
133 * It takes a 2 d array containing nitem arrays of length nframes
134 * for each item the ACF is calculated.
136 * A number of "modes" exist for computation of the ACF controlled
137 * by variable mode, with the following meaning.
139 * Mode | Function
140 * -----------|------------
141 * eacNormal | C(t) = < X (tau) * X (tau+t) >
142 * eacCos | C(t) = < cos (X(tau) - X(tau+t)) >
143 * eacIden | C(t) = < (X(tau) == X(tau+t)) > (not fully supported yet)
144 * eacVector | C(t) = < X(tau) * X(tau+t)
145 * eacP1 | C(t) = < cos (X(tau) * X(tau+t) >
146 * eacP2 | C(t) = 1/2 * < 3 cos (X(tau) * X(tau+t) - 1 >
147 * eacRcross | C(t) = < ( X(tau) * X(tau+t) )^2 >
149 * For modes eacVector, eacP1, eacP2 and eacRcross the input should be
150 * 3 x nframes long, where each triplet is taken as a 3D vector
152 * For mode eacCos inputdata must be in radians, not degrees!
154 * Other parameters are:
156 * \param[in] fn is output filename (.xvg) where the correlation function(s) are printed
157 * \param[in] oenv controls output file properties
158 * \param[in] title is the title in the output file
159 * \param[in] nframes is the number of frames in the time series
160 * \param[in] nitem is the number of items
161 * \param[in] nout
162 * \param[inout] c1 is an array of dimension [ 0 .. nitem-1 ] [ 0 .. nframes-1 ]
163 * on output, this array is filled with the correlation function
164 * to reduce storage
165 * \param[in] dt is the time between frames
166 * \param[in] mode Different types of ACF can be done, see above
167 * \param[in] nrestart is the number of steps between restarts for direct ACFs
168 * (i.e. without FFT) When set to 1 all points are used as
169 * time origin for averaging
170 * \param[in] bAver If set, all ndih C(t) functions are averaged into a single
171 * C(t)
172 * \param[in] bNormalize If set, all ACFs will be normalized to start at 0
173 * \param[in] bVerbose If set output to console will be generated
174 * \param[in] tbeginfit Time to start fitting to the ACF
175 * \param[in] tendfit Time to end fitting to the ACF
176 * \param[in] nfitparm Number of fitting parameters in a multi-exponential fit
178 void low_do_autocorr(const char *fn, const gmx_output_env_t *oenv,
179 const char *title, int nframes, int nitem,
180 int nout, real **c1, real dt, unsigned long mode,
181 int nrestart, gmx_bool bAver, gmx_bool bNormalize,
182 gmx_bool bVerbose, real tbeginfit, real tendfit,
183 int nfitparm);
185 #ifdef __cplusplus
187 #endif
189 #endif