NumeRe v1.1.4
NumeRe: Framework für Numerische Rechnungen
fitcontroller.hpp
Go to the documentation of this file.
1/*****************************************************************************
2 NumeRe: Framework fuer Numerische Rechnungen
3 Copyright (C) 2016 Erik Haenel et al.
4
5 This program is free software: you can redistribute it and/or modify
6 it under the terms of the GNU General Public License as published by
7 the Free Software Foundation, either version 3 of the License, or
8 (at your option) any later version.
9
10 This program is distributed in the hope that it will be useful,
11 but WITHOUT ANY WARRANTY; without even the implied warranty of
12 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 GNU General Public License for more details.
14
15 You should have received a copy of the GNU General Public License
16 along with this program. If not, see <http://www.gnu.org/licenses/>.
17******************************************************************************/
18
19
20
21
22
23// CLASS Fitcontroller
24
25#ifndef FITCONTROLLER_HPP
26#define FITCONTROLLER_HPP
27
28#include <iostream>
29#include <vector>
30#include <string>
31#include <gsl/gsl_multifit_nlin.h>
32#include <gsl/gsl_blas.h>
33#include <gsl/gsl_vector.h>
34
35#include "../ParserLib/muParser.h"
36#include "../utils/tools.hpp"
37#include "../ui/error.hpp"
38
39using namespace mu;
40
44typedef std::vector<double> FitVector;
45
49typedef std::vector<std::vector<double>> FitMatrix;
50
51
56struct FitData
57{
63 double dPrecision;
64};
65
66
73{
74 private:
75 static int fitfunction(const gsl_vector* params, void* data, gsl_vector* fvals);
76 static int fitjacobian(const gsl_vector* params, void* data, gsl_matrix* Jac);
77 static int fitfuncjac(const gsl_vector* params, void* data, gsl_vector* fvals, gsl_matrix* Jac);
78 static int fitfunctionrestricted(const gsl_vector* params, void* data, gsl_vector* fvals);
79 static int fitjacobianrestricted(const gsl_vector* params, void* data, gsl_matrix* Jac);
80 static int fitfuncjacrestricted(const gsl_vector* params, void* data, gsl_vector* fvals, gsl_matrix* Jac);
81 static double evalRestrictions(const mu::value_type* v, int nVals);
83 double dChiSqr;
84 std::string sExpr;
89
90 bool fitctrl(const std::string& __sExpr, const std::string& __sRestrictions, FitData& _fData, double __dPrecision, int nMaxIterations);
91 static void removeNANVals(gsl_vector* fvals, unsigned int nSize);
92 static void removeNANVals(gsl_matrix* Jac, unsigned int nLines, unsigned int nCols);
93
94 public:
96 static int nDimensions;
98
100 Fitcontroller(Parser* _parser);
102
110 inline void setParser(mu::Parser* _parser)
111 {
112 _fitParser = _parser;
113 mu::varmap_type mVars = _parser->GetVar();
114 xvar = mVars["x"];
115 yvar = mVars["y"];
116 zvar = mVars["z"];
117 return;
118 }
119
133 inline 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)
134 {
135 FitVector vy_w(vy.size(), 0.0);
136 return fit(vx,vy, vy_w, __sExpr, __sRestrictions, mParamsMap, __dPrecision, nMaxIterations);
137 }
138
139 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);
140
155 inline 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)
156 {
157 FitMatrix vz_w(vz.size(), std::vector<double>(vz[0].size(), 0.0));
158 return fit(vx, vy, vz, vz_w, __sExpr, __sRestrictions, mParamsMap, __dPrecision, nMaxIterations);
159 }
160
161 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);
162
170 inline double getFitChi() const
171 {return dChiSqr;}
172
180 inline int getIterations() const
181 {return nIterations;}
182
183 std::string getFitFunction();
185};
186
187
188#endif // FITCONTROLLER_HPP
189
This class contains the internal fit logic and the interface to the GSL fitting module.
static mu::value_type * xvar
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.
static int fitfunctionrestricted(const gsl_vector *params, void *data, gsl_vector *fvals)
Fit function respecting additional restrictions.
static mu::value_type * yvar
static mu::value_type * zvar
void setParser(mu::Parser *_parser)
Register the parser.
static mu::Parser * _fitParser
static double evalRestrictions(const mu::value_type *v, int nVals)
Evaluate additional restrictions.
int getIterations() const
Return the number of needed iterations.
static mu::varmap_type mParams
FitMatrix vCovarianceMatrix
FitMatrix getCovarianceMatrix() const
Returns the covariance matrix of the performed fit.
static void removeNANVals(gsl_vector *fvals, unsigned int nSize)
Replace NaNs with some very large value to avoid converging towards NaNs.
static int fitfuncjac(const gsl_vector *params, void *data, gsl_vector *fvals, gsl_matrix *Jac)
Combination of fit function and the corresponding jacobian.
static int fitjacobian(const gsl_vector *params, void *data, gsl_matrix *Jac)
Create the jacobian matrix without using restrictions.
static int fitfunction(const gsl_vector *params, void *data, gsl_vector *fvals)
Fit function for the usual case without any restrictions.
double getFitChi() const
Return the weighted sum of the residuals.
std::string sExpr
static int fitjacobianrestricted(const gsl_vector *params, void *data, gsl_matrix *Jac)
Jacobian matrix respecting additional restrictions.
static int nDimensions
~Fitcontroller()
Destructor.
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.
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.
Fitcontroller()
Default constructor.
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.
std::string getFitFunction()
Returns the fit function, where each parameter is replaced with its value converted to a string.
const varmap_type & GetVar() const
Return a map containing the used variables only.
Mathematical expressions parser.
Definition: muParser.h:51
std::vector< std::vector< double > > FitMatrix
Typedef for readability.
std::vector< double > FitVector
Typedef for readability.
Namespace for mathematical applications.
Definition: muParser.cpp:53
MUP_BASETYPE value_type
The numeric datatype used by the parser.
Definition: muParserDef.h:251
std::map< string_type, value_type * > varmap_type
Type used for storing variables.
Definition: muParserDef.h:273
Defines the datapoints, which shall be fitted.
FitMatrix vz_w
FitVector vy
double dPrecision
FitMatrix vz
FitVector vx
FitVector vy_w