Reorganizing analysis of correlation functions.
[alexxy/gromacs.git] / src / gromacs / correlationfunctions / expfit.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  * \file
37  * \brief
38  * Declares routine for fitting a data set to a curve
39  *
40  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
41  * \inlibraryapi
42  * \ingroup module_correlationfunctions
43  */
44 #ifndef GMX_EXPFIT_H
45 #define GMX_EXPFIT_H
46
47 #include "gromacs/legacyheaders/oenv.h"
48 #include "gromacs/utility/real.h"
49
50 #ifdef __cplusplus
51 extern "C" {
52 #endif
53
54 /*! \brief
55  * Enum to select fitting functions
56  */
57 enum {
58     effnNONE, effnEXP1, effnEXP2, effnEXP3,   effnVAC,
59     effnEXP5, effnEXP7, effnEXP9, effnERF, effnERREST, effnPRES, effnNR
60 };
61
62 /*! \brief
63  * Short name of each function type.
64  * This is exported for now in order to use when
65  * calling parse_common_args.
66  */
67 extern const char *s_ffn[effnNR+2];
68
69 /*! \brief
70  * Returns  description corresponding to the enum above, or NULL if out of range
71  * \param[in] effn Index
72  * \return Description or NULL
73  */
74 const char *effnDescription(int effn);
75
76 /*! \brief
77  * Returns  number of function parameters associated with a fitting function.
78  * \param[in] effn Index
79  * \return number or -1 if index out of range
80  */
81 int effnNparams(int effn);
82
83 /*! \brief
84  * Returns  corresponding to the selected enum option in sffn
85  * \param[in] sffn Two dimensional string array coming from parse_common_args
86  * \return the ffn enum
87  */
88 int sffn2effn(const char **sffn);
89
90 /*! \brief
91  * Returns the value of fit function eFitFn at x
92  * \param[in] eFitFn the index to the fitting function (0 .. effnNR)
93  * \param[in] parm Array of parameters, the length of which depends on eFitFn
94  * \param[in] x The value of x
95  * \return the value of the fit
96  */
97 double fit_function(int eFitFn, double *parm, double x);
98
99 /*! \brief
100  * Use Levenberg-Marquardt method to fit to a nfitparm parameter exponential
101  * or to a transverse current autocorrelation function.
102  *
103  * If x == NULL, the timestep dt will be used to create a time axis.
104  * \param[in] ndata Number of data points
105  * \param[in] c1 The data points
106  * \param[in] sig The standard deviation in the points
107  * \param[in] dt The time step
108  * \param[in] x The X-axis (may be NULL, see above)
109  * \param[in] begintimefit Starting time for fitting
110  * \param[in] endtimefit Ending time for fitting
111  * \param[in] oenv Output formatting information
112  * \param[in] bVerbose Should the routine write to console?
113  * \param[in] eFitFn Fitting function (0 .. effnNR)
114  * \param[out] fitparms[]
115  * \param[in] fix Constrains fit parameter i at it's starting value, when the i'th bit
116  * of fix is set.
117  * \return integral.
118  */
119 real do_lmfit(int ndata, real c1[], real sig[], real dt, real *x,
120               real begintimefit, real endtimefit, const output_env_t oenv,
121               gmx_bool bVerbose, int eFitFn, double fitparms[], int fix);
122
123 /*! \brief
124  * Fit an autocorrelation function to a pre-defined functional form
125  *
126  * \todo check parameters
127  * \param[in] ncorr
128  * \param[in] fitfn Fitting function (0 .. effnNR)
129  * \param[in] oenv Output formatting information
130  * \param[in] bVerbose Should the routine write to console?
131  * \param[in] tbeginfit Starting time for fitting
132  * \param[in] tendfit Ending time for fitting
133  * \param[in] dt The time step
134  * \param[in] c1 The data points
135  * \param[out] fit The fitting parameters
136  * \return the integral over the autocorrelation function?
137  */
138 real fit_acf(int ncorr, int fitfn, const output_env_t oenv, gmx_bool bVerbose,
139              real tbeginfit, real tendfit, real dt, real c1[], real *fit);
140
141 #ifdef __cplusplus
142 }
143 #endif
144
145 #endif