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