Gromacs  2024.3
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Enumerations | Functions | Variables
expfit.h File Reference
#include "gromacs/utility/basedefinitions.h"
#include "gromacs/utility/real.h"
+ Include dependency graph for expfit.h:
+ This graph shows which files directly or indirectly include this file:

Description

Declares routine for fitting a data set to a curve.

Author
David van der Spoel david.nosp@m..van.nosp@m.dersp.nosp@m.oel@.nosp@m.icm.u.nosp@m.u.se

Enumerations

enum  {
  effnNONE, effnEXP1, effnEXP2, effnEXPEXP,
  effnEXP5, effnEXP7, effnEXP9, effnVAC,
  effnERF, effnERREST, effnPRES, effnNR
}
 Enum to select fitting functions.
 

Functions

const char * effnDescription (int effn)
 Returns description corresponding to the enum above, or NULL if out of range. More...
 
int effnNparams (int effn)
 Returns number of function parameters associated with a fitting function. More...
 
int sffn2effn (const char **sffn)
 Returns corresponding to the selected enum option in sffn. More...
 
double fit_function (int eFitFn, const double parm[], double x)
 Returns the value of fit function eFitFn at x. More...
 
real do_lmfit (int ndata, const real c1[], real sig[], real dt, const real *x0, real begintimefit, real endtimefit, const gmx_output_env_t *oenv, gmx_bool bVerbose, int eFitFn, double fitparms[], int fix, const char *fn_fitted)
 Use Levenberg-Marquardt method to fit to a nfitparm parameter exponential or to a transverse current autocorrelation function. More...
 
real fit_acf (int ncorr, int fitfn, const gmx_output_env_t *oenv, gmx_bool bVerbose, real tbeginfit, real tendfit, real dt, real c1[], real *fit)
 Fit an autocorrelation function to a pre-defined functional form. More...
 

Variables

const char * s_ffn [effnNR+2]
 Short name of each function type. This is exported for now in order to use when calling parse_common_args.
 

Function Documentation

real do_lmfit ( int  ndata,
const real  c1[],
real  sig[],
real  dt,
const real x0,
real  begintimefit,
real  endtimefit,
const gmx_output_env_t *  oenv,
gmx_bool  bVerbose,
int  eFitFn,
double  fitparms[],
int  fix,
const char *  fn_fitted 
)

Use Levenberg-Marquardt method to fit to a nfitparm parameter exponential or to a transverse current autocorrelation function.

If x == NULL, the timestep dt will be used to create a time axis.

Parameters
[in]ndataNumber of data points
[in]c1The data points
[in]sigThe standard deviation in the points (can be NULL)
[in]dtThe time step
[in]x0The X-axis (may be NULL, see above)
[in]begintimefitStarting time for fitting
[in]endtimefitEnding time for fitting
[in]oenvOutput formatting information
[in]bVerboseShould the routine write to console?
[in]eFitFnFitting function (0 .. effnNR)
[in,out]fitparms[]Fitting parameters, see printed manual for a detailed description. Note that in all implemented cases the parameters corresponding to time constants will be generated with increasing values. Such input parameters should therefore be provided in increasing order. If this is not the case or if subsequent time parameters differ by less than a factor of 2, they will be modified to ensure tau_i+1 >= 2 tau_i.
[in]fixConstrains fit parameter i at it's starting value, when the i'th bit of fix is set. This works only when the N last parameters are fixed but not when a parameter somewhere in the middle needs to be fixed.
[in]fn_fittedIf not NULL file to print the data and fitted curve to
Returns
integral.
const char* effnDescription ( int  effn)

Returns description corresponding to the enum above, or NULL if out of range.

Parameters
[in]effnIndex
Returns
Description or NULL
int effnNparams ( int  effn)

Returns number of function parameters associated with a fitting function.

Parameters
[in]effnIndex
Returns
number or -1 if index out of range
real fit_acf ( int  ncorr,
int  fitfn,
const gmx_output_env_t *  oenv,
gmx_bool  bVerbose,
real  tbeginfit,
real  tendfit,
real  dt,
real  c1[],
real fit 
)

Fit an autocorrelation function to a pre-defined functional form.

Todo:
check parameters
Parameters
[in]ncorr
[in]fitfnFitting function (0 .. effnNR)
[in]oenvOutput formatting information
[in]bVerboseShould the routine write to console?
[in]tbeginfitStarting time for fitting
[in]tendfitEnding time for fitting
[in]dtThe time step
[in]c1The data points
[in,out]fitThe fitting parameters
Returns
the integral over the autocorrelation function?
double fit_function ( int  eFitFn,
const double  parm[],
double  x 
)

Returns the value of fit function eFitFn at x.

Parameters
[in]eFitFnthe index to the fitting function (0 .. effnNR)
[in]parmArray of parameters, the length of which depends on eFitFn
[in]xThe value of x
Returns
the value of the fit
int sffn2effn ( const char **  sffn)

Returns corresponding to the selected enum option in sffn.

Parameters
[in]sffnTwo dimensional string array coming from parse_common_args
Returns
the ffn enum