Compile nonbonded kernels as C++
[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,2015,2018, 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/utility/basedefinitions.h"
48 #include "gromacs/utility/real.h"
49
50 struct gmx_output_env_t;
51
52 /*! \brief
53  * Enum to select fitting functions
54  */
55 enum {
56     effnNONE, effnEXP1, effnEXP2, effnEXPEXP,
57     effnEXP5, effnEXP7, effnEXP9,
58     effnVAC,  effnERF,  effnERREST, effnPRES, effnNR
59 };
60
61 /*! \brief
62  * Short name of each function type.
63  * This is exported for now in order to use when
64  * calling parse_common_args.
65  */
66 extern const char *s_ffn[effnNR+2];
67
68 /*! \brief
69  * Returns  description corresponding to the enum above, or NULL if out of range
70  * \param[in] effn Index
71  * \return Description or NULL
72  */
73 const char *effnDescription(int effn);
74
75 /*! \brief
76  * Returns  number of function parameters associated with a fitting function.
77  * \param[in] effn Index
78  * \return number or -1 if index out of range
79  */
80 int effnNparams(int effn);
81
82 /*! \brief
83  * Returns  corresponding to the selected enum option in sffn
84  * \param[in] sffn Two dimensional string array coming from parse_common_args
85  * \return the ffn enum
86  */
87 int sffn2effn(const char **sffn);
88
89 /*! \brief
90  * Returns the value of fit function eFitFn at x
91  * \param[in] eFitFn the index to the fitting function (0 .. effnNR)
92  * \param[in] parm Array of parameters, the length of which depends on eFitFn
93  * \param[in] x The value of x
94  * \return the value of the fit
95  */
96 double fit_function(int eFitFn, const double parm[], double x);
97
98 /*! \brief
99  * Use Levenberg-Marquardt method to fit to a nfitparm parameter exponential
100  * or to a transverse current autocorrelation function.
101  *
102  * If x == NULL, the timestep dt will be used to create a time axis.
103  * \param[in] ndata Number of data points
104  * \param[in] c1 The data points
105  * \param[in] sig The standard deviation in the points (can be NULL)
106  * \param[in] dt The time step
107  * \param[in] x0 The X-axis (may be NULL, see above)
108  * \param[in] begintimefit Starting time for fitting
109  * \param[in] endtimefit Ending time for fitting
110  * \param[in] oenv Output formatting information
111  * \param[in] bVerbose Should the routine write to console?
112  * \param[in] eFitFn Fitting function (0 .. effnNR)
113  * \param[inout] fitparms[] Fitting parameters, see printed manual for a
114  * detailed description. Note that in all implemented cases the parameters
115  * corresponding to time constants will be generated with increasing values.
116  * Such input parameters should therefore be provided in increasing order.
117  * If this is not the case or if subsequent time parameters differ by less than
118  * a factor of 2, they will be modified to ensure tau_i+1 >= 2 tau_i.
119  * \param[in] fix Constrains fit parameter i at it's starting value, when the i'th bit
120  * of fix is set. This works only when the N last parameters are fixed
121  * but not when a parameter somewhere in the middle needs to be fixed.
122  * \param[in] fn_fitted If not NULL file to print the data and fitted curve to
123  * \return integral.
124  */
125 real do_lmfit(int ndata, const real c1[], real sig[], real dt, const real *x0,
126               real begintimefit, real endtimefit, const gmx_output_env_t *oenv,
127               gmx_bool bVerbose, int eFitFn, double fitparms[], int fix,
128               const char *fn_fitted);
129
130 /*! \brief
131  * Fit an autocorrelation function to a pre-defined functional form
132  *
133  * \todo check parameters
134  * \param[in] ncorr
135  * \param[in] fitfn Fitting function (0 .. effnNR)
136  * \param[in] oenv Output formatting information
137  * \param[in] bVerbose Should the routine write to console?
138  * \param[in] tbeginfit Starting time for fitting
139  * \param[in] tendfit Ending time for fitting
140  * \param[in] dt The time step
141  * \param[in] c1 The data points
142  * \param[inout] fit The fitting parameters
143  * \return the integral over the autocorrelation function?
144  */
145 real fit_acf(int ncorr, int fitfn, const gmx_output_env_t *oenv, gmx_bool bVerbose,
146              real tbeginfit, real tendfit, real dt, real c1[], real *fit);
147
148 #endif