NumeRe v1.1.4
NumeRe: Framework für Numerische Rechnungen
|
This class contains the internal fit logic and the interface to the GSL fitting module. More...
#include <fitcontroller.hpp>
Public Member Functions | |
Fitcontroller () | |
Default constructor. More... | |
Fitcontroller (Parser *_parser) | |
Constructor taking a pointer to the current parser. More... | |
~Fitcontroller () | |
Destructor. More... | |
void | setParser (mu::Parser *_parser) |
Register the parser. More... | |
bool | fit (FitVector &vx, FitVector &vy, const std::string &__sExpr, const std::string &__sRestrictions, mu::varmap_type &mParamsMap, double __dPrecision=1e-4, int nMaxIterations=500) |
Calculate an unweighted 1D-fit. More... | |
bool | fit (FitVector &vx, FitVector &vy, FitVector &vy_w, const std::string &__sExpr, const std::string &__sRestrictions, mu::varmap_type &mParamsMap, double __dPrecision=1e-4, int nMaxIterations=500) |
Calculate a weighted 1D-fit. More... | |
bool | fit (FitVector &vx, FitVector &vy, FitMatrix &vz, const std::string &__sExpr, const std::string &__sRestrictions, mu::varmap_type &mParamsMap, double __dPrecision=1e-4, int nMaxIterations=500) |
Calculate an unweighted 2D-fit. More... | |
bool | fit (FitVector &vx, FitVector &vy, FitMatrix &vz, FitMatrix &vz_w, const std::string &__sExpr, const std::string &__sRestrictions, mu::varmap_type &mParamsMap, double __dPrecision=1e-4, int nMaxIterations=500) |
Calculate a weighted 2D-fit. More... | |
double | getFitChi () const |
Return the weighted sum of the residuals. More... | |
int | getIterations () const |
Return the number of needed iterations. More... | |
std::string | getFitFunction () |
Returns the fit function, where each parameter is replaced with its value converted to a string. More... | |
FitMatrix | getCovarianceMatrix () const |
Returns the covariance matrix of the performed fit. More... | |
Static Public Attributes | |
static mu::Parser * | _fitParser = 0 |
static int | nDimensions = 0 |
static mu::varmap_type | mParams |
Private Member Functions | |
bool | fitctrl (const std::string &__sExpr, const std::string &__sRestrictions, FitData &_fData, double __dPrecision, int nMaxIterations) |
This is the central fitting function using the data passed from the interface functions. More... | |
Static Private Member Functions | |
static int | fitfunction (const gsl_vector *params, void *data, gsl_vector *fvals) |
Fit function for the usual case without any restrictions. More... | |
static int | fitjacobian (const gsl_vector *params, void *data, gsl_matrix *Jac) |
Create the jacobian matrix without using restrictions. More... | |
static int | fitfuncjac (const gsl_vector *params, void *data, gsl_vector *fvals, gsl_matrix *Jac) |
Combination of fit function and the corresponding jacobian. More... | |
static int | fitfunctionrestricted (const gsl_vector *params, void *data, gsl_vector *fvals) |
Fit function respecting additional restrictions. More... | |
static int | fitjacobianrestricted (const gsl_vector *params, void *data, gsl_matrix *Jac) |
Jacobian matrix respecting additional restrictions. More... | |
static int | fitfuncjacrestricted (const gsl_vector *params, void *data, gsl_vector *fvals, gsl_matrix *Jac) |
Combination of fit function and the corresponding jacobian respecting additional restrictions. More... | |
static double | evalRestrictions (const mu::value_type *v, int nVals) |
Evaluate additional restrictions. More... | |
static void | removeNANVals (gsl_vector *fvals, unsigned int nSize) |
Replace NaNs with some very large value to avoid converging towards NaNs. More... | |
static void | removeNANVals (gsl_matrix *Jac, unsigned int nLines, unsigned int nCols) |
Replace NaNs with some very large value to avoid converging towards NaNs. More... | |
Private Attributes | |
int | nIterations |
double | dChiSqr |
std::string | sExpr |
FitMatrix | vCovarianceMatrix |
Static Private Attributes | |
static mu::value_type * | xvar = 0 |
static mu::value_type * | yvar = 0 |
static mu::value_type * | zvar = 0 |
This class contains the internal fit logic and the interface to the GSL fitting module.
Definition at line 72 of file fitcontroller.hpp.
Fitcontroller::Fitcontroller | ( | ) |
Default constructor.
Definition at line 35 of file fitcontroller.cpp.
References dChiSqr, nIterations, sExpr, xvar, yvar, and zvar.
Fitcontroller::Fitcontroller | ( | Parser * | _parser | ) |
Constructor taking a pointer to the current parser.
_parser | Parser* |
Definition at line 54 of file fitcontroller.cpp.
References _fitParser, mu::ParserBase::GetVar(), xvar, yvar, and zvar.
Fitcontroller::~Fitcontroller | ( | ) |
Destructor.
Definition at line 67 of file fitcontroller.cpp.
References _fitParser, xvar, yvar, and zvar.
|
staticprivate |
Evaluate additional restrictions.
v | const value_type* |
nVals | int |
Definition at line 502 of file fitcontroller.cpp.
References mu::real().
Referenced by fitctrl(), fitfunctionrestricted(), and fitjacobianrestricted().
|
inline |
Calculate an unweighted 1D-fit.
vx | FitVector& |
vy | FitVector& |
__sExpr | const std::string& |
__sRestrictions | const std::string& |
mParamsMap | mu::varmap_type& |
__dPrecision | double |
nMaxIterations | int |
Definition at line 133 of file fitcontroller.hpp.
References fit().
Referenced by applyFitAlgorithm(), calculateChiMap(), and fit().
|
inline |
Calculate an unweighted 2D-fit.
vx | FitVector& |
vy | FitVector& |
vz | FitMatrix& |
__sExpr | const std::string& |
__sRestrictions | const std::string& |
mParamsMap | mu::varmap_type& |
__dPrecision | double |
nMaxIterations | int |
Definition at line 155 of file fitcontroller.hpp.
References fit().
bool Fitcontroller::fit | ( | FitVector & | vx, |
FitVector & | vy, | ||
FitMatrix & | vz, | ||
FitMatrix & | vz_w, | ||
const std::string & | __sExpr, | ||
const std::string & | __sRestrictions, | ||
mu::varmap_type & | mParamsMap, | ||
double | __dPrecision = 1e-4 , |
||
int | nMaxIterations = 500 |
||
) |
Calculate a weighted 2D-fit.
vx | FitVector& |
vy | FitVector& |
vz | FitMatrix& |
vz_w | FitMatrix& |
__sExpr | const string& |
__sRestrictions | const string& |
mParamsMap | mu::varmap_type& |
__dPrecision | double |
nMaxIterations | int |
Definition at line 844 of file fitcontroller.cpp.
References fitctrl(), mParams, FitData::vx, FitData::vy, FitData::vz, and FitData::vz_w.
bool Fitcontroller::fit | ( | FitVector & | vx, |
FitVector & | vy, | ||
FitVector & | vy_w, | ||
const std::string & | __sExpr, | ||
const std::string & | __sRestrictions, | ||
mu::varmap_type & | mParamsMap, | ||
double | __dPrecision = 1e-4 , |
||
int | nMaxIterations = 500 |
||
) |
Calculate a weighted 1D-fit.
vx | FitVector& |
vy | FitVector& |
vy_w | FitVector& |
__sExpr | const string& |
__sRestrictions | const string& |
mParamsMap | mu::varmap_type& |
__dPrecision | double |
nMaxIterations | int |
Definition at line 817 of file fitcontroller.cpp.
References fitctrl(), mParams, FitData::vx, FitData::vy, and FitData::vy_w.
|
private |
This is the central fitting function using the data passed from the interface functions.
__sExpr | const string& |
__sRestrictions | const string& |
_fData | FitData& |
__dPrecision | double |
nMaxIterations | int |
Definition at line 615 of file fitcontroller.cpp.
References _fitParser, dChiSqr, FitData::dPrecision, mu::ParserBase::Eval(), evalRestrictions(), fitfuncjac(), fitfuncjacrestricted(), fitfunction(), fitfunctionrestricted(), fitjacobian(), fitjacobianrestricted(), func, NumeReKernel::GetAsyncCancelState(), SyntaxError::invalid_position, mu::isinf(), mu::isnan(), mParams, nIterations, SyntaxError::PROCESS_ABORTED_BY_USER, mu::ParserBase::SetExpr(), sExpr, vCovarianceMatrix, FitData::vx, FitData::vy, FitData::vy_w, FitData::vz, and FitData::vz_w.
Referenced by fit().
|
staticprivate |
Combination of fit function and the corresponding jacobian.
params | const gsl_vector* |
data | void* |
fvals | gsl_vector* |
Jac | gsl_matrix* |
Definition at line 272 of file fitcontroller.cpp.
References fitfunction(), and fitjacobian().
Referenced by fitctrl().
|
staticprivate |
Combination of fit function and the corresponding jacobian respecting additional restrictions.
params | const gsl_vector* |
data | void* |
fvals | gsl_vector* |
Jac | gsl_matrix* |
Definition at line 486 of file fitcontroller.cpp.
References fitfunctionrestricted(), and fitjacobianrestricted().
Referenced by fitctrl().
|
staticprivate |
Fit function for the usual case without any restrictions.
params | const gsl_vector* |
data | void* |
fvals | gsl_vector* |
Definition at line 87 of file fitcontroller.cpp.
References _fitParser, mu::ParserBase::Eval(), mu::isinf(), mu::isnan(), mParams, removeNANVals(), FitData::vx, FitData::vy, FitData::vy_w, FitData::vz, FitData::vz_w, xvar, and yvar.
Referenced by fitctrl(), and fitfuncjac().
|
staticprivate |
Fit function respecting additional restrictions.
params | const gsl_vector* |
data | void* |
fvals | gsl_vector* |
Definition at line 290 of file fitcontroller.cpp.
References _fitParser, mu::ParserBase::Eval(), evalRestrictions(), mu::isinf(), mu::isnan(), mParams, removeNANVals(), FitData::vx, FitData::vy, FitData::vy_w, FitData::vz, FitData::vz_w, xvar, and yvar.
Referenced by fitctrl(), and fitfuncjacrestricted().
|
staticprivate |
Create the jacobian matrix without using restrictions.
params | const gsl_vector* |
data | void* |
Jac | gsl_matrix* |
Definition at line 152 of file fitcontroller.cpp.
References _fitParser, FitData::dPrecision, mu::ParserBase::Eval(), mu::isinf(), mu::isnan(), mParams, removeNANVals(), FitData::vx, FitData::vy, FitData::vy_w, FitData::vz, FitData::vz_w, xvar, and yvar.
Referenced by fitctrl(), and fitfuncjac().
|
staticprivate |
Jacobian matrix respecting additional restrictions.
params | const gsl_vector* |
data | void* |
Jac | gsl_matrix* |
Definition at line 359 of file fitcontroller.cpp.
References _fitParser, FitData::dPrecision, mu::ParserBase::Eval(), evalRestrictions(), mu::isinf(), mu::isnan(), mParams, removeNANVals(), FitData::vx, FitData::vy, FitData::vy_w, FitData::vz, FitData::vz_w, xvar, and yvar.
Referenced by fitctrl(), and fitfuncjacrestricted().
FitMatrix Fitcontroller::getCovarianceMatrix | ( | ) | const |
Returns the covariance matrix of the performed fit.
Definition at line 918 of file fitcontroller.cpp.
References vCovarianceMatrix.
Referenced by fitDataSet().
|
inline |
Return the weighted sum of the residuals.
Definition at line 170 of file fitcontroller.hpp.
References dChiSqr.
Referenced by calculateChiMap(), createTeXExport(), fitDataSet(), and getFitOptionsTable().
string Fitcontroller::getFitFunction | ( | ) |
Returns the fit function, where each parameter is replaced with its value converted to a string.
Definition at line 865 of file fitcontroller.cpp.
References checkDelimiter(), mParams, sExpr, StripSpaces(), and toString().
Referenced by fitDataSet().
|
inline |
Return the number of needed iterations.
Definition at line 180 of file fitcontroller.hpp.
References nIterations.
Referenced by createTeXExport(), getFitAnalysis(), and getFitOptionsTable().
|
staticprivate |
Replace NaNs with some very large value to avoid converging towards NaNs.
Jac | gsl_matrix* |
nLines | unsigned int |
nCols | unsigned int |
Definition at line 568 of file fitcontroller.cpp.
References mu::isinf(), and mu::isnan().
|
staticprivate |
Replace NaNs with some very large value to avoid converging towards NaNs.
fvals | gsl_vector* |
nSize | unsigned int |
Definition at line 530 of file fitcontroller.cpp.
References mu::isinf(), and mu::isnan().
Referenced by fitfunction(), fitfunctionrestricted(), fitjacobian(), and fitjacobianrestricted().
|
inline |
Register the parser.
_parser | mu::Parser* |
Definition at line 110 of file fitcontroller.hpp.
References _fitParser, mu::ParserBase::GetVar(), xvar, yvar, and zvar.
|
static |
Definition at line 95 of file fitcontroller.hpp.
Referenced by Fitcontroller(), fitctrl(), fitfunction(), fitfunctionrestricted(), fitjacobian(), fitjacobianrestricted(), setParser(), and ~Fitcontroller().
|
private |
Definition at line 83 of file fitcontroller.hpp.
Referenced by Fitcontroller(), fitctrl(), and getFitChi().
|
static |
Definition at line 97 of file fitcontroller.hpp.
Referenced by fit(), fitctrl(), fitfunction(), fitfunctionrestricted(), fitjacobian(), fitjacobianrestricted(), and getFitFunction().
|
static |
Definition at line 96 of file fitcontroller.hpp.
|
private |
Definition at line 82 of file fitcontroller.hpp.
Referenced by Fitcontroller(), fitctrl(), and getIterations().
|
private |
Definition at line 84 of file fitcontroller.hpp.
Referenced by Fitcontroller(), fitctrl(), and getFitFunction().
|
private |
Definition at line 85 of file fitcontroller.hpp.
Referenced by fitctrl(), and getCovarianceMatrix().
|
staticprivate |
Definition at line 86 of file fitcontroller.hpp.
Referenced by Fitcontroller(), fitfunction(), fitfunctionrestricted(), fitjacobian(), fitjacobianrestricted(), setParser(), and ~Fitcontroller().
|
staticprivate |
Definition at line 87 of file fitcontroller.hpp.
Referenced by Fitcontroller(), fitfunction(), fitfunctionrestricted(), fitjacobian(), fitjacobianrestricted(), setParser(), and ~Fitcontroller().
|
staticprivate |
Definition at line 88 of file fitcontroller.hpp.
Referenced by Fitcontroller(), setParser(), and ~Fitcontroller().