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